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 | Copyright (C) 2011-2015 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 "multiLevelDecomp.H"
28 #include "IFstream.H"
29 #include "globalIndex.H"
30 #include "mapDistribute.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  defineTypeNameAndDebug(multiLevelDecomp, 0);
37 
39  (
40  decompositionMethod,
41  multiLevelDecomp,
42  dictionary
43  );
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 // Given a subset of cells determine the new global indices. The problem
50 // is in the cells from neighbouring processors which need to be renumbered.
52 (
53  const label nDomains,
54  const label domainI,
55  const labelList& dist,
56 
57  const labelListList& cellCells,
58  const labelList& set,
59  labelListList& subCellCells,
60  labelList& cutConnections
61 ) const
62 {
63  // Determine new index for cells by inverting subset
64  labelList oldToNew(invert(cellCells.size(), set));
65 
66  globalIndex globalCells(cellCells.size());
67 
68  // Subset locally the elements for which I have data
69  subCellCells = UIndirectList<labelList>(cellCells, set);
70 
71  // Get new indices for neighbouring processors
72  List<Map<label> > compactMap;
73  mapDistribute map(globalCells, subCellCells, compactMap);
74  map.distribute(oldToNew);
75  labelList allDist(dist);
76  map.distribute(allDist);
77 
78  // Now we have:
79  // oldToNew : the locally-compact numbering of all our cellCells. -1 if
80  // cellCell is not in set.
81  // allDist : destination domain for all our cellCells
82  // subCellCells : indexes into oldToNew and allDist
83 
84  // Globally compact numbering for cells in set.
85  globalIndex globalSubCells(set.size());
86 
87  // Now subCellCells contains indices into oldToNew which are the
88  // new locations of the neighbouring cells.
89 
90  cutConnections.setSize(nDomains);
91  cutConnections = 0;
92 
93  forAll(subCellCells, subCellI)
94  {
95  labelList& cCells = subCellCells[subCellI];
96 
97  // Keep the connections to valid mapped cells
98  label newI = 0;
99  forAll(cCells, i)
100  {
101  // Get locally-compact cell index of neighbouring cell
102  label nbrCellI = oldToNew[cCells[i]];
103  if (nbrCellI == -1)
104  {
105  cutConnections[allDist[cCells[i]]]++;
106  }
107  else
108  {
109  // Reconvert local cell index into global one
110 
111  // Get original neighbour
112  label cellI = set[subCellI];
113  label oldNbrCellI = cellCells[cellI][i];
114  // Get processor from original neighbour
115  label procI = globalCells.whichProcID(oldNbrCellI);
116  // Convert into global compact numbering
117  cCells[newI++] = globalSubCells.toGlobal(procI, nbrCellI);
118  }
119  }
120  cCells.setSize(newI);
121  }
122 }
123 
124 
126 (
127  const labelListList& pointPoints,
128  const pointField& points,
129  const scalarField& pointWeights,
130  const labelList& pointMap, // map back to original points
131  const label levelI,
132 
133  labelField& finalDecomp
134 )
135 {
136  labelList dist
137  (
138  methods_[levelI].decompose
139  (
140  pointPoints,
141  points,
142  pointWeights
143  )
144  );
145 
146  forAll(pointMap, i)
147  {
148  label orig = pointMap[i];
149  finalDecomp[orig] += dist[i];
150  }
151 
152  if (levelI != methods_.size()-1)
153  {
154  // Recurse
155 
156  // Determine points per domain
157  label n = methods_[levelI].nDomains();
158  labelListList domainToPoints(invertOneToMany(n, dist));
159 
160  // 'Make space' for new levels of decomposition
161  finalDecomp *= methods_[levelI+1].nDomains();
162 
163  // Extract processor+local index from point-point addressing
164  if (debug && Pstream::master())
165  {
166  Pout<< "Decomposition at level " << levelI << " :" << endl;
167  }
168 
169  for (label domainI = 0; domainI < n; domainI++)
170  {
171  // Extract elements for current domain
172  const labelList domainPoints(findIndices(dist, domainI));
173 
174  // Subset point-wise data.
175  pointField subPoints(points, domainPoints);
176  scalarField subWeights(pointWeights, domainPoints);
177  labelList subPointMap(UIndirectList<label>(pointMap, domainPoints));
178  // Subset point-point addressing (adapt global numbering)
179  labelListList subPointPoints;
180  labelList nOutsideConnections;
181  subsetGlobalCellCells
182  (
183  n,
184  domainI,
185  dist,
186 
187  pointPoints,
188  domainPoints,
189 
190  subPointPoints,
191  nOutsideConnections
192  );
193 
194  label nPoints = returnReduce(domainPoints.size(), plusOp<label>());
195  Pstream::listCombineGather(nOutsideConnections, plusEqOp<label>());
196  Pstream::listCombineScatter(nOutsideConnections);
197  label nPatches = 0;
198  label nFaces = 0;
199  forAll(nOutsideConnections, i)
200  {
201  if (nOutsideConnections[i] > 0)
202  {
203  nPatches++;
204  nFaces += nOutsideConnections[i];
205  }
206  }
207 
208  string oldPrefix;
209  if (debug && Pstream::master())
210  {
211  Pout<< " Domain " << domainI << nl
212  << " Number of cells = " << nPoints << nl
213  << " Number of inter-domain patches = " << nPatches
214  << nl
215  << " Number of inter-domain faces = " << nFaces << nl
216  << endl;
217  oldPrefix = Pout.prefix();
218  Pout.prefix() = " " + oldPrefix;
219  }
220 
221  decompose
222  (
223  subPointPoints,
224  subPoints,
225  subWeights,
226  subPointMap,
227  levelI+1,
228 
229  finalDecomp
230  );
231  if (debug && Pstream::master())
232  {
233  Pout.prefix() = oldPrefix;
234  }
235  }
236 
237 
238  if (debug)
239  {
240  // Do straight decompose of two levels
241  label nNext = methods_[levelI+1].nDomains();
242  label nTotal = n*nNext;
243 
244  // Retrieve original level0 dictionary and modify number of domains
245  dictionary::const_iterator iter =
246  decompositionDict_.subDict(typeName + "Coeffs").begin();
247  dictionary myDict = iter().dict();
248  myDict.set("numberOfSubdomains", nTotal);
249 
250  if (debug && Pstream::master())
251  {
252  Pout<< "Reference decomposition with " << myDict << " :"
253  << endl;
254  }
255 
257  (
258  myDict
259  );
260  labelList dist
261  (
262  method0().decompose
263  (
264  pointPoints,
265  points,
266  pointWeights
267  )
268  );
269 
270  for (label blockI = 0; blockI < n; blockI++)
271  {
272  // Count the number inbetween blocks of nNext size
273 
274  label nPoints = 0;
275  labelList nOutsideConnections(n, 0);
276  forAll(pointPoints, pointI)
277  {
278  if ((dist[pointI] / nNext) == blockI)
279  {
280  nPoints++;
281 
282  const labelList& pPoints = pointPoints[pointI];
283 
284  forAll(pPoints, i)
285  {
286  label distBlockI = dist[pPoints[i]] / nNext;
287  if (distBlockI != blockI)
288  {
289  nOutsideConnections[distBlockI]++;
290  }
291  }
292  }
293  }
294 
296  Pstream::listCombineGather
297  (
298  nOutsideConnections,
300  );
301  Pstream::listCombineScatter(nOutsideConnections);
302  label nPatches = 0;
303  label nFaces = 0;
304  forAll(nOutsideConnections, i)
305  {
306  if (nOutsideConnections[i] > 0)
307  {
308  nPatches++;
309  nFaces += nOutsideConnections[i];
310  }
311  }
312 
313  if (debug && Pstream::master())
314  {
315  Pout<< " Domain " << blockI << nl
316  << " Number of cells = " << nPoints << nl
317  << " Number of inter-domain patches = "
318  << nPatches << nl
319  << " Number of inter-domain faces = " << nFaces
320  << nl << endl;
321  }
322  }
323  }
324  }
325 }
326 
327 
328 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
329 
331 :
332  decompositionMethod(decompositionDict),
333  methodsDict_(decompositionDict_.subDict(typeName + "Coeffs"))
334 {
335  methods_.setSize(methodsDict_.size());
336  label i = 0;
338  {
339  methods_.set(i++, decompositionMethod::New(iter().dict()));
340  }
341 
342  label n = 1;
343  Info<< "decompositionMethod " << type() << " :" << endl;
344  forAll(methods_, i)
345  {
346  Info<< " level " << i << " decomposing with " << methods_[i].type()
347  << " into " << methods_[i].nDomains() << " subdomains." << endl;
348 
349  n *= methods_[i].nDomains();
350  }
351 
352  if (n != nDomains())
353  {
355  << "Top level decomposition specifies " << nDomains()
356  << " domains which is not equal to the product of"
357  << " all sub domains " << n
358  << exit(FatalError);
359  }
360 }
361 
362 
363 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
364 
366 {
367  forAll(methods_, i)
368  {
369  if (!methods_[i].parallelAware())
370  {
371  return false;
372  }
373  }
374  return true;
375 }
376 
377 
379 (
380  const polyMesh& mesh,
381  const pointField& cc,
382  const scalarField& cWeights
383 )
384 {
385  CompactListList<label> cellCells;
386  calcCellCells(mesh, identity(cc.size()), cc.size(), true, cellCells);
387 
388  labelField finalDecomp(cc.size(), 0);
389  labelList cellMap(identity(cc.size()));
390 
391  decompose
392  (
393  cellCells(),
394  cc,
395  cWeights,
396  cellMap, // map back to original cells
397  0,
398 
399  finalDecomp
400  );
401 
402  return finalDecomp;
403 }
404 
405 
407 (
408  const labelListList& globalPointPoints,
409  const pointField& points,
410  const scalarField& pointWeights
411 )
412 {
413  labelField finalDecomp(points.size(), 0);
414  labelList pointMap(identity(points.size()));
415 
416  decompose
417  (
418  globalPointPoints,
419  points,
420  pointWeights,
421  pointMap, // map back to original points
422  0,
423 
424  finalDecomp
425  );
426 
427  return finalDecomp;
428 }
429 
430 
431 // ************************************************************************* //
multiLevelDecomp.H
Foam::multiLevelDecomp::multiLevelDecomp
multiLevelDecomp(const multiLevelDecomp &)
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
nPatches
label nPatches
Definition: readKivaGrid.H:402
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Foam::decompositionMethod::New
static autoPtr< decompositionMethod > New(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
Definition: decompositionMethod.C:179
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::CompactListList
A packed storage unstructured matrix of objects of type <T> using an offset table for access.
Definition: polyTopoChange.H:91
globalIndex.H
Foam::plusOp
Definition: ops.H:164
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::multiLevelDecomp::parallelAware
virtual bool parallelAware() const
Is method parallel aware (i.e. does it synchronize domains across.
Definition: multiLevelDecomp.C:365
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::multiLevelDecomp::decompose
void decompose(const labelListList &pointPoints, const pointField &points, const scalarField &pointWeights, const labelList &pointMap, const label levelI, labelField &finalDecomp)
Decompose level methodI without addressing.
Definition: multiLevelDecomp.C:126
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
IFstream.H
Foam::mapDistribute::distribute
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
Definition: mapDistributeTemplates.C:155
Foam::identity
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::decompositionMethod::nDomains
label nDomains() const
Definition: decompositionMethod.H:112
Foam::multiLevelDecomp::subsetGlobalCellCells
void subsetGlobalCellCells(const label nDomains, const label domainI, const labelList &dist, const labelListList &cellCells, const labelList &set, labelListList &subCellCells, labelList &cutConnections) const
Given connectivity across processors work out connectivity.
Definition: multiLevelDecomp.C:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::plusEqOp
Definition: ops.H:71
Foam::decompositionMethod
Abstract base class for decomposition.
Definition: decompositionMethod.H:48
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:117
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
mapDistribute.H
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::multiLevelDecomp::methodsDict_
dictionary methodsDict_
Definition: multiLevelDecomp.H:52
Foam::globalIndex::whichProcID
label whichProcID(const label i) const
Which processor does global come from? Binary search.
Definition: globalIndexI.H:123
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::invert
labelList invert(const label len, const labelUList &)
Invert one-to-one map. Unmapped elements will be -1.
Definition: ListOps.C:37
Foam::multiLevelDecomp::methods_
PtrList< decompositionMethod > methods_
Definition: multiLevelDecomp.H:54
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::findIndices
labelList findIndices(const ListType &, typename ListType::const_reference, const label start=0)
Find all occurences of given element. Linear search.
Foam::prefixOSstream::prefix
const string & prefix() const
Return the prefix of the stream.
Definition: prefixOSstream.H:86
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
Foam::dictionary::set
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:856