metis.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 "metis.H"
27 #include "Time.H"
29 
30 extern "C"
31 {
32  #define OMPI_SKIP_MPICXX
33  #include "metis.h"
34 }
35 
36 
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace decompositionMethods
42 {
43  defineTypeNameAndDebug(metis, 0);
44  addToRunTimeSelectionTable(decompositionMethod, metis, decomposer);
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 (
53  const labelList& adjncy,
54  const labelList& xadj,
55  const scalarField& cellWeights,
56  labelList& decomp
57 )
58 {
59  // Method of decomposition
60  // recursive: multi-level recursive bisection (default)
61  // kWay: multi-level k-way
62  word method("recursive");
63 
64  label numCells = xadj.size()-1;
65 
66  // Decomposition options
67  labelList options(METIS_NOPTIONS);
68  METIS_SetDefaultOptions(options.begin());
69 
70  // Processor weights initialised with no size, only used if specified in
71  // a file
72  Field<real_t> processorWeights;
73 
74  if (cellWeights.size() > 0 && cellWeights.size() != numCells)
75  {
77  << "Number of cell weights " << cellWeights.size()
78  << " does not equal number of cells " << numCells
79  << exit(FatalError);
80  }
81 
82  label nWeights = 1;
83 
84  // Cell weights (so on the vertices of the dual)
85  labelList intWeights(scaleWeights(cellWeights, nWeights, false));
86 
87  // Face weights (so on the edges of the dual)
89 
90 
91  // Check for user supplied weights and decomp options
92  if (decompositionDict_.found("metisCoeffs"))
93  {
94  const dictionary& metisCoeffs =
95  decompositionDict_.subDict("metisCoeffs");
96 
97  if (metisCoeffs.readIfPresent("method", method))
98  {
99  if (method != "recursive" && method != "kWay")
100  {
101  FatalIOErrorInFunction(metisCoeffs)
102  << "Method " << method << " in metisCoeffs in dictionary : "
104  << " should be 'recursive' or 'kWay'"
105  << exit(FatalIOError);
106  }
107 
108  Info<< "metis : Using Metis method " << method
109  << nl << endl;
110  }
111 
112  if (metisCoeffs.readIfPresent("options", options))
113  {
114  if (options.size() != METIS_NOPTIONS)
115  {
116  FatalIOErrorInFunction(metisCoeffs)
117  << "Number of options in metisCoeffs dictionary : "
119  << " should be " << METIS_NOPTIONS << " found " << options
120  << exit(FatalIOError);
121  }
122 
123  Info<< "Using Metis options " << options << nl << endl;
124  }
125 
126  if (metisCoeffs.readIfPresent("processorWeights", processorWeights))
127  {
128  processorWeights /= sum(processorWeights);
129 
130  if (processorWeights.size() != nProcessors_)
131  {
132  FatalIOErrorInFunction(metisCoeffs)
133  << "Number of processor weights "
134  << processorWeights.size()
135  << " does not equal number of domains " << nProcessors_
136  << exit(FatalIOError);
137  }
138  }
139  }
140 
141  label ncon = 1;
142  label nProcs = nProcessors_;
143 
144  // Output: cell -> processor addressing
145  decomp.setSize(numCells);
146 
147  // Output: number of cut edges
148  label edgeCut = 0;
149 
150  if (method == "recursive")
151  {
152  METIS_PartGraphRecursive
153  (
154  &numCells, // num vertices in graph
155  &ncon, // num balancing constraints
156  const_cast<labelList&>(xadj).begin(), // indexing into adjncy
157  const_cast<labelList&>(adjncy).begin(), // neighbour info
158  intWeights.begin(), // vertexweights
159  nullptr, // vsize: total communication vol
160  faceWeights.begin(), // edgeweights
161  &nProcs, // nParts
162  processorWeights.begin(), // tpwgts
163  nullptr, // ubvec: processor imbalance (default)
164  options.begin(),
165  &edgeCut,
166  decomp.begin()
167  );
168  }
169  else
170  {
171  METIS_PartGraphKway
172  (
173  &numCells, // num vertices in graph
174  &ncon, // num balancing constraints
175  const_cast<labelList&>(xadj).begin(), // indexing into adjncy
176  const_cast<labelList&>(adjncy).begin(), // neighbour info
177  intWeights.begin(), // vertexweights
178  nullptr, // vsize: total communication vol
179  faceWeights.begin(), // edgeweights
180  &nProcs, // nParts
181  processorWeights.begin(), // tpwgts
182  nullptr, // ubvec: processor imbalance (default)
183  options.begin(),
184  &edgeCut,
185  decomp.begin()
186  );
187  }
188 
189  return edgeCut;
190 }
191 
192 
193 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
194 
195 Foam::decompositionMethods::metis::metis(const dictionary& decompositionDict)
196 :
197  decompositionMethod(decompositionDict)
198 {}
199 
200 
201 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
202 
204 (
205  const polyMesh& mesh,
206  const pointField& points,
207  const scalarField& pointWeights
208 )
209 {
210  if (points.size() != mesh.nCells())
211  {
213  << "Can use this decomposition method only for the whole mesh"
214  << endl
215  << "and supply one coordinate (cellCentre) for every cell." << endl
216  << "The number of coordinates " << points.size() << endl
217  << "The number of cells in the mesh " << mesh.nCells()
218  << exit(FatalError);
219  }
220 
221  checkWeights(points, pointWeights);
222 
223  // Make Metis CSR (Compressed Storage Format) storage
224  // adjncy : contains neighbours (= edges in graph)
225  // xadj(celli) : start of information in adjncy for celli
226  CompactListList<label> cellCells;
228  (
229  mesh,
230  identityMap(mesh.nCells()),
231  mesh.nCells(),
232  false,
233  cellCells
234  );
235 
236  // Decompose using default weights
237  labelList decomp;
238  decompose
239  (
240  cellCells.m(),
241  cellCells.offsets(),
242  pointWeights,
243  decomp
244  );
245 
246  return decomp;
247 }
248 
249 
251 (
252  const polyMesh& mesh,
253  const labelList& cellToRegion,
254  const pointField& regionPoints,
255  const scalarField& regionWeights
256 )
257 {
258  if (cellToRegion.size() != mesh.nCells())
259  {
261  << "Size of cell-to-coarse map " << cellToRegion.size()
262  << " differs from number of cells in mesh " << mesh.nCells()
263  << exit(FatalError);
264  }
265 
266  checkWeights(regionPoints, regionWeights);
267 
268  // Make Metis CSR (Compressed Storage Format) storage
269  // adjncy : contains neighbours (= edges in graph)
270  // xadj(celli) : start of information in adjncy for celli
271 
272  CompactListList<label> cellCells;
273  calcCellCells(mesh, cellToRegion, regionPoints.size(), false, cellCells);
274 
275  // Decompose using default weights
276  labelList decomp;
277  decompose(cellCells.m(), cellCells.offsets(), regionWeights, decomp);
278 
279 
280  // Rework back into decomposition for original mesh
281  labelList fineDistribution(cellToRegion.size());
282 
283  forAll(fineDistribution, i)
284  {
285  fineDistribution[i] = decomp[cellToRegion[i]];
286  }
287 
288  return fineDistribution;
289 }
290 
291 
293 (
294  const labelListList& globalCellCells,
295  const pointField& cellCentres,
296  const scalarField& cellWeights
297 )
298 {
299  if (cellCentres.size() != globalCellCells.size())
300  {
302  << "Inconsistent number of cells (" << globalCellCells.size()
303  << ") and number of cell centres (" << cellCentres.size()
304  << ")." << exit(FatalError);
305  }
306 
307  checkWeights(cellCentres, cellWeights);
308 
309  // Make Metis CSR (Compressed Storage Format) storage
310  // adjncy : contains neighbours (= edges in graph)
311  // xadj(celli) : start of information in adjncy for celli
312  CompactListList<label> cellCells(globalCellCells);
313 
314  // Decompose using default weights
315  labelList decomp;
316  decompose(cellCells.m(), cellCells.offsets(), cellWeights, decomp);
317 
318  return decomp;
319 }
320 
321 
322 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Macros for easy insertion into run-time selection tables.
iterator begin()
Return an iterator to begin traversing the UList.
Definition: UListI.H:216
label nWeights(const pointField &points, const scalarField &pointWeights) const
Return the number of weights per point.
label checkWeights(const pointField &points, const scalarField &pointWeights) const
Check the weights against the points.
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.
static labelList scaleWeights(const scalarField &weights, label &nWeights, const bool distributed=true)
Convert the given scalar weights to labels.
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Inherit decompose from decompositionMethod.
metis(const dictionary &)
Construct given the decomposition dictionary.
Definition: dummyMetis.C:73
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:111
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField scalarField(fieldObject, mesh)
const pointField & points
addToRunTimeSelectionTable(decompositionMethod, metis, decomposer)
tmp< scalarField > faceWeights(const polyMesh &mesh, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs)
Generate interpolation factors field.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
IOerror FatalIOError
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
static const char nl
Definition: Ostream.H:266