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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "metisDecomp.H"
31 #include "Time.H"
33 
34 // Probably not needed...
35 #ifndef MPICH_SKIP_MPICXX
36 #define MPICH_SKIP_MPICXX
37 #endif
38 #ifndef OMPI_SKIP_MPICXX
39 #define OMPI_SKIP_MPICXX
40 #endif
41 
42 #include "metis.h"
43 
44 // Provide a clear error message if we have a severe size mismatch
45 // Allow widening, but not narrowing
46 //
47 // Metis has an 'idx_t' type, but the IDXTYPEWIDTH define is perhaps
48 // more future-proof?
49 //#ifdef IDXTYPEWIDTH
50 //static_assert
51 //(
52 // sizeof(Foam::label) > (IDXTYPEWIDTH/8),
53 // "sizeof(Foam::label) > (IDXTYPEWIDTH/8), check your metis headers"
54 //);
55 //#endif
56 
57 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
58 
59 namespace Foam
60 {
61  defineTypeNameAndDebug(metisDecomp, 0);
63  (
64  decompositionMethod,
65  metisDecomp,
67  );
68 }
69 
70 
71 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
72 
74 (
75  const labelList& adjncy,
76  const labelList& xadj,
77  const List<scalar>& cWeights,
78  labelList& decomp
79 ) const
80 {
81  // Method of decomposition
82  // recursive: multi-level recursive bisection (default)
83  // k-way: multi-level k-way
84  word method("recursive");
85 
86  const dictionary* coeffsDictPtr = decompDict_.findDict("metisCoeffs");
87 
88  idx_t numCells = max(0, (xadj.size()-1));
89 
90  // Decomposition options
91  List<idx_t> options(METIS_NOPTIONS);
92  METIS_SetDefaultOptions(options.data());
93 
94  // Processor weights initialised with no size, only used if specified in
95  // a file
96  Field<real_t> procWeights;
97 
98  // Cell weights (so on the vertices of the dual)
99  List<idx_t> cellWeights;
100 
101  // Face weights (so on the edges of the dual)
102  List<idx_t> faceWeights;
103 
104  bool hasWeights = !cWeights.empty();
105 
106  // Note: min, not gMin since routine runs on master only.
107  const scalar minWeights = hasWeights ? min(cWeights) : scalar(1);
108 
109  if (minWeights <= 0)
110  {
111  hasWeights = false;
113  << "Illegal minimum weight " << minWeights
114  << " ... ignoring"
115  << endl;
116  }
117  else if (hasWeights && (cWeights.size() != numCells))
118  {
120  << "Number of weights (" << cWeights.size()
121  << ") != number of cells (" << numCells << ")"
122  << exit(FatalError);
123  }
124 
125  if (hasWeights)
126  {
127  // Convert to integers.
128  cellWeights.resize_nocopy(cWeights.size());
129  forAll(cellWeights, i)
130  {
131  cellWeights[i] = idx_t(cWeights[i]/minWeights);
132  }
133  }
134 
135 
136  // Check for user supplied weights and decomp options
137  if (coeffsDictPtr)
138  {
139  const dictionary& coeffDict = *coeffsDictPtr;
140 
141  word weightsFile;
142 
143  if (coeffDict.readIfPresent("method", method))
144  {
145  if (method != "recursive" && method != "k-way")
146  {
148  << "Method " << method << " in metisCoeffs in dictionary : "
149  << decompDict_.name()
150  << " should be 'recursive' or 'k-way'"
151  << exit(FatalError);
152  }
153 
154  Info<< "metisDecomp : Using Metis method " << method
155  << nl << endl;
156  }
157 
158  if (coeffDict.readIfPresent("options", options))
159  {
160  if (options.size() != METIS_NOPTIONS)
161  {
163  << "Number of options in metisCoeffs in dictionary : "
164  << decompDict_.name()
165  << " should be " << METIS_NOPTIONS
166  << exit(FatalError);
167  }
168 
169  Info<< "metisDecomp : Using Metis options " << options
170  << nl << endl;
171  }
172 
173  if (coeffDict.readIfPresent("processorWeights", procWeights))
174  {
175  if (procWeights.size() != nDomains_)
176  {
177  FatalIOErrorInFunction(coeffDict)
178  << "processorWeights (" << procWeights.size()
179  << ") != number of domains (" << nDomains_ << ")" << nl
180  << exit(FatalIOError);
181  }
182 
183  procWeights /= sum(procWeights);
184  }
185  }
186 
187  idx_t ncon = 1;
188  idx_t nProcs = nDomains_;
189 
190  // Addressing
191  ConstPrecisionAdaptor<idx_t, label, List> xadj_param(xadj);
192  ConstPrecisionAdaptor<idx_t, label, List> adjncy_param(adjncy);
193 
194  // Output: cell -> processor addressing
195  decomp.resize(numCells);
196  PrecisionAdaptor<idx_t, label, List> decomp_param(decomp, false);
197 
198  // Avoid potential nullptr issues with zero-sized arrays
199  labelList adjncy_dummy, xadj_dummy, decomp_dummy;
200  if (!numCells)
201  {
202  adjncy_dummy.resize(1, 0);
203  adjncy_param.set(adjncy_dummy);
204 
205  xadj_dummy.resize(2, 0);
206  xadj_param.set(xadj_dummy);
207 
208  decomp_dummy.resize(1, 0);
209  decomp_param.clear(); // Avoid propagating spurious values
210  decomp_param.set(decomp_dummy);
211  }
212 
213 
214  //
215  // Decompose
216  //
217 
218  // Output: number of cut edges
219  idx_t edgeCut = 0;
220 
221  if (method == "recursive")
222  {
223  METIS_PartGraphRecursive
224  (
225  &numCells, // num vertices in graph
226  &ncon, // num balancing constraints
227  xadj_param.constCast().data(), // indexing into adjncy
228  adjncy_param.constCast().data(), // neighbour info
229  cellWeights.data(), // vertex wts
230  nullptr, // vsize: total communication vol
231  faceWeights.data(), // edge wts
232  &nProcs, // nParts
233  procWeights.data(), // tpwgts
234  nullptr, // ubvec: processor imbalance (default)
235  options.data(),
236  &edgeCut,
237  decomp_param.ref().data()
238  );
239  }
240  else
241  {
242  METIS_PartGraphKway
243  (
244  &numCells, // num vertices in graph
245  &ncon, // num balancing constraints
246  xadj_param.constCast().data(), // indexing into adjncy
247  adjncy_param.constCast().data(), // neighbour info
248  cellWeights.data(), // vertex wts
249  nullptr, // vsize: total communication vol
250  faceWeights.data(), // edge wts
251  &nProcs, // nParts
252  procWeights.data(), // tpwgts
253  nullptr, // ubvec: processor imbalance (default)
254  options.data(),
255  &edgeCut,
256  decomp_param.ref().data()
257  );
258  }
260  return edgeCut;
261 }
262 
263 
264 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
265 
266 Foam::metisDecomp::metisDecomp(const label numDomains)
267 :
268  metisLikeDecomp(numDomains)
269 {}
270 
271 
273 (
274  const dictionary& decompDict,
275  const word& regionName
276 )
277 :
278  metisLikeDecomp(typeName, decompDict, regionName, selectionType::NULL_DICT)
279 {}
280 
281 
282 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:160
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
label nDomains_
Number of domains for the decomposition.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
bool empty() const noexcept
True if List is empty (ie, size() is zero)
Definition: UList.H:666
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:272
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:175
metisDecomp(const metisDecomp &)=delete
No copy construct.
const fileName & name() const noexcept
The dictionary name.
Definition: dictionaryI.H:41
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
Generic templated field type.
Definition: Field.H:62
A class for handling words, derived from Foam::string.
Definition: word.H:63
Domain decomposition using METIS-like data structures.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
defineTypeNameAndDebug(combustionModel, 0)
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:627
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
const dictionary & decompDict_
Top-level decomposition dictionary (eg, decomposeParDict)
messageStream Info
Information stream (stdout output on master, null elsewhere)
List< label > labelList
A List of labels.
Definition: List.H:62
virtual label decomposeSerial(const labelList &adjncy, const labelList &xadj, const List< scalar > &cellWeights, labelList &decomp) const
Decompose non-parallel.
Definition: metisDecomp.C:67
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and a sub-dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...