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-2023 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
54 struct vectorLessOp
55 {
56  const UList<vector>& values;
58 
59  vectorLessOp(const UList<vector>& list, direction cmpt = vector::X)
60  :
61  values(list),
62  sortCmpt(cmpt)
63  {}
64 
65  void setComponent(direction cmpt)
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 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
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 namespace Foam
88 {
89 
90 static void assignToProcessorGroup
91 (
92  labelList& processorGroup,
93  const label nProcGroup
94 )
95 {
96  const label jump = processorGroup.size()/nProcGroup;
97  const label jumpb = jump + 1;
98  const label fstProcessorGroup = processorGroup.size() - jump*nProcGroup;
99 
100  label ind = 0;
101  label j = 0;
102 
103  // assign cells to the first few processor groups (those with
104  // one extra cell each
105  for (j=0; j<fstProcessorGroup; j++)
106  {
107  for (label k=0; k<jumpb; k++)
108  {
109  processorGroup[ind++] = j;
110  }
111  }
112 
113  // and now to the `normal' processor groups
114  for (; j<nProcGroup; j++)
115  {
116  for (label k=0; k<jump; k++)
117  {
118  processorGroup[ind++] = j;
119  }
120  }
121 }
122 
123 
124 static void assignToProcessorGroup
125 (
126  labelList& processorGroup,
127  const label nProcGroup,
128  const labelList& indices,
129  const scalarField& weights,
130  const scalar summedWeights
131 )
132 {
133  // This routine gets the sorted points.
134  // Easiest to explain with an example.
135  // E.g. 400 points, sum of weights : 513.
136  // Now with number of divisions in this direction (nProcGroup) : 4
137  // gives the split at 513/4 = 128
138  // So summed weight from 0..128 goes into bin 0,
139  // ,, 128..256 goes into bin 1
140  // etc.
141  // Finally any remaining ones go into the last bin (3).
142 
143  const scalar jump = summedWeights/nProcGroup;
144  const label nProcGroupM1 = nProcGroup - 1;
145 
146  scalar sumWeights = 0;
147  label ind = 0;
148  label j = 0;
149 
150  // assign cells to all except last group.
151  for (j=0; j<nProcGroupM1; j++)
152  {
153  const scalar limit = jump*scalar(j + 1);
154  while (sumWeights < limit)
155  {
156  sumWeights += weights[indices[ind]];
157  processorGroup[ind++] = j;
158  }
159  }
160  // Ensure last included.
161  while (ind < processorGroup.size())
162  {
163  processorGroup[ind++] = nProcGroupM1;
164  }
165 }
166 
167 } // End namespace Foam
168 
169 
170 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
171 
172 Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
173 (
174  const pointField& points
175 ) const
176 {
177  // construct a list for the final result
178  labelList finalDecomp(points.size());
179 
180  labelList processorGroups(points.size());
181 
182  labelList pointIndices(identity(points.size()));
183 
184  const pointField rotatedPoints(adjustPoints(points));
185 
186  vectorLessOp sorter(rotatedPoints);
187 
188  // and one to take the processor group id's. For each direction.
189  // we assign the processors to groups of processors labelled
190  // 0..nX to give a banded structure on the mesh. Then we
191  // construct the actual processor number by treating this as
192  // the units part of the processor number.
193 
194  sorter.setComponent(vector::X);
195  Foam::sort(pointIndices, sorter);
196 
197  assignToProcessorGroup(processorGroups, n_.x());
198 
199  forAll(points, i)
200  {
201  finalDecomp[pointIndices[i]] = processorGroups[i];
202  }
203 
204 
205  // now do the same thing in the Y direction. These processor group
206  // numbers add multiples of nX to the proc. number (columns)
207 
208  sorter.setComponent(vector::Y);
209  Foam::sort(pointIndices, sorter);
210 
211  assignToProcessorGroup(processorGroups, n_.y());
212 
213  forAll(points, i)
214  {
215  finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
216  }
217 
218 
219  // finally in the Z direction. Now we add multiples of nX*nY to give
220  // layers
221 
222  sorter.setComponent(vector::Z);
223  Foam::sort(pointIndices, sorter);
224 
225  assignToProcessorGroup(processorGroups, n_.z());
226 
227  forAll(points, i)
228  {
229  finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
230  }
231 
232  return finalDecomp;
233 }
234 
235 
236 Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
237 (
238  const pointField& points,
239  const scalarField& weights
240 ) const
241 {
242  if (weights.empty())
243  {
244  return decomposeOneProc(points);
245  }
246 
247  // construct a list for the final result
248  labelList finalDecomp(points.size());
249 
250  labelList processorGroups(points.size());
251 
252  labelList pointIndices(identity(points.size()));
253 
254  const pointField rotatedPoints(adjustPoints(points));
255 
256  vectorLessOp sorter(rotatedPoints);
257 
258  // and one to take the processor group id's. For each direction.
259  // we assign the processors to groups of processors labelled
260  // 0..nX to give a banded structure on the mesh. Then we
261  // construct the actual processor number by treating this as
262  // the units part of the processor number.
263 
264  sorter.setComponent(vector::X);
265  Foam::sort(pointIndices, sorter);
266 
267  const scalar summedWeights = sum(weights);
269  (
270  processorGroups,
271  n_.x(),
272  pointIndices,
273  weights,
274  summedWeights
275  );
276 
277  forAll(points, i)
278  {
279  finalDecomp[pointIndices[i]] = processorGroups[i];
280  }
281 
282 
283  // now do the same thing in the Y direction. These processor group
284  // numbers add multiples of nX to the proc. number (columns)
285 
286  sorter.setComponent(vector::Y);
287  Foam::sort(pointIndices, sorter);
288 
290  (
291  processorGroups,
292  n_.y(),
293  pointIndices,
294  weights,
295  summedWeights
296  );
297 
298  forAll(points, i)
299  {
300  finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
301  }
302 
303 
304  // finally in the Z direction. Now we add multiples of nX*nY to give
305  // layers
306 
307  sorter.setComponent(vector::Z);
308  Foam::sort(pointIndices, sorter);
309 
311  (
312  processorGroups,
313  n_.z(),
314  pointIndices,
315  weights,
316  summedWeights
317  );
318 
319  forAll(points, i)
320  {
321  finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
322  }
324  return finalDecomp;
325 }
326 
327 
328 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
329 
331 :
332  geomDecomp(divisions)
333 {}
334 
335 
337 (
338  const dictionary& decompDict,
339  const word& regionName
340 )
341 :
342  geomDecomp(typeName, decompDict, regionName)
343 {}
344 
345 
346 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
347 
349 (
350  const pointField& points,
351  const scalarField& weights
352 ) const
353 {
354  // Uniform weighting if weights are empty or poorly sized
355  const bool hasWeights = returnReduceAnd(points.size() == weights.size());
356 
357  if (!UPstream::parRun())
358  {
359  if (hasWeights)
360  {
361  return decomposeOneProc(points, weights);
362  }
363  else
364  {
365  return decomposeOneProc(points);
366  }
367  }
368  else
369  {
370  // Parallel
371  const globalIndex globalNumbers(points.size());
372 
373  labelList allDecomp;
374  pointField allPoints(globalNumbers.gather(points));
375  scalarField allWeights;
376 
377  if (hasWeights)
378  {
379  // Non-uniform weighting
380  allWeights = globalNumbers.gather(weights);
381  }
382 
383  if (UPstream::master())
384  {
385  if (hasWeights)
386  {
387  allDecomp = decomposeOneProc(allPoints, allWeights);
388  }
389  else
390  {
391  allDecomp = decomposeOneProc(allPoints);
392  }
393  allPoints.clear(); // Not needed anymore
394  allWeights.clear(); // ...
395  }
396 
397  return globalNumbers.scatter(allDecomp);
398  }
399 }
400 
401 
402 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
simpleGeomDecomp(const simpleGeomDecomp &)=delete
No copy construct.
uint8_t direction
Definition: direction.H:46
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
const UList< vector > & values
bool operator()(const label a, const label b) const
Vector< label > n_
The divisions.
Definition: geomDecomp.H:122
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
virtual labelList decompose(const pointField &points, const scalarField &weights=scalarField::null()) const
Return for every coordinate the wanted processor number. using uniform or specified point weights...
const Cmpt & y() const noexcept
Access to the vector y component.
Definition: Vector.H:140
label k
Boltzmann constant.
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
Base for geometrical domain decomposition methods.
Definition: geomDecomp.H:82
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
const dimensionedScalar b
Wien displacement law constant: default SI units: [m.K].
Definition: createFields.H:27
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
complex limit(const complex &c1, const complex &c2)
Definition: complexI.H:235
static void assignToProcessorGroup(labelList &processorGroup, const label nProcGroup)
void setComponent(direction cmpt)
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:105
const Cmpt & x() const noexcept
Access to the vector x component.
Definition: Vector.H:135
defineTypeNameAndDebug(combustionModel, 0)
tmp< pointField > adjustPoints(const pointField &) const
Apply delta (jitter) or rotation to coordinates.
Definition: geomDecomp.C:115
const Cmpt & z() const noexcept
Access to the vector z component.
Definition: Vector.H:145
vectorLessOp(const UList< vector > &list, direction cmpt=vector::X)
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
PtrList< volScalarField > & Y
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1082
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
List< label > labelList
A List of labels.
Definition: List.H:62
Namespace for OpenFOAM.
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)