multiLevelDecomp.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-2017 OpenFOAM Foundation
9  Copyright (C) 2017-2024 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 "multiLevelDecomp.H"
31 #include "IFstream.H"
32 #include "globalIndex.H"
33 #include "globalMeshData.H"
34 #include "mapDistribute.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(multiLevelDecomp, 0);
42  (
43  decompositionMethod,
44  multiLevelDecomp,
45  dictionary
46  );
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::multiLevelDecomp::createMethodsDict()
53 {
54  methodsDict_.clear();
55 
56  word defaultMethod;
57  labelList domains;
58 
59  label nTotal = 0;
60  label nLevels = 0;
61 
62  // Found (non-recursive, no patterns) "method" and "domains" ?
63  // Allow as quick short-cut entry
64  if
65  (
66  // non-recursive, no patterns
67  coeffsDict_.readIfPresent("method", defaultMethod, keyType::LITERAL)
68  // non-recursive, no patterns
69  && coeffsDict_.readIfPresent("domains", domains, keyType::LITERAL)
70  )
71  {
72  // Short-cut version specified by method, domains only
73 
74  nTotal = (domains.empty() ? 0 : 1);
75 
76  for (const label n : domains)
77  {
78  nTotal *= n;
79  ++nLevels;
80  }
81 
82  if (nTotal == 1)
83  {
84  // Emit Warning
85  nTotal = nDomains();
86  nLevels = 1;
87 
88  domains.setSize(1);
89  domains[0] = nTotal;
90  }
91  else if (nTotal > 0 && nTotal < nDomains() && !(nDomains() % nTotal))
92  {
93  // nTotal < nDomains, but with an integral factor,
94  // which we insert as level 0
95  ++nLevels;
96 
97  labelList old(std::move(domains));
98 
99  domains.setSize(old.size()+1);
100 
101  domains[0] = nDomains() / nTotal;
102  forAll(old, i)
103  {
104  domains[i+1] = old[i];
105  }
106  nTotal *= domains[0];
107 
108  Info<<" inferred level0 with " << domains[0]
109  << " domains" << nl << nl;
110  }
111 
112  if (!nLevels || nTotal != nDomains())
113  {
115  << "Top level decomposition specifies " << nDomains()
116  << " domains which is not equal to the product of"
117  << " all sub domains " << nTotal
118  << exit(FatalError);
119  }
120 
121  // Create editable methods dictionaries
122  nLevels = 0;
123 
124  // Common coeffs dictionary
125  const dictionary& subMethodCoeffsDict
126  (
128  (
129  coeffsDict_,
130  defaultMethod + "Coeffs",
131  selectionType::NULL_DICT
132  )
133  );
134 
135  for (const label n : domains)
136  {
137  const word levelName("level" + Foam::name(nLevels++));
138 
139  entry* dictptr = methodsDict_.set(levelName, dictionary());
140 
141  dictionary& dict = dictptr->dict();
142  dict.add("method", defaultMethod);
143  dict.add("numberOfSubdomains", n);
144 
145  // Inject coeffs dictionary too
146  if (subMethodCoeffsDict.size())
147  {
148  dict.add(subMethodCoeffsDict.dictName(), subMethodCoeffsDict);
149  }
150  }
151  }
152  else
153  {
154  // Specified by full dictionaries
155 
156  // Create editable methods dictionaries
157  // - Only consider sub-dictionaries with a "numberOfSubdomains" entry
158  // This automatically filters out any coeffs dictionaries
159 
160  for (const entry& dEntry : coeffsDict_)
161  {
162  word methodName;
163 
164  if
165  (
166  dEntry.isDict()
167  // non-recursive, no patterns
168  && dEntry.dict().found("numberOfSubdomains", keyType::LITERAL)
169  )
170  {
171  // No method specified? can use a default method?
172 
173  const bool addDefaultMethod
174  (
175  !(dEntry.dict().found("method", keyType::LITERAL))
176  && !defaultMethod.empty()
177  );
178 
179  entry* e = methodsDict_.add(dEntry);
180 
181  if (addDefaultMethod && e && e->isDict())
182  {
183  e->dict().add("method", defaultMethod);
184  }
185  }
186  }
187  }
188 }
189 
190 
191 void Foam::multiLevelDecomp::setMethods()
192 {
193  // Assuming methodsDict_ has be properly created, convert the method
194  // dictionaries to actual methods
195 
196  label nLevels = 0;
197 
198  methods_.resize_null(methodsDict_.size());
199  for (const entry& dEntry : methodsDict_)
200  {
201  // Dictionary entries only
202  // - these method dictionaries are non-regional
203  if (dEntry.isDict())
204  {
205  methods_.set
206  (
207  nLevels++,
208  // non-verbose would be nicer
209  decompositionMethod::New(dEntry.dict())
210  );
211  }
212  }
213 
214  methods_.resize(nLevels);
215 
216  // Verify that nTotal is correct based on what each method delivers
217 
218  Info<< nl
219  << "Decompose " << type() << " [" << nDomains() << "] in "
220  << nLevels << " levels:" << endl;
221 
222  label nTotal = 1;
223  forAll(methods_, i)
224  {
225  Info<< " level " << i << " : " << methods_[i].type()
226  << " [" << methods_[i].nDomains() << "]" << endl;
227 
228  nTotal *= methods_[i].nDomains();
229  }
230 
231  if (nTotal != nDomains())
232  {
234  << "Top level decomposition specifies " << nDomains()
235  << " domains which is not equal to the product of"
236  << " all sub domains " << nTotal
237  << exit(FatalError);
238  }
239 }
240 
241 
242 // Given a subset of cells determine the new global indices. The problem
243 // is in the cells from neighbouring processors which need to be renumbered.
244 void Foam::multiLevelDecomp::subsetGlobalCellCells
245 (
246  const label nDomains,
247  const label domainI,
248  const labelList& dist,
249 
250  const labelListList& cellCells,
251  const labelList& set,
252  labelListList& subCellCells,
253  labelList& cutConnections
254 ) const
255 {
256  const globalIndex globalCells(cellCells.size());
257 
258  // Determine new index for cells by inverting subset
259  labelList oldToNew(invert(cellCells.size(), set));
260 
261  // Subset locally the elements for which I have data
262  subCellCells = UIndirectList<labelList>(cellCells, set);
263 
264  // Get new indices for neighbouring processors
265  List<Map<label>> compactMap;
266  mapDistribute map(globalCells, subCellCells, compactMap);
267  map.distribute(oldToNew);
268  labelList allDist(dist);
269  map.distribute(allDist);
270 
271  // Now we have:
272  // oldToNew : the locally-compact numbering of all our cellCells. -1 if
273  // cellCell is not in set.
274  // allDist : destination domain for all our cellCells
275  // subCellCells : indexes into oldToNew and allDist
276 
277  // Globally compact numbering for cells in set.
278  const globalIndex globalSubCells(set.size());
279 
280  // Now subCellCells contains indices into oldToNew which are the
281  // new locations of the neighbouring cells.
282 
283  cutConnections.resize_nocopy(nDomains);
284  cutConnections = 0;
285 
286  forAll(subCellCells, subCelli)
287  {
288  labelList& cCells = subCellCells[subCelli];
289 
290  // Keep the connections to valid mapped cells
291  label newI = 0;
292  forAll(cCells, i)
293  {
294  // Get locally-compact cell index of neighbouring cell
295  const label nbrCelli = oldToNew[cCells[i]];
296  if (nbrCelli == -1)
297  {
298  cutConnections[allDist[cCells[i]]]++;
299  }
300  else
301  {
302  // Reconvert local cell index into global one
303 
304  // Get original neighbour
305  const label celli = set[subCelli];
306  const label oldNbrCelli = cellCells[celli][i];
307  // Get processor from original neighbour
308  const label proci = globalCells.whichProcID(oldNbrCelli);
309  // Convert into global compact numbering
310  cCells[newI++] = globalSubCells.toGlobal(proci, nbrCelli);
311  }
312  }
313  cCells.setSize(newI);
314  }
315 }
316 
317 
318 void Foam::multiLevelDecomp::decompose
319 (
320  const labelListList& pointPoints,
321  const pointField& points,
322  const scalarField& pointWeights,
323  const labelUList& pointMap, // map back to original points
324  const label currLevel,
325  const label leafOffset,
326 
327  labelList& finalDecomp
328 ) const
329 {
330  labelList dist
331  (
332  methods_[currLevel].decompose
333  (
334  pointPoints,
335  points,
336  pointWeights
337  )
338  );
339 
340  // The next recursion level
341  const label nextLevel = currLevel+1;
342 
343  // Number of domains at this current level
344  const label nCurrDomains = methods_[currLevel].nDomains();
345 
346  // Calculate the domain remapping.
347  // The decompose() method delivers a distribution of [0..nDomains-1]
348  // which we map to the final location according to the decomposition
349  // leaf we are on.
350 
351  labelList domainLookup(nCurrDomains);
352  {
353  label sizes = 1; // Cumulative number of domains
354  for (label i = 0; i <= currLevel; ++i)
355  {
356  sizes *= methods_[i].nDomains();
357  }
358 
359  // Distribution of domains at this level
360  sizes = this->nDomains() / sizes;
361 
362  forAll(domainLookup, i)
363  {
364  domainLookup[i] = i * sizes + leafOffset;
365  }
366  }
367 
368  if (debug)
369  {
370  Info<< "Distribute at level " << currLevel
371  << " to domains" << nl
372  << flatOutput(domainLookup) << endl;
373  }
374 
375  // Extract processor+local index from point-point addressing
376  forAll(pointMap, i)
377  {
378  const label orig = pointMap[i];
379  finalDecomp[orig] = domainLookup[dist[i]];
380  }
381 
382  if (nextLevel < methods_.size())
383  {
384  // Recurse
385 
386  // Extract processor+local index from point-point addressing
387  if (debug && Pstream::master())
388  {
389  Pout<< "Decomposition at level " << currLevel << " :" << endl;
390  }
391 
392  for (label domainI = 0; domainI < nCurrDomains; ++domainI)
393  {
394  // Extract elements for current domain
395  const labelList domainPoints(findIndices(dist, domainI));
396 
397  // Subset point-wise data.
398  pointField subPoints(points, domainPoints);
399  scalarField subWeights;
400  if (pointWeights.size() == points.size())
401  {
402  subWeights = scalarField(pointWeights, domainPoints);
403  }
404 
405  labelList subPointMap(labelUIndList(pointMap, domainPoints));
406  // Subset point-point addressing (adapt global numbering)
407  labelListList subPointPoints;
408  labelList nOutsideConnections;
409  subsetGlobalCellCells
410  (
411  nCurrDomains,
412  domainI,
413  dist,
414 
415  pointPoints,
416  domainPoints,
417 
418  subPointPoints,
419  nOutsideConnections
420  );
421 
422  label nPoints = returnReduce(domainPoints.size(), sumOp<label>());
423 
424  Pstream::listCombineReduce(nOutsideConnections, plusEqOp<label>());
425  label nPatches = 0;
426  label nFaces = 0;
427  for (const label nConnect : nOutsideConnections)
428  {
429  if (nConnect > 0)
430  {
431  ++nPatches;
432  nFaces += nConnect;
433  }
434  }
435 
436  string oldPrefix;
437  if (debug && Pstream::master())
438  {
439  Pout<< " Domain " << domainI << nl
440  << " Number of cells = " << nPoints << nl
441  << " Number of inter-domain patches = " << nPatches
442  << nl
443  << " Number of inter-domain faces = " << nFaces << nl
444  << endl;
445  oldPrefix = Pout.prefix();
446  Pout.prefix() = " " + oldPrefix;
447  }
448 
449  decompose
450  (
451  subPointPoints,
452  subPoints,
453  subWeights,
454  subPointMap,
455  nextLevel,
456  domainLookup[domainI], // The offset for this level and leaf
457 
458  finalDecomp
459  );
460  if (debug && Pstream::master())
461  {
462  Pout.prefix() = oldPrefix;
463  }
464  }
465 
466 
467  if (debug)
468  {
469  // Do straight decompose of two levels
470  const label nNext = methods_[nextLevel].nDomains();
471  const label nTotal = nCurrDomains * nNext;
472 
473  // Get original level0 dictionary and modify numberOfSubdomains
474  dictionary level0Dict;
475  for (const entry& dEntry : methodsDict_)
476  {
477  if (dEntry.isDict())
478  {
479  level0Dict = dEntry.dict();
480  break;
481  }
482  }
483  level0Dict.set("numberOfSubdomains", nTotal);
484 
485  if (debug && Pstream::master())
486  {
487  Pout<< "Reference decomposition with " << level0Dict << " :"
488  << endl;
489  }
490 
491  autoPtr<decompositionMethod> method0 = decompositionMethod::New
492  (
493  level0Dict
494  );
495  labelList dist
496  (
497  method0().decompose
498  (
499  pointPoints,
500  points,
501  pointWeights
502  )
503  );
504 
505  for (label blockI = 0; blockI < nCurrDomains; ++blockI)
506  {
507  // Count the number inbetween blocks of nNext size
508 
509  label nPoints = 0;
510  labelList nOutsideConnections(nCurrDomains, Foam::zero{});
511  forAll(pointPoints, pointi)
512  {
513  if ((dist[pointi] / nNext) == blockI)
514  {
515  nPoints++;
516 
517  const labelList& pPoints = pointPoints[pointi];
518 
519  forAll(pPoints, i)
520  {
521  const label distBlockI = dist[pPoints[i]] / nNext;
522  if (distBlockI != blockI)
523  {
524  nOutsideConnections[distBlockI]++;
525  }
526  }
527  }
528  }
529 
530  reduce(nPoints, sumOp<label>());
532  (
533  nOutsideConnections,
534  plusEqOp<label>()
535  );
536 
537  label nPatches = 0;
538  label nFaces = 0;
539  for (const label nConnect : nOutsideConnections)
540  {
541  if (nConnect > 0)
542  {
543  ++nPatches;
544  nFaces += nConnect;
545  }
546  }
547 
548  if (debug && Pstream::master())
549  {
550  Pout<< " Domain " << blockI << nl
551  << " Number of cells = " << nPoints << nl
552  << " Number of inter-domain patches = "
553  << nPatches << nl
554  << " Number of inter-domain faces = " << nFaces
555  << nl << endl;
556  }
557  }
558  }
559  }
560 }
561 
562 
563 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
564 
566 (
567  const dictionary& decompDict,
568  const word& regionName
569 )
570 :
571  decompositionMethod(decompDict, regionName),
572  coeffsDict_
573  (
574  findCoeffsDict
575  (
576  typeName + "Coeffs",
577  (selectionType::EXACT | selectionType::MANDATORY)
578  )
579  ),
580  methodsDict_(),
581  methods_()
582 {
583  createMethodsDict();
584  setMethods();
585 }
586 
587 
588 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
589 
591 {
592  for (const decompositionMethod& meth : methods_)
593  {
594  if (!meth.parallelAware())
595  {
596  return false;
597  }
598  }
599 
600  return true;
601 }
602 
603 
604 Foam::labelList Foam::multiLevelDecomp::decompose
605 (
606  const polyMesh& mesh,
607  const pointField& cc,
608  const scalarField& cWeights
609 ) const
610 {
611  labelList finalDecomp(cc.size(), Foam::zero{});
612  labelList cellMap(Foam::identity(cc.size()));
613 
614  CompactListList<label> cellCells;
616  (
617  mesh,
618  cellMap,
619  cellMap.size(),
620  true, // Global mesh connectivity
621  cellCells
622  );
623 
624  decompose
625  (
626  cellCells.unpack(),
627  cc,
628  cWeights,
629  cellMap, // map back to original cell points
630  0, // currLevel
631  0, // leafOffset
632 
633  finalDecomp
634  );
635 
636  return finalDecomp;
637 }
638 
639 
640 Foam::labelList Foam::multiLevelDecomp::decompose
641 (
642  const CompactListList<label>& globalCellCells,
643  const pointField& cc,
644  const scalarField& cWeights
645 ) const
646 {
647  labelList finalDecomp(cc.size(), Foam::zero{});
648  labelList cellMap(Foam::identity(cc.size()));
649 
650  decompose
651  (
652  globalCellCells.unpack(),
653  cc,
654  cWeights,
655  cellMap, // map back to original cell points
656  0, // currLevel
657  0, // leafOffset
658 
659  finalDecomp
660  );
661 
662  return finalDecomp;
663 }
664 
665 
666 Foam::labelList Foam::multiLevelDecomp::decompose
667 (
668  const labelListList& globalCellCells,
669  const pointField& cc,
670  const scalarField& cWeights
671 ) const
672 {
673  labelList finalDecomp(cc.size(), Foam::zero{});
674  labelList cellMap(Foam::identity(cc.size()));
675 
676  decompose
677  (
678  globalCellCells,
679  cc,
680  cWeights,
681  cellMap, // map back to original cell points
682  0, // currLevel
683  0, // leafOffset
684 
685  finalDecomp
686  );
687 
688  return finalDecomp;
689 }
690 
691 
692 // ************************************************************************* //
label nPatches
Definition: readKivaGrid.H:396
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
Macros for easy insertion into run-time selection tables.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
label size() const noexcept
The number of elements in table.
Definition: HashTable.H:358
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const pointField & points
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
Definition: labelLists.C:44
List< SubListType > unpack() const
Return non-compact list of lists.
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const string & prefix() const noexcept
Return the stream prefix.
Abstract base class for domain decomposition.
String literal.
Definition: keyType.H:82
static const dictionary & findCoeffsDict(const dictionary &dict, const word &coeffsName, int select=selectionType::DEFAULT)
Locate coeffsName dictionary or the fallback "coeffs" dictionary within an enclosing dictionary...
label nDomains() const noexcept
Number of domains.
A packed storage of objects of type <T> using an offset table for access.
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition: ListOps.C:30
static autoPtr< decompositionMethod > New(const dictionary &decompDict, const word &regionName="")
Return a reference to the selected decomposition method, optionally region-specific.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
multiLevelDecomp(const multiLevelDecomp &)=delete
No copy construct.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
label n
static void calcCellCells(const polyMesh &mesh, const labelUList &agglom, const label nLocalCoarse, const bool parallel, CompactListList< label > &cellCells)
Determine (local or global) cellCells from mesh agglomeration.
virtual bool parallelAware() const
Is parallel aware when all sub-methods are also parallel-aware.
List< label > labelList
A List of labels.
Definition: List.H:62
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
static void listCombineReduce(List< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Combines List elements. After completion all processors have the same data.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225