metisDecomp.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "metisDecomp.H"
28 #include "floatScalar.H"
29 #include "Time.H"
30 
31 extern "C"
32 {
33  #define OMPI_SKIP_MPICXX
34  #include "metis.h"
35 }
36 
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(metisDecomp, 0);
43  addToRunTimeSelectionTable(decompositionMethod, metisDecomp, dictionary);
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
50 (
51  const List<label>& adjncy,
52  const List<label>& xadj,
53  const scalarField& cWeights,
54 
55  List<label>& finalDecomp
56 )
57 {
58  // Method of decomposition
59  // recursive: multi-level recursive bisection (default)
60  // k-way: multi-level k-way
61  word method("recursive");
62 
63  label numCells = xadj.size()-1;
64 
65  // decomposition options
66  List<label> options(METIS_NOPTIONS);
67  METIS_SetDefaultOptions(options.begin());
68 
69  // processor weights initialised with no size, only used if specified in
70  // a file
71  Field<floatScalar> processorWeights;
72 
73  // cell weights (so on the vertices of the dual)
74  List<label> cellWeights;
75 
76  // face weights (so on the edges of the dual)
77  List<label> faceWeights;
78 
79 
80  // Check for externally provided cellweights and if so initialise weights
81  scalar minWeights = gMin(cWeights);
82  if (cWeights.size() > 0)
83  {
84  if (minWeights <= 0)
85  {
87  << "Illegal minimum weight " << minWeights
88  << endl;
89  }
90 
91  if (cWeights.size() != numCells)
92  {
94  << "Number of cell weights " << cWeights.size()
95  << " does not equal number of cells " << numCells
96  << exit(FatalError);
97  }
98  // Convert to integers.
99  cellWeights.setSize(cWeights.size());
100  forAll(cellWeights, i)
101  {
102  cellWeights[i] = int(cWeights[i]/minWeights);
103  }
104  }
105 
106 
107  // Check for user supplied weights and decomp options
108  if (decompositionDict_.found("metisCoeffs"))
109  {
110  const dictionary& metisCoeffs =
111  decompositionDict_.subDict("metisCoeffs");
112  word weightsFile;
113 
114  if (metisCoeffs.readIfPresent("method", method))
115  {
116  if (method != "recursive" && method != "k-way")
117  {
119  << "Method " << method << " in metisCoeffs in dictionary : "
121  << " should be 'recursive' or 'k-way'"
122  << exit(FatalError);
123  }
124 
125  Info<< "metisDecomp : Using Metis method " << method
126  << nl << endl;
127  }
128 
129  if (metisCoeffs.readIfPresent("options", options))
130  {
131  if (options.size() != METIS_NOPTIONS)
132  {
134  << "Number of options in metisCoeffs in dictionary : "
136  << " should be " << METIS_NOPTIONS
137  << exit(FatalError);
138  }
139 
140  Info<< "metisDecomp : Using Metis options " << options
141  << nl << endl;
142  }
143 
144  if (metisCoeffs.readIfPresent("processorWeights", processorWeights))
145  {
146  processorWeights /= sum(processorWeights);
147 
148  if (processorWeights.size() != nProcessors_)
149  {
151  << "Number of processor weights "
152  << processorWeights.size()
153  << " does not equal number of domains " << nProcessors_
154  << exit(FatalError);
155  }
156  }
157 
158  //if (metisCoeffs.readIfPresent("cellWeightsFile", weightsFile))
159  //{
160  // Info<< "metisDecomp : Using cell-based weights." << endl;
161  //
162  // IOList<label> cellIOWeights
163  // (
164  // IOobject
165  // (
166  // weightsFile,
167  // mesh_.time().timeName(),
168  // mesh_,
169  // IOobject::MUST_READ,
170  // IOobject::AUTO_WRITE
171  // )
172  // );
173  // cellWeights.transfer(cellIOWeights);
174  //
175  // if (cellWeights.size() != xadj.size()-1)
176  // {
177  // FatalErrorInFunction
178  // << "Number of cell weights " << cellWeights.size()
179  // << " does not equal number of cells " << xadj.size()-1
180  // << exit(FatalError);
181  // }
182  //}
183  }
184 
185  label ncon = 1;
186 
187  label nProcs = nProcessors_;
188 
189  // output: cell -> processor addressing
190  finalDecomp.setSize(numCells);
191 
192  // output: number of cut edges
193  label edgeCut = 0;
194 
195  if (method == "recursive")
196  {
197  METIS_PartGraphRecursive
198  (
199  &numCells, // num vertices in graph
200  &ncon, // num balancing constraints
201  const_cast<List<label>&>(xadj).begin(), // indexing into adjncy
202  const_cast<List<label>&>(adjncy).begin(), // neighbour info
203  cellWeights.begin(),// vertexweights
204  NULL, // vsize: total communication vol
205  faceWeights.begin(),// edgeweights
206  &nProcs, // nParts
207  processorWeights.begin(), // tpwgts
208  NULL, // ubvec: processor imbalance (default)
209  options.begin(),
210  &edgeCut,
211  finalDecomp.begin()
212  );
213  }
214  else
215  {
216  METIS_PartGraphKway
217  (
218  &numCells, // num vertices in graph
219  &ncon, // num balancing constraints
220  const_cast<List<label>&>(xadj).begin(), // indexing into adjncy
221  const_cast<List<label>&>(adjncy).begin(), // neighbour info
222  cellWeights.begin(),// vertexweights
223  NULL, // vsize: total communication vol
224  faceWeights.begin(),// edgeweights
225  &nProcs, // nParts
226  processorWeights.begin(), // tpwgts
227  NULL, // ubvec: processor imbalance (default)
228  options.begin(),
229  &edgeCut,
230  finalDecomp.begin()
231  );
232  }
233 
234  return edgeCut;
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
239 
240 Foam::metisDecomp::metisDecomp(const dictionary& decompositionDict)
241 :
242  decompositionMethod(decompositionDict)
243 {}
244 
245 
246 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
247 
249 (
250  const polyMesh& mesh,
251  const pointField& points,
252  const scalarField& pointWeights
253 )
254 {
255  if (points.size() != mesh.nCells())
256  {
258  << "Can use this decomposition method only for the whole mesh"
259  << endl
260  << "and supply one coordinate (cellCentre) for every cell." << endl
261  << "The number of coordinates " << points.size() << endl
262  << "The number of cells in the mesh " << mesh.nCells()
263  << exit(FatalError);
264  }
265 
266  CompactListList<label> cellCells;
268  (
269  mesh,
270  identity(mesh.nCells()),
271  mesh.nCells(),
272  false,
273  cellCells
274  );
275 
276  // Decompose using default weights
277  labelList decomp;
278  decompose(cellCells.m(), cellCells.offsets(), pointWeights, decomp);
279 
280  return decomp;
281 }
282 
283 
285 (
286  const polyMesh& mesh,
287  const labelList& agglom,
288  const pointField& agglomPoints,
289  const scalarField& agglomWeights
290 )
291 {
292  if (agglom.size() != mesh.nCells())
293  {
295  << "Size of cell-to-coarse map " << agglom.size()
296  << " differs from number of cells in mesh " << mesh.nCells()
297  << exit(FatalError);
298  }
299 
300  // Make Metis CSR (Compressed Storage Format) storage
301  // adjncy : contains neighbours (= edges in graph)
302  // xadj(celli) : start of information in adjncy for celli
303 
304  CompactListList<label> cellCells;
305  calcCellCells(mesh, agglom, agglomPoints.size(), false, cellCells);
306 
307  // Decompose using default weights
308  labelList finalDecomp;
309  decompose(cellCells.m(), cellCells.offsets(), agglomWeights, finalDecomp);
310 
311 
312  // Rework back into decomposition for original mesh
313  labelList fineDistribution(agglom.size());
314 
315  forAll(fineDistribution, i)
316  {
317  fineDistribution[i] = finalDecomp[agglom[i]];
318  }
319 
320  return finalDecomp;
321 }
322 
323 
325 (
326  const labelListList& globalCellCells,
327  const pointField& cellCentres,
328  const scalarField& cellWeights
329 )
330 {
331  if (cellCentres.size() != globalCellCells.size())
332  {
334  << "Inconsistent number of cells (" << globalCellCells.size()
335  << ") and number of cell centres (" << cellCentres.size()
336  << ")." << exit(FatalError);
337  }
338 
339 
340  // Make Metis CSR (Compressed Storage Format) storage
341  // adjncy : contains neighbours (= edges in graph)
342  // xadj(celli) : start of information in adjncy for celli
343 
344  CompactListList<label> cellCells(globalCellCells);
345 
346  // Decompose using default weights
347  labelList decomp;
348  decompose(cellCells.m(), cellCells.offsets(), cellWeights, decomp);
349 
350  return decomp;
351 }
352 
353 
354 // ************************************************************************* //
metisDecomp.H
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::decompositionMethod::calcCellCells
static void calcCellCells(const polyMesh &mesh, const labelList &agglom, const label nLocalCoarse, const bool global, CompactListList< label > &cellCells)
Helper: determine (local or global) cellCells from mesh.
Definition: decompositionMethod.C:289
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::dictionaryName::name
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:103
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
List::size
int size() const
Definition: ListI.H:83
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::dictionary::found
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:304
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::decompositionMethod::decompositionDict_
const dictionary & decompositionDict_
Definition: decompositionMethod.H:55
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
Foam::FatalError
error FatalError
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::metisDecomp::metisDecomp
metisDecomp(const metisDecomp &)
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
Foam::decompositionMethod::nProcessors_
label nProcessors_
Definition: decompositionMethod.H:56
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
scalarField
volScalarField scalarField(fieldObject, mesh)
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
floatScalar.H
Foam::metisDecomp::decompose
label decompose(const List< label > &adjncy, const List< label > &xadj, const scalarField &cellWeights, List< label > &finalDecomp)
Call Metis with options from dictionary.
Definition: dummyMetisDecomp.C:59
List
Definition: Test.C:19
Foam::gMin
Type gMin(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:563
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259