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-2021 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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"
32#include "PrecisionAdaptor.H"
33
34// Probably not needed...
35#define MPICH_SKIP_MPICXX
36#define OMPI_SKIP_MPICXX
37
38#include "metis.h"
39
40// Provide a clear error message if we have a severe size mismatch
41// Allow widening, but not narrowing
42//
43// Metis has an 'idx_t' type, but the IDXTYPEWIDTH define is perhaps
44// more future-proof?
45//#ifdef IDXTYPEWIDTH
46//static_assert
47//(
48// sizeof(Foam::label) > (IDXTYPEWIDTH/8),
49// "sizeof(Foam::label) > (IDXTYPEWIDTH/8), check your metis headers"
50//);
51//#endif
52
53// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
54
55namespace Foam
56{
57 defineTypeNameAndDebug(metisDecomp, 0);
59 (
60 decompositionMethod,
61 metisDecomp,
62 dictionary
63 );
64}
65
66
67// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
68
70(
71 const labelList& adjncy,
72 const labelList& xadj,
73 const List<scalar>& cWeights,
74 labelList& decomp
75) const
76{
77 // Method of decomposition
78 // recursive: multi-level recursive bisection (default)
79 // k-way: multi-level k-way
80 word method("recursive");
81
82 const dictionary* coeffsDictPtr = decompDict_.findDict("metisCoeffs");
83
84 idx_t numCells = max(0, (xadj.size()-1));
85
86 // Decomposition options
87 List<idx_t> options(METIS_NOPTIONS);
88 METIS_SetDefaultOptions(options.data());
89
90 // Processor weights initialised with no size, only used if specified in
91 // a file
92 Field<real_t> procWeights;
93
94 // Cell weights (so on the vertices of the dual)
95 List<idx_t> cellWeights;
96
97 // Face weights (so on the edges of the dual)
98 List<idx_t> faceWeights;
99
100 // Check for externally provided cellweights and if so initialise weights
101 // Note: min, not gMin since routine runs on master only.
102 const scalar minWeights = min(cWeights);
103
104 if (!cWeights.empty())
105 {
106 if (minWeights <= 0)
107 {
109 << "Illegal minimum weight " << minWeights
110 << endl;
111 }
112
113 if (cWeights.size() != numCells)
114 {
116 << "Number of cell weights " << cWeights.size()
117 << " does not equal number of cells " << numCells
118 << exit(FatalError);
119 }
120
121 // Convert to integers.
122 cellWeights.setSize(cWeights.size());
123 forAll(cellWeights, i)
124 {
125 cellWeights[i] = idx_t(cWeights[i]/minWeights);
126 }
127 }
128
129
130 // Check for user supplied weights and decomp options
131 if (coeffsDictPtr)
132 {
133 const dictionary& coeffDict = *coeffsDictPtr;
134
135 word weightsFile;
136
137 if (coeffDict.readIfPresent("method", method))
138 {
139 if (method != "recursive" && method != "k-way")
140 {
142 << "Method " << method << " in metisCoeffs in dictionary : "
143 << decompDict_.name()
144 << " should be 'recursive' or 'k-way'"
145 << exit(FatalError);
146 }
147
148 Info<< "metisDecomp : Using Metis method " << method
149 << nl << endl;
150 }
151
152 if (coeffDict.readIfPresent("options", options))
153 {
154 if (options.size() != METIS_NOPTIONS)
155 {
157 << "Number of options in metisCoeffs in dictionary : "
158 << decompDict_.name()
159 << " should be " << METIS_NOPTIONS
160 << exit(FatalError);
161 }
162
163 Info<< "metisDecomp : Using Metis options " << options
164 << nl << endl;
165 }
166
167 if (coeffDict.readIfPresent("processorWeights", procWeights))
168 {
169 if (procWeights.size() != nDomains_)
170 {
171 FatalIOErrorInFunction(coeffDict)
172 << "processorWeights (" << procWeights.size()
173 << ") != number of domains (" << nDomains_ << ")" << nl
174 << exit(FatalIOError);
175 }
176
177 procWeights /= sum(procWeights);
178 }
179 }
180
181 idx_t ncon = 1;
182 idx_t nProcs = nDomains_;
183
184 // Addressing
185 ConstPrecisionAdaptor<idx_t, label, List> xadj_param(xadj);
186 ConstPrecisionAdaptor<idx_t, label, List> adjncy_param(adjncy);
187
188 // Output: cell -> processor addressing
189 decomp.resize(numCells);
190 PrecisionAdaptor<idx_t, label, List> decomp_param(decomp, false);
191
192 // Avoid potential nullptr issues with zero-sized arrays
193 labelList adjncy_dummy, xadj_dummy, decomp_dummy;
194 if (!numCells)
195 {
196 adjncy_dummy.resize(1, 0);
197 adjncy_param.set(adjncy_dummy);
198
199 xadj_dummy.resize(2, 0);
200 xadj_param.set(xadj_dummy);
201
202 decomp_dummy.resize(1, 0);
203 decomp_param.clear(); // Avoid propagating spurious values
204 decomp_param.set(decomp_dummy);
205 }
206
207
208 //
209 // Decompose
210 //
211
212 // Output: number of cut edges
213 idx_t edgeCut = 0;
214
215 if (method == "recursive")
216 {
217 METIS_PartGraphRecursive
218 (
219 &numCells, // num vertices in graph
220 &ncon, // num balancing constraints
221 xadj_param.constCast().data(), // indexing into adjncy
222 adjncy_param.constCast().data(), // neighbour info
223 cellWeights.data(), // vertex wts
224 nullptr, // vsize: total communication vol
225 faceWeights.data(), // edge wts
226 &nProcs, // nParts
227 procWeights.data(), // tpwgts
228 nullptr, // ubvec: processor imbalance (default)
229 options.data(),
230 &edgeCut,
231 decomp_param.ref().data()
232 );
233 }
234 else
235 {
236 METIS_PartGraphKway
237 (
238 &numCells, // num vertices in graph
239 &ncon, // num balancing constraints
240 xadj_param.constCast().data(), // indexing into adjncy
241 adjncy_param.constCast().data(), // neighbour info
242 cellWeights.data(), // vertex wts
243 nullptr, // vsize: total communication vol
244 faceWeights.data(), // edge wts
245 &nProcs, // nParts
246 procWeights.data(), // tpwgts
247 nullptr, // ubvec: processor imbalance (default)
248 options.data(),
249 &edgeCut,
250 decomp_param.ref().data()
251 );
252 }
253
254 return edgeCut;
255}
256
257
258// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
259
261(
262 const dictionary& decompDict,
263 const word& regionName
264)
265:
266 metisLikeDecomp(typeName, decompDict, regionName, selectionType::NULL_DICT)
267{}
268
269
270// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
label nDomains_
Number of domains for the decomposition.
const dictionary & decompDict_
Top-level decomposition dictionary (eg, decomposeParDict)
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
const fileName & name() const noexcept
The dictionary name.
Definition: dictionaryI.H:48
Metis domain decomposition.
Definition: metisDecomp.H:95
virtual label decomposeSerial(const labelList &adjncy, const labelList &xadj, const List< scalar > &cellWeights, labelList &decomp) const
Decompose non-parallel.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
Foam::word regionName(Foam::polyMesh::defaultRegion)
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
#define WarningInFunction
Report a warning using Foam::Warning.
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
List< label > labelList
A List of labels.
Definition: List.H:66
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
IOerror FatalIOError
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333