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-2020 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 "SortableList.H"
32 #include "globalIndex.H"
33 #include "SubField.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(simpleGeomDecomp, 0);
40 
42  (
43  decompositionMethod,
44  simpleGeomDecomp,
45  dictionary
46  );
47 
49  (
50  decompositionMethod,
51  simpleGeomDecomp,
52  dictionaryRegion
53  );
54 }
55 
56 
57 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59 // assignToProcessorGroup : given nCells cells and nProcGroup processor
60 // groups to share them, how do we share them out? Answer : each group
61 // gets nCells/nProcGroup cells, and the first few get one
62 // extra to make up the numbers. This should produce almost
63 // perfect load balancing
64 
65 void Foam::simpleGeomDecomp::assignToProcessorGroup
66 (
67  labelList& processorGroup,
68  const label nProcGroup
69 )
70 {
71  label jump = processorGroup.size()/nProcGroup;
72  label jumpb = jump + 1;
73  label fstProcessorGroup = processorGroup.size() - jump*nProcGroup;
74 
75  label ind = 0;
76  label j = 0;
77 
78  // assign cells to the first few processor groups (those with
79  // one extra cell each
80  for (j=0; j<fstProcessorGroup; j++)
81  {
82  for (label k=0; k<jumpb; k++)
83  {
84  processorGroup[ind++] = j;
85  }
86  }
87 
88  // and now to the `normal' processor groups
89  for (; j<nProcGroup; j++)
90  {
91  for (label k=0; k<jump; k++)
92  {
93  processorGroup[ind++] = j;
94  }
95  }
96 }
97 
98 
99 void Foam::simpleGeomDecomp::assignToProcessorGroup
100 (
101  labelList& processorGroup,
102  const label nProcGroup,
103  const labelList& indices,
104  const scalarField& weights,
105  const scalar summedWeights
106 )
107 {
108  // This routine gets the sorted points.
109  // Easiest to explain with an example.
110  // E.g. 400 points, sum of weights : 513.
111  // Now with number of divisions in this direction (nProcGroup) : 4
112  // gives the split at 513/4 = 128
113  // So summed weight from 0..128 goes into bin 0,
114  // ,, 128..256 goes into bin 1
115  // etc.
116  // Finally any remaining ones go into the last bin (3).
117 
118  const scalar jump = summedWeights/nProcGroup;
119  const label nProcGroupM1 = nProcGroup - 1;
120  scalar sumWeights = 0;
121  label ind = 0;
122  label j = 0;
123 
124  // assign cells to all except last group.
125  for (j=0; j<nProcGroupM1; j++)
126  {
127  const scalar limit = jump*scalar(j + 1);
128  while (sumWeights < limit)
129  {
130  sumWeights += weights[indices[ind]];
131  processorGroup[ind++] = j;
132  }
133  }
134  // Ensure last included.
135  while (ind < processorGroup.size())
136  {
137  processorGroup[ind++] = nProcGroupM1;
138  }
139 }
140 
141 
142 Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
143 (
144  const pointField& points
145 ) const
146 {
147  // construct a list for the final result
148  labelList finalDecomp(points.size());
149 
150  labelList processorGroups(points.size());
151 
152  labelList pointIndices(points.size());
153  forAll(pointIndices, i)
154  {
155  pointIndices[i] = i;
156  }
157 
158  const pointField rotatedPoints(rotDelta_ & points);
159 
160  // and one to take the processor group id's. For each direction.
161  // we assign the processors to groups of processors labelled
162  // 0..nX to give a banded structure on the mesh. Then we
163  // construct the actual processor number by treating this as
164  // the units part of the processor number.
165  sort
166  (
167  pointIndices,
168  UList<scalar>::less(rotatedPoints.component(vector::X))
169  );
170 
171  assignToProcessorGroup(processorGroups, n_.x());
172 
173  forAll(points, i)
174  {
175  finalDecomp[pointIndices[i]] = processorGroups[i];
176  }
177 
178 
179  // now do the same thing in the Y direction. These processor group
180  // numbers add multiples of nX to the proc. number (columns)
181  sort
182  (
183  pointIndices,
184  UList<scalar>::less(rotatedPoints.component(vector::Y))
185  );
186 
187  assignToProcessorGroup(processorGroups, n_.y());
188 
189  forAll(points, i)
190  {
191  finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
192  }
193 
194 
195  // finally in the Z direction. Now we add multiples of nX*nY to give
196  // layers
197  sort
198  (
199  pointIndices,
200  UList<scalar>::less(rotatedPoints.component(vector::Z))
201  );
202 
203  assignToProcessorGroup(processorGroups, n_.z());
204 
205  forAll(points, i)
206  {
207  finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
208  }
209 
210  return finalDecomp;
211 }
212 
213 
214 Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
215 (
216  const pointField& points,
217  const scalarField& weights
218 ) const
219 {
220  // construct a list for the final result
221  labelList finalDecomp(points.size());
222 
223  labelList processorGroups(points.size());
224 
225  labelList pointIndices(points.size());
226  forAll(pointIndices, i)
227  {
228  pointIndices[i] = i;
229  }
230 
231  const pointField rotatedPoints(rotDelta_ & points);
232 
233  // and one to take the processor group id's. For each direction.
234  // we assign the processors to groups of processors labelled
235  // 0..nX to give a banded structure on the mesh. Then we
236  // construct the actual processor number by treating this as
237  // the units part of the processor number.
238  sort
239  (
240  pointIndices,
241  UList<scalar>::less(rotatedPoints.component(vector::X))
242  );
243 
244  const scalar summedWeights = sum(weights);
245  assignToProcessorGroup
246  (
247  processorGroups,
248  n_.x(),
249  pointIndices,
250  weights,
251  summedWeights
252  );
253 
254  forAll(points, i)
255  {
256  finalDecomp[pointIndices[i]] = processorGroups[i];
257  }
258 
259 
260  // now do the same thing in the Y direction. These processor group
261  // numbers add multiples of nX to the proc. number (columns)
262  sort
263  (
264  pointIndices,
265  UList<scalar>::less(rotatedPoints.component(vector::Y))
266  );
267 
268  assignToProcessorGroup
269  (
270  processorGroups,
271  n_.y(),
272  pointIndices,
273  weights,
274  summedWeights
275  );
276 
277  forAll(points, i)
278  {
279  finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
280  }
281 
282 
283  // finally in the Z direction. Now we add multiples of nX*nY to give
284  // layers
285  sort
286  (
287  pointIndices,
288  UList<scalar>::less(rotatedPoints.component(vector::Z))
289  );
290 
291  assignToProcessorGroup
292  (
293  processorGroups,
294  n_.z(),
295  pointIndices,
296  weights,
297  summedWeights
298  );
299 
300  forAll(points, i)
301  {
302  finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
303  }
304 
305  return finalDecomp;
306 }
307 
308 
309 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
310 
311 Foam::simpleGeomDecomp::simpleGeomDecomp(const dictionary& decompDict)
312 :
313  geomDecomp(typeName, decompDict)
314 {}
315 
316 
317 Foam::simpleGeomDecomp::simpleGeomDecomp
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 slaves
353  for (const int slave : Pstream::subProcs())
354  {
355  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
356  pointField nbrPoints(fromSlave);
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 slave : Pstream::subProcs())
371  {
373  toSlave << SubField<label>
374  (
375  finalDecomp,
376  globalNumbers.localSize(slave),
377  globalNumbers.offset(slave)
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 slaves
438  for (const int slave : Pstream::subProcs())
439  {
440  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
441  pointField nbrPoints(fromSlave);
442  scalarField nbrWeights(fromSlave);
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 slave : Pstream::subProcs())
463  {
465  toSlave << SubField<label>
466  (
467  finalDecomp,
468  globalNumbers.localSize(slave),
469  globalNumbers.offset(slave)
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:71
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::UPstream::masterNo
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:452
Foam::Vector::Z
Definition: Vector.H:81
Foam::Vector::Y
Definition: Vector.H:81
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
SubField.H
Foam::less
static bool less(const vector &x, const vector &y)
To compare normals.
Definition: meshRefinementRefine.C:57
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
Foam::simpleGeomDecomp::decompose
virtual labelList decompose(const pointField &points) const
Decompose with uniform weights.
Definition: simpleGeomDecomp.C:330
globalIndex.H
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
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.
Definition: Field.H:64
SortableList.H
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
Foam::geomDecomp
Base for geometrical domain decomposition methods.
Definition: geomDecomp.H:71
Foam::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:516
Foam::Field< vector >
Foam::Vector::X
Definition: Vector.H:81
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:254
Foam::UPstream::commsTypes::scheduled
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
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::List< label >
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::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
simpleGeomDecomp.H