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-2019 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"
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 
55 namespace Foam
56 {
57  defineTypeNameAndDebug(metisDecomp, 0);
58 
60  (
61  decompositionMethod,
62  metisDecomp,
63  dictionary
64  );
65 
67  (
68  decompositionMethod,
69  metisDecomp,
70  dictionaryRegion
71  );
72 }
73 
74 
75 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
76 
78 (
79  const labelList& adjncy,
80  const labelList& xadj,
81  const List<scalar>& cWeights,
82  labelList& decomp
83 ) const
84 {
85  // Method of decomposition
86  // recursive: multi-level recursive bisection (default)
87  // k-way: multi-level k-way
88  word method("recursive");
89 
90  const dictionary* coeffsDictPtr = decompDict_.findDict("metisCoeffs");
91 
92  idx_t numCells = xadj.size()-1;
93 
94  // Decomposition options
95  List<idx_t> options(METIS_NOPTIONS);
96  METIS_SetDefaultOptions(options.data());
97 
98  // Processor weights initialised with no size, only used if specified in
99  // a file
100  Field<real_t> processorWeights;
101 
102  // Cell weights (so on the vertices of the dual)
103  List<idx_t> cellWeights;
104 
105  // Face weights (so on the edges of the dual)
106  List<idx_t> faceWeights;
107 
108  // Check for externally provided cellweights and if so initialise weights
109  // Note: min, not gMin since routine runs on master only.
110  const scalar minWeights = min(cWeights);
111 
112  if (!cWeights.empty())
113  {
114  if (minWeights <= 0)
115  {
117  << "Illegal minimum weight " << minWeights
118  << endl;
119  }
120 
121  if (cWeights.size() != numCells)
122  {
124  << "Number of cell weights " << cWeights.size()
125  << " does not equal number of cells " << numCells
126  << exit(FatalError);
127  }
128 
129  // Convert to integers.
130  cellWeights.setSize(cWeights.size());
131  forAll(cellWeights, i)
132  {
133  cellWeights[i] = idx_t(cWeights[i]/minWeights);
134  }
135  }
136 
137 
138  // Check for user supplied weights and decomp options
139  if (coeffsDictPtr)
140  {
141  const dictionary& coeffDict = *coeffsDictPtr;
142 
143  word weightsFile;
144 
145  if (coeffDict.readIfPresent("method", method))
146  {
147  if (method != "recursive" && method != "k-way")
148  {
150  << "Method " << method << " in metisCoeffs in dictionary : "
151  << decompDict_.name()
152  << " should be 'recursive' or 'k-way'"
153  << exit(FatalError);
154  }
155 
156  Info<< "metisDecomp : Using Metis method " << method
157  << nl << endl;
158  }
159 
160  if (coeffDict.readIfPresent("options", options))
161  {
162  if (options.size() != METIS_NOPTIONS)
163  {
165  << "Number of options in metisCoeffs in dictionary : "
166  << decompDict_.name()
167  << " should be " << METIS_NOPTIONS
168  << exit(FatalError);
169  }
170 
171  Info<< "metisDecomp : Using Metis options " << options
172  << nl << endl;
173  }
174 
175  if (coeffDict.readIfPresent("processorWeights", processorWeights))
176  {
177  processorWeights /= sum(processorWeights);
178 
179  if (processorWeights.size() != nDomains_)
180  {
182  << "Number of processor weights "
183  << processorWeights.size()
184  << " does not equal number of domains " << nDomains_
185  << exit(FatalError);
186  }
187  }
188  }
189 
190  idx_t ncon = 1;
191  idx_t nProcs = nDomains_;
192 
193  // Addressing
194  ConstPrecisionAdaptor<idx_t, label, List> xadj_param(xadj);
195  ConstPrecisionAdaptor<idx_t, label, List> adjncy_param(adjncy);
196 
197  // Output: cell -> processor addressing
198  decomp.resize(numCells);
199  PrecisionAdaptor<idx_t, label, List> decomp_param(decomp);
200 
201  // Output: number of cut edges
202  idx_t edgeCut = 0;
203 
204  if (method == "recursive")
205  {
206  METIS_PartGraphRecursive
207  (
208  &numCells, // num vertices in graph
209  &ncon, // num balancing constraints
210  xadj_param.constCast().data(), // indexing into adjncy
211  adjncy_param.constCast().data(), // neighbour info
212  cellWeights.data(), // vertex wts
213  nullptr, // vsize: total communication vol
214  faceWeights.data(), // edge wts
215  &nProcs, // nParts
216  processorWeights.data(), // tpwgts
217  nullptr, // ubvec: processor imbalance (default)
218  options.data(),
219  &edgeCut,
220  decomp_param.ref().data()
221  );
222  }
223  else
224  {
225  METIS_PartGraphKway
226  (
227  &numCells, // num vertices in graph
228  &ncon, // num balancing constraints
229  xadj_param.constCast().data(), // indexing into adjncy
230  adjncy_param.constCast().data(), // neighbour info
231  cellWeights.data(), // vertex wts
232  nullptr, // vsize: total communication vol
233  faceWeights.data(), // edge wts
234  &nProcs, // nParts
235  processorWeights.data(), // tpwgts
236  nullptr, // ubvec: processor imbalance (default)
237  options.data(),
238  &edgeCut,
239  decomp_param.ref().data()
240  );
241  }
242 
243  return edgeCut;
244 }
245 
246 
247 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
248 
249 Foam::metisDecomp::metisDecomp(const dictionary& decompDict)
250 :
251  metisLikeDecomp(typeName, decompDict, selectionType::NULL_DICT)
252 {}
253 
254 
256 (
257  const dictionary& decompDict,
258  const word& regionName
259 )
260 :
261  metisLikeDecomp(typeName, decompDict, regionName, selectionType::NULL_DICT)
262 {}
263 
264 
265 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::decompositionMethod::nDomains_
label nDomains_
Number of domains for the decomposition.
Definition: decompositionMethod.H:95
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionary.C:503
metisDecomp.H
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::decompositionMethod::decompDict_
const dictionary & decompDict_
Top-level decomposition dictionary (eg, decomposeParDict)
Definition: decompositionMethod.H:89
PrecisionAdaptor.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::metisDecomp::decomposeSerial
virtual label decomposeSerial(const labelList &adjncy, const labelList &xadj, const List< scalar > &cellWeights, labelList &decomp) const
Decompose non-parallel.
Definition: dummyMetisDecomp.C:59
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::dictionary::name
const fileName & name() const
The dictionary name.
Definition: dictionary.H:446
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:62
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::metisDecomp::metisDecomp
metisDecomp(const metisDecomp &)=delete
No copy construct.
Foam::FatalError
error FatalError
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Time.H
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294