Loading...
Searching...
No Matches
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-2025 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 "multiLevelDecomp.H"
31#include "IFstream.H"
32#include "globalIndex.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
42 (
46 );
47}
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
52void 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",
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
191void 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.
244void 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
318void 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 Pstream::listReduce(nOutsideConnections, sumOp<label>());
424
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
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
531 Pstream::listReduce(nOutsideConnections, sumOp<label>());
532
533 label nPatches = 0;
534 label nFaces = 0;
535 for (const label nConnect : nOutsideConnections)
536 {
537 if (nConnect > 0)
538 {
539 ++nPatches;
540 nFaces += nConnect;
541 }
542 }
543
544 if (debug && Pstream::master())
545 {
546 Pout<< " Domain " << blockI << nl
547 << " Number of cells = " << nPoints << nl
548 << " Number of inter-domain patches = "
549 << nPatches << nl
550 << " Number of inter-domain faces = " << nFaces
551 << nl << endl;
552 }
553 }
555 }
556}
557
558
559// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
560
562(
563 const dictionary& decompDict,
564 const word& regionName
565)
566:
567 decompositionMethod(decompDict, regionName),
568 coeffsDict_
569 (
570 findCoeffsDict
571 (
572 typeName + "Coeffs",
573 (selectionType::EXACT | selectionType::MANDATORY)
574 )
575 ),
576 methodsDict_(),
577 methods_()
578{
579 createMethodsDict();
580 setMethods();
581}
582
583
584// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
585
587{
588 for (const decompositionMethod& meth : methods_)
589 {
590 if (!meth.parallelAware())
591 {
592 return false;
594 }
595
596 return true;
597}
598
599
600Foam::labelList Foam::multiLevelDecomp::decompose
601(
602 const polyMesh& mesh,
603 const pointField& cc,
604 const scalarField& cWeights
605) const
606{
607 labelList finalDecomp(cc.size(), Foam::zero{});
608 labelList cellMap(Foam::identity(cc.size()));
609
610 CompactListList<label> cellCells;
612 (
613 mesh,
614 cellMap,
615 cellMap.size(),
616 true, // Global mesh connectivity
617 cellCells
618 );
619
620 decompose
621 (
622 cellCells.unpack(),
623 cc,
624 cWeights,
625 cellMap, // map back to original cell points
626 0, // currLevel
627 0, // leafOffset
628
629 finalDecomp
630 );
631
632 return finalDecomp;
633}
634
635
636Foam::labelList Foam::multiLevelDecomp::decompose
637(
638 const CompactListList<label>& globalCellCells,
639 const pointField& cc,
640 const scalarField& cWeights
641) const
642{
643 labelList finalDecomp(cc.size(), Foam::zero{});
644 labelList cellMap(Foam::identity(cc.size()));
645
646 decompose
647 (
648 globalCellCells.unpack(),
649 cc,
650 cWeights,
651 cellMap, // map back to original cell points
652 0, // currLevel
653 0, // leafOffset
654
655 finalDecomp
656 );
657
658 return finalDecomp;
659}
660
661
662Foam::labelList Foam::multiLevelDecomp::decompose
663(
664 const labelListList& globalCellCells,
665 const pointField& cc,
666 const scalarField& cWeights
667) const
668{
669 labelList finalDecomp(cc.size(), Foam::zero{});
670 labelList cellMap(Foam::identity(cc.size()));
671
672 decompose
673 (
674 globalCellCells,
675 cc,
676 cWeights,
677 cellMap, // map back to original cell points
678 0, // currLevel
679 0, // leafOffset
680
681 finalDecomp
682 );
683
684 return finalDecomp;
685}
686
687
688// ************************************************************************* //
label n
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A packed storage of objects of type <T> using an offset table for access.
List< SubListType > unpack() const
Return non-compact list of lists.
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
A List with indirect addressing. Like IndirectList but does not store addressing.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Abstract base class for domain decomposition.
selectionType
Selection type when handling the coefficients dictionary.
@ EXACT
No fallback to "coeffs" if main name not found.
@ MANDATORY
Fatal if dictionary could not be found.
decompositionMethod(const label numDomains)
Construct with specified number of domains, no coefficients or constraints.
static autoPtr< decompositionMethod > New(const dictionary &decompDict, const word &regionName="")
Return a reference to the selected decomposition method, optionally region-specific.
static FOAM_NO_DANGLING_REFERENCE 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 list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
void clear()
Clear the dictionary.
Definition dictionary.C:862
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition dictionary.C:765
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
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.
@ LITERAL
String literal.
Definition keyType.H:82
Class containing processor-to-processor mapping information.
Decompose given using consecutive application of decomposers.
multiLevelDecomp(const multiLevelDecomp &)=delete
No copy construct.
virtual bool parallelAware() const
Is parallel aware when all sub-methods are also parallel-aware.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
dynamicFvMesh & mesh
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const pointField & points
label nPoints
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
UIndirectList< label > labelUIndList
UIndirectList of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
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...
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
labelList findIndices(const ListType &input, typename ListType::const_reference val, label start=0)
Linear search to find all occurrences of given element.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
label nPatches
dictionary dict
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299