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 (!methodDict_.empty())
93  {
94  if (methodDict_.readIfPresent("method", method))
95  {
96  if (method != "recursive" && method != "kWay")
97  {
98  FatalIOErrorInFunction(methodDict_)
99  << "Method " << method
100  << " should be 'recursive' or 'kWay'"
101  << exit(FatalIOError);
102  }
103 
104  Info<< "metis : Using Metis method " << method
105  << nl << endl;
106  }
107 
108  if (methodDict_.readIfPresent("options", options))
109  {
110  if (options.size() != METIS_NOPTIONS)
111  {
112  FatalIOErrorInFunction(methodDict_)
113  << "Number of options should be " << METIS_NOPTIONS
114  << " found " << options
115  << exit(FatalIOError);
116  }
117 
118  Info<< "Using Metis options " << options << nl << endl;
119  }
120 
121  if (methodDict_.readIfPresent("processorWeights", processorWeights))
122  {
123  processorWeights /= sum(processorWeights);
124 
125  if (processorWeights.size() != nProcessors_)
126  {
127  FatalIOErrorInFunction(methodDict_)
128  << "Number of processor weights "
129  << processorWeights.size()
130  << " does not equal number of domains " << nProcessors_
131  << exit(FatalIOError);
132  }
133  }
134  }
135 
136  label ncon = 1;
137  label nProcs = nProcessors_;
138 
139  // Output: cell -> processor addressing
140  decomp.setSize(numCells);
141 
142  // Output: number of cut edges
143  label edgeCut = 0;
144 
145  if (method == "recursive")
146  {
147  METIS_PartGraphRecursive
148  (
149  &numCells, // num vertices in graph
150  &ncon, // num balancing constraints
151  const_cast<labelList&>(xadj).begin(), // indexing into adjncy
152  const_cast<labelList&>(adjncy).begin(), // neighbour info
153  intWeights.begin(), // vertexweights
154  nullptr, // vsize: total communication vol
155  faceWeights.begin(), // edgeweights
156  &nProcs, // nParts
157  processorWeights.begin(), // tpwgts
158  nullptr, // ubvec: processor imbalance (default)
159  options.begin(),
160  &edgeCut,
161  decomp.begin()
162  );
163  }
164  else
165  {
166  METIS_PartGraphKway
167  (
168  &numCells, // num vertices in graph
169  &ncon, // num balancing constraints
170  const_cast<labelList&>(xadj).begin(), // indexing into adjncy
171  const_cast<labelList&>(adjncy).begin(), // neighbour info
172  intWeights.begin(), // vertexweights
173  nullptr, // vsize: total communication vol
174  faceWeights.begin(), // edgeweights
175  &nProcs, // nParts
176  processorWeights.begin(), // tpwgts
177  nullptr, // ubvec: processor imbalance (default)
178  options.begin(),
179  &edgeCut,
180  decomp.begin()
181  );
182  }
183 
184  return edgeCut;
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
189 
191 (
192  const dictionary& decompositionDict,
193  const dictionary& methodDict
194 )
195 :
196  decompositionMethod(decompositionDict),
197  methodDict_(methodDict)
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,
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:433
Macros for easy insertion into run-time selection tables.
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
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 &decompositionDict, const dictionary &methodDict)
Construct given the decomposition dictionary.
Definition: dummyMetis.C:73
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
label nCells() const
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#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:258
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
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
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:267