simpleGeomDecomp.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-2021 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 "simpleGeomDecomp.H"
31 #include "globalIndex.H"
32 #include "SubField.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(simpleGeomDecomp, 0);
40  (
41  decompositionMethod,
42  simpleGeomDecomp,
43  dictionary
44  );
45 }
46 
47 
48 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
49 
50 namespace Foam
51 {
52 
53 // A list compare binary predicate for normal sort by vector component
55 {
58 
60  :
61  values(list),
62  sortCmpt(cmpt)
63  {}
64 
66  {
67  sortCmpt = cmpt;
68  }
69 
70  bool operator()(const label a, const label b) const
71  {
72  return values[a][sortCmpt] < values[b][sortCmpt];
73  }
74 };
75 
76 } // End namespace Foam
77 
78 
79 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
80 
81 // assignToProcessorGroup : given nCells cells and nProcGroup processor
82 // groups to share them, how do we share them out? Answer : each group
83 // gets nCells/nProcGroup cells, and the first few get one
84 // extra to make up the numbers. This should produce almost
85 // perfect load balancing
86 
87 void Foam::simpleGeomDecomp::assignToProcessorGroup
88 (
89  labelList& processorGroup,
90  const label nProcGroup
91 )
92 {
93  label jump = processorGroup.size()/nProcGroup;
94  label jumpb = jump + 1;
95  label fstProcessorGroup = processorGroup.size() - jump*nProcGroup;
96 
97  label ind = 0;
98  label j = 0;
99 
100  // assign cells to the first few processor groups (those with
101  // one extra cell each
102  for (j=0; j<fstProcessorGroup; j++)
103  {
104  for (label k=0; k<jumpb; k++)
105  {
106  processorGroup[ind++] = j;
107  }
108  }
109 
110  // and now to the `normal' processor groups
111  for (; j<nProcGroup; j++)
112  {
113  for (label k=0; k<jump; k++)
114  {
115  processorGroup[ind++] = j;
116  }
117  }
118 }
119 
120 
121 void Foam::simpleGeomDecomp::assignToProcessorGroup
122 (
123  labelList& processorGroup,
124  const label nProcGroup,
125  const labelList& indices,
126  const scalarField& weights,
127  const scalar summedWeights
128 )
129 {
130  // This routine gets the sorted points.
131  // Easiest to explain with an example.
132  // E.g. 400 points, sum of weights : 513.
133  // Now with number of divisions in this direction (nProcGroup) : 4
134  // gives the split at 513/4 = 128
135  // So summed weight from 0..128 goes into bin 0,
136  // ,, 128..256 goes into bin 1
137  // etc.
138  // Finally any remaining ones go into the last bin (3).
139 
140  const scalar jump = summedWeights/nProcGroup;
141  const label nProcGroupM1 = nProcGroup - 1;
142  scalar sumWeights = 0;
143  label ind = 0;
144  label j = 0;
145 
146  // assign cells to all except last group.
147  for (j=0; j<nProcGroupM1; j++)
148  {
149  const scalar limit = jump*scalar(j + 1);
150  while (sumWeights < limit)
151  {
152  sumWeights += weights[indices[ind]];
153  processorGroup[ind++] = j;
154  }
155  }
156  // Ensure last included.
157  while (ind < processorGroup.size())
158  {
159  processorGroup[ind++] = nProcGroupM1;
160  }
161 }
162 
163 
164 Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
165 (
166  const pointField& points
167 ) const
168 {
169  // construct a list for the final result
170  labelList finalDecomp(points.size());
171 
172  labelList processorGroups(points.size());
173 
174  labelList pointIndices(identity(points.size()));
175 
176  const pointField rotatedPoints(adjustPoints(points));
177 
178  vectorLessOp sorter(rotatedPoints);
179 
180  // and one to take the processor group id's. For each direction.
181  // we assign the processors to groups of processors labelled
182  // 0..nX to give a banded structure on the mesh. Then we
183  // construct the actual processor number by treating this as
184  // the units part of the processor number.
185 
186  sorter.setComponent(vector::X);
187  Foam::sort(pointIndices, sorter);
188 
189  assignToProcessorGroup(processorGroups, n_.x());
190 
191  forAll(points, i)
192  {
193  finalDecomp[pointIndices[i]] = processorGroups[i];
194  }
195 
196 
197  // now do the same thing in the Y direction. These processor group
198  // numbers add multiples of nX to the proc. number (columns)
199 
200  sorter.setComponent(vector::Y);
201  Foam::sort(pointIndices, sorter);
202 
203  assignToProcessorGroup(processorGroups, n_.y());
204 
205  forAll(points, i)
206  {
207  finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
208  }
209 
210 
211  // finally in the Z direction. Now we add multiples of nX*nY to give
212  // layers
213 
214  sorter.setComponent(vector::Z);
215  Foam::sort(pointIndices, sorter);
216 
217  assignToProcessorGroup(processorGroups, n_.z());
218 
219  forAll(points, i)
220  {
221  finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
222  }
223 
224  return finalDecomp;
225 }
226 
227 
228 Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
229 (
230  const pointField& points,
231  const scalarField& weights
232 ) const
233 {
234  // construct a list for the final result
235  labelList finalDecomp(points.size());
236 
237  labelList processorGroups(points.size());
238 
239  labelList pointIndices(identity(points.size()));
240 
241  const pointField rotatedPoints(adjustPoints(points));
242 
243  vectorLessOp sorter(rotatedPoints);
244 
245  // and one to take the processor group id's. For each direction.
246  // we assign the processors to groups of processors labelled
247  // 0..nX to give a banded structure on the mesh. Then we
248  // construct the actual processor number by treating this as
249  // the units part of the processor number.
250 
251  sorter.setComponent(vector::X);
252  Foam::sort(pointIndices, sorter);
253 
254  const scalar summedWeights = sum(weights);
255  assignToProcessorGroup
256  (
257  processorGroups,
258  n_.x(),
259  pointIndices,
260  weights,
261  summedWeights
262  );
263 
264  forAll(points, i)
265  {
266  finalDecomp[pointIndices[i]] = processorGroups[i];
267  }
268 
269 
270  // now do the same thing in the Y direction. These processor group
271  // numbers add multiples of nX to the proc. number (columns)
272 
273  sorter.setComponent(vector::Y);
274  Foam::sort(pointIndices, sorter);
275 
276  assignToProcessorGroup
277  (
278  processorGroups,
279  n_.y(),
280  pointIndices,
281  weights,
282  summedWeights
283  );
284 
285  forAll(points, i)
286  {
287  finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
288  }
289 
290 
291  // finally in the Z direction. Now we add multiples of nX*nY to give
292  // layers
293 
294  sorter.setComponent(vector::Z);
295  Foam::sort(pointIndices, sorter);
296 
297  assignToProcessorGroup
298  (
299  processorGroups,
300  n_.z(),
301  pointIndices,
302  weights,
303  summedWeights
304  );
305 
306  forAll(points, i)
307  {
308  finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
309  }
310 
311  return finalDecomp;
312 }
313 
314 
315 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
316 
318 (
319  const dictionary& decompDict,
320  const word& regionName
321 )
322 :
323  geomDecomp(typeName, decompDict, regionName)
324 {}
325 
326 
327 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
328 
330 (
331  const pointField& points
332 ) const
333 {
334  if (!Pstream::parRun())
335  {
336  return decomposeOneProc(points);
337  }
338  else
339  {
340  globalIndex globalNumbers(points.size());
341 
342  // Collect all points on master
343  if (Pstream::master())
344  {
345  pointField allPoints(globalNumbers.size());
346 
347  label nTotalPoints = 0;
348  // Master first
350  nTotalPoints += points.size();
351 
352  // Add received
353  for (const int subproci : Pstream::subProcs())
354  {
355  IPstream fromProc(Pstream::commsTypes::scheduled, subproci);
356  pointField nbrPoints(fromProc);
358  (
359  allPoints,
360  nbrPoints.size(),
361  nTotalPoints
362  ) = nbrPoints;
363  nTotalPoints += nbrPoints.size();
364  }
365 
366  // Decompose
367  labelList finalDecomp(decomposeOneProc(allPoints));
368 
369  // Send back
370  for (const int subproci : Pstream::subProcs())
371  {
372  OPstream toProc(Pstream::commsTypes::scheduled, subproci);
373  toProc << SubField<label>
374  (
375  finalDecomp,
376  globalNumbers.localSize(subproci),
377  globalNumbers.offset(subproci)
378  );
379  }
380  // Get my own part
381  finalDecomp.setSize(points.size());
382 
383  return finalDecomp;
384  }
385  else
386  {
387  // Send my points
388  {
389  OPstream toMaster
390  (
393  );
394  toMaster<< points;
395  }
396 
397  // Receive back decomposition
398  IPstream fromMaster
399  (
402  );
403  labelList finalDecomp(fromMaster);
404 
405  return finalDecomp;
406  }
407  }
408 }
409 
410 
412 (
413  const pointField& points,
414  const scalarField& weights
415 ) const
416 {
417  if (!Pstream::parRun())
418  {
419  return decomposeOneProc(points, weights);
420  }
421  else
422  {
423  globalIndex globalNumbers(points.size());
424 
425  // Collect all points on master
426  if (Pstream::master())
427  {
428  pointField allPoints(globalNumbers.size());
429  scalarField allWeights(allPoints.size());
430 
431  label nTotalPoints = 0;
432  // Master first
434  SubField<scalar>(allWeights, points.size()) = weights;
435  nTotalPoints += points.size();
436 
437  // Add received
438  for (const int subproci : Pstream::subProcs())
439  {
440  IPstream fromProc(Pstream::commsTypes::scheduled, subproci);
441  pointField nbrPoints(fromProc);
442  scalarField nbrWeights(fromProc);
444  (
445  allPoints,
446  nbrPoints.size(),
447  nTotalPoints
448  ) = nbrPoints;
450  (
451  allWeights,
452  nbrWeights.size(),
453  nTotalPoints
454  ) = nbrWeights;
455  nTotalPoints += nbrPoints.size();
456  }
457 
458  // Decompose
459  labelList finalDecomp(decomposeOneProc(allPoints, allWeights));
460 
461  // Send back
462  for (const int subproci : Pstream::subProcs())
463  {
464  OPstream toProc(Pstream::commsTypes::scheduled, subproci);
465  toProc << SubField<label>
466  (
467  finalDecomp,
468  globalNumbers.localSize(subproci),
469  globalNumbers.offset(subproci)
470  );
471  }
472  // Get my own part
473  finalDecomp.setSize(points.size());
474 
475  return finalDecomp;
476  }
477  else
478  {
479  // Send my points
480  {
481  OPstream toMaster
482  (
485  );
486  toMaster<< points << weights;
487  }
488 
489  // Receive back decomposition
490  IPstream fromMaster
491  (
494  );
495  labelList finalDecomp(fromMaster);
496 
497  return finalDecomp;
498  }
499  }
500 }
501 
502 
503 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::vectorLessOp::sortCmpt
direction sortCmpt
Definition: simpleGeomDecomp.C:57
Foam::UPstream::masterNo
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:451
Foam::Vector< scalar >::Z
Definition: Vector.H:81
Foam::Vector< scalar >::Y
Definition: Vector.H:81
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
SubField.H
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:53
Foam::simpleGeomDecomp::decompose
virtual labelList decompose(const pointField &points) const
Decompose with uniform weights.
Definition: simpleGeomDecomp.C:330
globalIndex.H
Foam::vectorLessOp::values
const UList< vector > & values
Definition: simpleGeomDecomp.C:56
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::vectorLessOp::vectorLessOp
vectorLessOp(const UList< vector > &list, direction cmpt=vector::X)
Definition: simpleGeomDecomp.C:59
Foam::simpleGeomDecomp::simpleGeomDecomp
simpleGeomDecomp(const simpleGeomDecomp &)=delete
No copy assignment.
Foam::vectorLessOp::setComponent
void setComponent(direction cmpt)
Definition: simpleGeomDecomp.C:65
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::SubField
SubField is a Field obtained as a section of another Field, without its own allocation....
Definition: Field.H:64
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::geomDecomp
Base for geometrical domain decomposition methods.
Definition: geomDecomp.H:86
Foam::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:515
Foam::constant::physicoChemical::b
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
Foam::Field< vector >
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::Vector< scalar >::X
Definition: Vector.H:81
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:261
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::vectorLessOp::operator()
bool operator()(const label a, const label b) const
Definition: simpleGeomDecomp.C:70
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::UPstream::commsTypes::scheduled
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:103
points
const pointField & points
Definition: gmvOutputHeader.H:1
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::limit
complex limit(const complex &, const complex &)
Definition: complexI.H:263
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::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::direction
uint8_t direction
Definition: direction.H:52
Foam::vectorLessOp
Definition: simpleGeomDecomp.C:54
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:53
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
simpleGeomDecomp.H