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