multiLevel.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 "multiLevel.H"
27 #include "IFstream.H"
28 #include "globalIndex.H"
29 #include "distributionMap.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace decompositionMethods
37 {
39 
41  (
43  multiLevel,
44  decomposer
45  );
46 
48  (
50  multiLevel,
51  distributor
52  );
53 }
54 }
55 
56 
57 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
58 
59 void Foam::decompositionMethods::multiLevel::subsetGlobalCellCells
60 (
61  const label nDomains,
62  const label domainI,
63  const labelList& dist,
64 
65  const labelListList& cellCells,
66  const labelList& set,
67  labelListList& subCellCells,
68  labelList& cutConnections
69 ) const
70 {
71  // Determine new index for cells by inverting subset
72  labelList oldToNew(invert(cellCells.size(), set));
73 
74  globalIndex globalCells(cellCells.size());
75 
76  // Subset locally the elements for which I have data
77  subCellCells = UIndirectList<labelList>(cellCells, set);
78 
79  // Get new indices for neighbouring processors
80  List<Map<label>> compactMap;
81  distributionMap map(globalCells, subCellCells, compactMap);
82  map.distribute(oldToNew);
83  labelList allDist(dist);
84  map.distribute(allDist);
85 
86  // Now we have:
87  // oldToNew : the locally-compact numbering of all our cellCells. -1 if
88  // cellCell is not in set.
89  // allDist : destination domain for all our cellCells
90  // subCellCells : indexes into oldToNew and allDist
91 
92  // Globally compact numbering for cells in set.
93  globalIndex globalSubCells(set.size());
94 
95  // Now subCellCells contains indices into oldToNew which are the
96  // new locations of the neighbouring cells.
97 
98  cutConnections.setSize(nDomains);
99  cutConnections = 0;
100 
101  forAll(subCellCells, subCelli)
102  {
103  labelList& cCells = subCellCells[subCelli];
104 
105  // Keep the connections to valid mapped cells
106  label newI = 0;
107  forAll(cCells, i)
108  {
109  // Get locally-compact cell index of neighbouring cell
110  label nbrCelli = oldToNew[cCells[i]];
111  if (nbrCelli == -1)
112  {
113  cutConnections[allDist[cCells[i]]]++;
114  }
115  else
116  {
117  // Reconvert local cell index into global one
118 
119  // Get original neighbour
120  label celli = set[subCelli];
121  label oldNbrCelli = cellCells[celli][i];
122  // Get processor from original neighbour
123  label proci = globalCells.whichProcID(oldNbrCelli);
124  // Convert into global compact numbering
125  cCells[newI++] = globalSubCells.toGlobal(proci, nbrCelli);
126  }
127  }
128  cCells.setSize(newI);
129  }
130 }
131 
132 
134 (
135  const labelListList& pointPoints,
136  const pointField& points,
137  const scalarField& pointWeights,
138  const labelList& pointMap, // map back to original points
139  const label levelI,
140  labelField& decomp
141 )
142 {
143  labelList dist
144  (
145  methods_[levelI].decompose
146  (
147  pointPoints,
148  points,
149  pointWeights
150  )
151  );
152 
153  forAll(pointMap, i)
154  {
155  label orig = pointMap[i];
156  decomp[orig] += dist[i];
157  }
158 
159  if (levelI != methods_.size()-1)
160  {
161  // Recurse
162 
163  // Determine points per domain
164  label n = methods_[levelI].nDomains();
165  labelListList domainToPoints(invertOneToMany(n, dist));
166 
167  // 'Make space' for new levels of decomposition
168  decomp *= methods_[levelI+1].nDomains();
169 
170  // Extract processor+local index from point-point addressing
171  if (debug && Pstream::master())
172  {
173  Pout<< "Decomposition at level " << levelI << " :" << endl;
174  }
175 
176  for (label domainI = 0; domainI < n; domainI++)
177  {
178  // Extract elements for current domain
179  const labelList domainPoints(findIndices(dist, domainI));
180 
181  // Subset point-wise data.
182  pointField subPoints(points, domainPoints);
183  scalarField subWeights(pointWeights, domainPoints);
184  labelList subPointMap(UIndirectList<label>(pointMap, domainPoints));
185  // Subset point-point addressing (adapt global numbering)
186  labelListList subPointPoints;
187  labelList nOutsideConnections;
188  subsetGlobalCellCells
189  (
190  n,
191  domainI,
192  dist,
193 
194  pointPoints,
195  domainPoints,
196 
197  subPointPoints,
198  nOutsideConnections
199  );
200 
201  label nPoints = returnReduce(domainPoints.size(), plusOp<label>());
202  Pstream::listCombineGather(nOutsideConnections, plusEqOp<label>());
203  Pstream::listCombineScatter(nOutsideConnections);
204  label nPatches = 0;
205  label nFaces = 0;
206  forAll(nOutsideConnections, i)
207  {
208  if (nOutsideConnections[i] > 0)
209  {
210  nPatches++;
211  nFaces += nOutsideConnections[i];
212  }
213  }
214 
215  string oldPrefix;
216  if (debug && Pstream::master())
217  {
218  Pout<< " Domain " << domainI << nl
219  << " Number of cells = " << nPoints << nl
220  << " Number of inter-domain patches = " << nPatches
221  << nl
222  << " Number of inter-domain faces = " << nFaces << nl
223  << endl;
224  oldPrefix = Pout.prefix();
225  Pout.prefix() = " " + oldPrefix;
226  }
227 
228  decompose
229  (
230  subPointPoints,
231  subPoints,
232  subWeights,
233  subPointMap,
234  levelI+1,
235 
236  decomp
237  );
238  if (debug && Pstream::master())
239  {
240  Pout.prefix() = oldPrefix;
241  }
242  }
243 
244 
245  if (debug)
246  {
247  // Do straight decompose of two levels
248  label nNext = methods_[levelI+1].nDomains();
249  label nTotal = n*nNext;
250 
251  // Retrieve original level0 dictionary and modify number of domains
253  decompositionDict_.optionalSubDict(typeName + "Coeffs").begin();
254  dictionary myDict = iter().dict();
255  myDict.set("numberOfSubdomains", nTotal);
256 
257  if (debug && Pstream::master())
258  {
259  Pout<< "Reference decomposition with " << myDict << " :"
260  << endl;
261  }
262 
263  autoPtr<decompositionMethod> method0 =
265  (
266  myDict
267  );
268 
269  labelList dist
270  (
271  method0().decompose
272  (
273  pointPoints,
274  points,
275  pointWeights
276  )
277  );
278 
279  for (label blockI = 0; blockI < n; blockI++)
280  {
281  // Count the number in between blocks of nNext size
282 
283  label nPoints = 0;
284  labelList nOutsideConnections(n, 0);
285  forAll(pointPoints, pointi)
286  {
287  if ((dist[pointi] / nNext) == blockI)
288  {
289  nPoints++;
290 
291  const labelList& pPoints = pointPoints[pointi];
292 
293  forAll(pPoints, i)
294  {
295  label distBlockI = dist[pPoints[i]] / nNext;
296  if (distBlockI != blockI)
297  {
298  nOutsideConnections[distBlockI]++;
299  }
300  }
301  }
302  }
303 
304  reduce(nPoints, plusOp<label>());
306  (
307  nOutsideConnections,
308  plusEqOp<label>()
309  );
310  Pstream::listCombineScatter(nOutsideConnections);
311  label nPatches = 0;
312  label nFaces = 0;
313  forAll(nOutsideConnections, i)
314  {
315  if (nOutsideConnections[i] > 0)
316  {
317  nPatches++;
318  nFaces += nOutsideConnections[i];
319  }
320  }
321 
322  if (debug && Pstream::master())
323  {
324  Pout<< " Domain " << blockI << nl
325  << " Number of cells = " << nPoints << nl
326  << " Number of inter-domain patches = "
327  << nPatches << nl
328  << " Number of inter-domain faces = " << nFaces
329  << nl << endl;
330  }
331  }
332  }
333  }
334 }
335 
336 
337 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
338 
340 (
341  const dictionary& decompositionDict
342 )
343 :
344  decompositionMethod(decompositionDict),
345  methodsDict_(decompositionDict_.optionalSubDict(typeName + "Coeffs"))
346 {
347  methods_.setSize(methodsDict_.size());
348  label i = 0;
349  forAllConstIter(dictionary, methodsDict_, iter)
350  {
351  methods_.set(i++, decompositionMethod::NewDecomposer(iter().dict()));
352  }
353 
354  label n = 1;
355  Info<< "decompositionMethod " << type() << " :" << endl;
356  forAll(methods_, i)
357  {
358  Info<< " level " << i << " decomposing with " << methods_[i].type()
359  << " into " << methods_[i].nDomains() << " subdomains." << endl;
360 
361  n *= methods_[i].nDomains();
362  }
363 
364  if (n != nDomains())
365  {
367  << "Top level decomposition specifies " << nDomains()
368  << " domains which is not equal to the product of"
369  << " all sub domains " << n
370  << exit(FatalError);
371  }
372 }
373 
374 
375 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
376 
378 (
379  const polyMesh& mesh,
380  const pointField& cellCentres,
381  const scalarField& cellWeights
382 )
383 {
384  CompactListList<label> cellCells;
385  calcCellCells
386  (
387  mesh,
388  identityMap(cellCentres.size()),
389  cellCentres.size(),
390  true,
391  cellCells
392  );
393 
394  labelField decomp(cellCentres.size(), 0);
395  labelList cellMap(identityMap(cellCentres.size()));
396 
397  decompose
398  (
399  cellCells.list(),
400  cellCentres,
401  cellWeights,
402  cellMap, // map back to original cells
403  0,
404 
405  decomp
406  );
407 
408  return decomp;
409 }
410 
411 
413 (
414  const labelListList& globalPointPoints,
415  const pointField& points,
416  const scalarField& pointWeights
417 )
418 {
419  labelField decomp(points.size(), 0);
420  labelList pointMap(identityMap(points.size()));
421 
422  decompose
423  (
424  globalPointPoints,
425  points,
426  pointWeights,
427  pointMap, // map back to original points
428  0,
429 
430  decomp
431  );
432 
433  return decomp;
434 }
435 
436 
437 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
List< Container > list() const
Convert to List<Container>
friend class const_iterator
Definition: UILList.H:78
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Abstract base class for decomposition.
static autoPtr< decompositionMethod > NewDecomposer(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
Decomposition given using consecutive application of decomposers.
Definition: multiLevel.H:51
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Inherit decompose from decompositionMethod.
multiLevel(const dictionary &decompositionDict)
Construct given the decomposition dictionary.
Definition: multiLevel.C:340
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const string & prefix() const
Return the prefix of the stream.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const pointField & points
label nPoints
addToRunTimeSelectionTable(decompositionMethod, metis, decomposer)
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
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurrences of given element. Linear search.
messageStream Info
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Field< label > labelField
Specialisation of Field<T> for label.
Definition: labelField.H:49
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label nPatches
Definition: readKivaGrid.H:396
dictionary dict