hierarchGeomDecomp.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 "hierarchGeomDecomp.H"
28 #include "PstreamReduceOps.H"
29 #include "SortableList.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  defineTypeNameAndDebug(hierarchGeomDecomp, 0);
36 
38  (
39  decompositionMethod,
40  hierarchGeomDecomp,
41  decomposer
42  );
43 
45  (
46  decompositionMethod,
47  hierarchGeomDecomp,
48  distributor
49  );
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 void Foam::hierarchGeomDecomp::setDecompOrder()
56 {
57  const word order(geomDecomDict_.lookup("order"));
58 
59  if (order.size() != 3)
60  {
62  (
63  decompositionDict_
64  ) << "number of characters in order (" << order << ") != 3"
65  << exit(FatalIOError);
66  }
67 
68  for (label i = 0; i < 3; ++i)
69  {
70  if (order[i] == 'x')
71  {
72  decompOrder_[i] = 0;
73  }
74  else if (order[i] == 'y')
75  {
76  decompOrder_[i] = 1;
77  }
78  else if (order[i] == 'z')
79  {
80  decompOrder_[i] = 2;
81  }
82  else
83  {
85  (
86  decompositionDict_
87  ) << "Illegal decomposition order " << order << endl
88  << "It should only contain x, y or z" << exit(FatalError);
89  }
90  }
91 }
92 
93 
94 Foam::label Foam::hierarchGeomDecomp::findLower
95 (
96  const List<scalar>& l,
97  const scalar t,
98  const label initLow,
99  const label initHigh
100 )
101 {
102  if (initHigh <= initLow)
103  {
104  return initLow;
105  }
106 
107  label low = initLow;
108  label high = initHigh;
109 
110  while ((high - low) > 1)
111  {
112  label mid = (low + high)/2;
113 
114  if (l[mid] < t)
115  {
116  low = mid;
117  }
118  else
119  {
120  high = mid;
121  }
122  }
123 
124  // high and low can still differ by one. Choose best.
125 
126  label tIndex = -1;
127 
128  if (l[high-1] < t)
129  {
130  tIndex = high;
131  }
132  else
133  {
134  tIndex = low;
135  }
136 
137  return tIndex;
138 }
139 
140 
141 // Create a mapping between the index and the weighted size.
142 // For convenience, sortedWeightedSize is one size bigger than current. This
143 // avoids extra tests.
144 void Foam::hierarchGeomDecomp::calculateSortedWeightedSizes
145 (
146  const labelList& current,
147  const labelList& indices,
148  const scalarField& weights,
149  const label globalCurrentSize,
150 
151  scalarField& sortedWeightedSizes
152 )
153 {
154  // Evaluate cumulative weights.
155  sortedWeightedSizes[0] = 0;
156  forAll(current, i)
157  {
158  label pointi = current[indices[i]];
159  sortedWeightedSizes[i + 1] = sortedWeightedSizes[i] + weights[pointi];
160  }
161  // Non-dimensionalise and multiply by size.
162  scalar globalCurrentLength = returnReduce
163  (
164  sortedWeightedSizes[current.size()],
165  sumOp<scalar>()
166  );
167  // Normalise weights by global sum of weights and multiply through
168  // by global size.
169  sortedWeightedSizes *= (globalCurrentSize/globalCurrentLength);
170 }
171 
172 
173 // Find position in values so between minIndex and this position there
174 // are wantedSize elements.
175 void Foam::hierarchGeomDecomp::findBinary
176 (
177  const label sizeTol,
178  const List<scalar>& values,
179  const label minIndex, // index of previous value
180  const scalar minValue, // value at minIndex
181  const scalar maxValue, // global max of values
182  const scalar wantedSize, // wanted size
183 
184  label& mid, // index where size of bin is
185  // wantedSize (to within sizeTol)
186  scalar& midValue // value at mid
187 )
188 {
189  label low = minIndex;
190  scalar lowValue = minValue;
191 
192  scalar highValue = maxValue;
193  // (one beyond) index of highValue
194  label high = values.size();
195 
196  // Safeguards to avoid infinite loop.
197  scalar midValuePrev = vGreat;
198 
199  while (true)
200  {
201  label size = returnReduce(mid-minIndex, sumOp<label>());
202 
203  if (debug)
204  {
205  Pout<< " low:" << low << " lowValue:" << lowValue
206  << " high:" << high << " highValue:" << highValue
207  << " mid:" << mid << " midValue:" << midValue << endl
208  << " globalSize:" << size << " wantedSize:" << wantedSize
209  << " sizeTol:" << sizeTol << endl;
210  }
211 
212  if (wantedSize < size - sizeTol)
213  {
214  high = mid;
215  highValue = midValue;
216  }
217  else if (wantedSize > size + sizeTol)
218  {
219  low = mid;
220  lowValue = midValue;
221  }
222  else
223  {
224  break;
225  }
226 
227  // Update mid, midValue
228  midValue = 0.5*(lowValue+highValue);
229  mid = findLower(values, midValue, low, high);
230 
231  // Safeguard if same as previous.
232  bool hasNotChanged = (mag(midValue-midValuePrev) < small);
233 
234  if (returnReduce(hasNotChanged, andOp<bool>()))
235  {
237  << "unable to find desired decomposition split, making do!"
238  << endl;
239  break;
240  }
241 
242  midValuePrev = midValue;
243  }
244 }
245 
246 
247 // Find position in values so between minIndex and this position there
248 // are wantedSize elements.
249 void Foam::hierarchGeomDecomp::findBinary
250 (
251  const label sizeTol,
252  const List<scalar>& sortedWeightedSizes,
253  const List<scalar>& values,
254  const label minIndex, // index of previous value
255  const scalar minValue, // value at minIndex
256  const scalar maxValue, // global max of values
257  const scalar wantedSize, // wanted size
258 
259  label& mid, // index where size of bin is
260  // wantedSize (to within sizeTol)
261  scalar& midValue // value at mid
262 )
263 {
264  label low = minIndex;
265  scalar lowValue = minValue;
266 
267  scalar highValue = maxValue;
268  // (one beyond) index of highValue
269  label high = values.size();
270 
271  // Safeguards to avoid infinite loop.
272  scalar midValuePrev = vGreat;
273 
274  while (true)
275  {
276  scalar weightedSize = returnReduce
277  (
278  sortedWeightedSizes[mid] - sortedWeightedSizes[minIndex],
279  sumOp<scalar>()
280  );
281 
282  if (debug)
283  {
284  Pout<< " low:" << low << " lowValue:" << lowValue
285  << " high:" << high << " highValue:" << highValue
286  << " mid:" << mid << " midValue:" << midValue << endl
287  << " globalSize:" << weightedSize
288  << " wantedSize:" << wantedSize
289  << " sizeTol:" << sizeTol << endl;
290  }
291 
292  if (wantedSize < weightedSize - sizeTol)
293  {
294  high = mid;
295  highValue = midValue;
296  }
297  else if (wantedSize > weightedSize + sizeTol)
298  {
299  low = mid;
300  lowValue = midValue;
301  }
302  else
303  {
304  break;
305  }
306 
307  // Update mid, midValue
308  midValue = 0.5*(lowValue+highValue);
309  mid = findLower(values, midValue, low, high);
310 
311  // Safeguard if same as previous.
312  bool hasNotChanged = (mag(midValue-midValuePrev) < small);
313 
314  if (returnReduce(hasNotChanged, andOp<bool>()))
315  {
317  << "unable to find desired decomposition split, making do!"
318  << endl;
319  break;
320  }
321 
322  midValuePrev = midValue;
323  }
324 }
325 
326 
327 // Sort points into bins according to one component. Recurses to next component.
328 void Foam::hierarchGeomDecomp::sortComponent
329 (
330  const label sizeTol,
331  const pointField& points,
332  const labelList& current, // slice of points to decompose
333  const direction componentIndex, // index in decompOrder_
334  const label mult, // multiplication factor for finalDecomp
335  labelList& finalDecomp
336 )
337 {
338  // Current component
339  label compI = decompOrder_[componentIndex];
340 
341  if (debug)
342  {
343  Pout<< "sortComponent : Sorting slice of size " << current.size()
344  << " in component " << compI << endl;
345  }
346 
347  // Storage for sorted component compI
348  SortableList<scalar> sortedCoord(current.size());
349 
350  forAll(current, i)
351  {
352  label pointi = current[i];
353 
354  sortedCoord[i] = points[pointi][compI];
355  }
356  sortedCoord.sort();
357 
358  label globalCurrentSize = returnReduce(current.size(), sumOp<label>());
359 
360  scalar minCoord = returnReduce
361  (
362  (
363  sortedCoord.size()
364  ? sortedCoord[0]
365  : great
366  ),
367  minOp<scalar>()
368  );
369 
370  scalar maxCoord = returnReduce
371  (
372  (
373  sortedCoord.size()
374  ? sortedCoord.last()
375  : -great
376  ),
377  maxOp<scalar>()
378  );
379 
380  if (debug)
381  {
382  Pout<< "sortComponent : minCoord:" << minCoord
383  << " maxCoord:" << maxCoord << endl;
384  }
385 
386  // starting index (in sortedCoord) of bin (= local)
387  label leftIndex = 0;
388  // starting value of bin (= global since coordinate)
389  scalar leftCoord = minCoord;
390 
391  // Sort bins of size n
392  for (label bin = 0; bin < n_[compI]; bin++)
393  {
394  // Now we need to determine the size of the bin (dx). This is
395  // determined by the 'pivot' values - everything to the left of this
396  // value goes in the current bin, everything to the right into the next
397  // bins.
398 
399  // Local number of elements
400  label localSize = -1; // offset from leftOffset
401 
402  // Value at right of bin (leftIndex+localSize-1)
403  scalar rightCoord = -great;
404 
405  if (bin == n_[compI]-1)
406  {
407  // Last bin. Copy all.
408  localSize = current.size()-leftIndex;
409  rightCoord = maxCoord; // note: not used anymore
410  }
411  else if (Pstream::nProcs() == 1)
412  {
413  // No need for binary searching of bin size
414  localSize = label(current.size()/n_[compI]);
415  rightCoord = sortedCoord[leftIndex+localSize];
416  }
417  else
418  {
419  // For the current bin (starting at leftCoord) we want a rightCoord
420  // such that the sum of all sizes are globalCurrentSize/n_[compI].
421  // We have to iterate to obtain this.
422 
423  label rightIndex = current.size();
424  rightCoord = maxCoord;
425 
426  // Calculate rightIndex/rightCoord to have wanted size
427  findBinary
428  (
429  sizeTol,
430  sortedCoord,
431  leftIndex,
432  leftCoord,
433  maxCoord,
434  globalCurrentSize/n_[compI], // wanted size
435 
436  rightIndex,
437  rightCoord
438  );
439  localSize = rightIndex - leftIndex;
440  }
441 
442  if (debug)
443  {
444  Pout<< "For component " << compI << ", bin " << bin
445  << " copying" << endl
446  << "from " << leftCoord << " at local index "
447  << leftIndex << endl
448  << "to " << rightCoord << " localSize:"
449  << localSize << endl
450  << endl;
451  }
452 
453 
454  // Copy localSize elements starting from leftIndex.
455  labelList slice(localSize);
456 
457  forAll(slice, i)
458  {
459  label pointi = current[sortedCoord.indices()[leftIndex+i]];
460 
461  // Mark point into correct bin
462  finalDecomp[pointi] += bin*mult;
463 
464  // And collect for next sorting action
465  slice[i] = pointi;
466  }
467 
468  // Sort slice in next component
469  if (componentIndex < 2)
470  {
471  string oldPrefix;
472  if (debug)
473  {
474  oldPrefix = Pout.prefix();
475  Pout.prefix() = " " + oldPrefix;
476  }
477 
478  sortComponent
479  (
480  sizeTol,
481  points,
482  slice,
483  componentIndex+1,
484  mult*n_[compI], // Multiplier to apply to decomposition.
485  finalDecomp
486  );
487 
488  if (debug)
489  {
490  Pout.prefix() = oldPrefix;
491  }
492  }
493 
494  // Step to next bin.
495  leftIndex += localSize;
496  leftCoord = rightCoord;
497  }
498 }
499 
500 
501 // Sort points into bins according to one component. Recurses to next component.
502 void Foam::hierarchGeomDecomp::sortComponent
503 (
504  const label sizeTol,
505  const scalarField& weights,
506  const pointField& points,
507  const labelList& current, // slice of points to decompose
508  const direction componentIndex, // index in decompOrder_
509  const label mult, // multiplication factor for finalDecomp
510  labelList& finalDecomp
511 )
512 {
513  // Current component
514  label compI = decompOrder_[componentIndex];
515 
516  if (debug)
517  {
518  Pout<< "sortComponent : Sorting slice of size " << current.size()
519  << " in component " << compI << endl;
520  }
521 
522  // Storage for sorted component compI
523  SortableList<scalar> sortedCoord(current.size());
524 
525  forAll(current, i)
526  {
527  label pointi = current[i];
528 
529  sortedCoord[i] = points[pointi][compI];
530  }
531  sortedCoord.sort();
532 
533  label globalCurrentSize = returnReduce(current.size(), sumOp<label>());
534 
535  // Now evaluate local cumulative weights, based on the sorting.
536  // Make one bigger than the nodes.
537  scalarField sortedWeightedSizes(current.size()+1, 0);
538  calculateSortedWeightedSizes
539  (
540  current,
541  sortedCoord.indices(),
542  weights,
543  globalCurrentSize,
544  sortedWeightedSizes
545  );
546 
547  scalar minCoord = returnReduce
548  (
549  (
550  sortedCoord.size()
551  ? sortedCoord[0]
552  : great
553  ),
554  minOp<scalar>()
555  );
556 
557  scalar maxCoord = returnReduce
558  (
559  (
560  sortedCoord.size()
561  ? sortedCoord.last()
562  : -great
563  ),
564  maxOp<scalar>()
565  );
566 
567  if (debug)
568  {
569  Pout<< "sortComponent : minCoord:" << minCoord
570  << " maxCoord:" << maxCoord << endl;
571  }
572 
573  // starting index (in sortedCoord) of bin (= local)
574  label leftIndex = 0;
575  // starting value of bin (= global since coordinate)
576  scalar leftCoord = minCoord;
577 
578  // Sort bins of size n
579  for (label bin = 0; bin < n_[compI]; bin++)
580  {
581  // Now we need to determine the size of the bin (dx). This is
582  // determined by the 'pivot' values - everything to the left of this
583  // value goes in the current bin, everything to the right into the next
584  // bins.
585 
586  // Local number of elements
587  label localSize = -1; // offset from leftOffset
588 
589  // Value at right of bin (leftIndex+localSize-1)
590  scalar rightCoord = -great;
591 
592  if (bin == n_[compI]-1)
593  {
594  // Last bin. Copy all.
595  localSize = current.size()-leftIndex;
596  rightCoord = maxCoord; // note: not used anymore
597  }
598  else
599  {
600  // For the current bin (starting at leftCoord) we want a rightCoord
601  // such that the sum of all weighted sizes are
602  // globalCurrentLength/n_[compI].
603  // We have to iterate to obtain this.
604 
605  label rightIndex = current.size();
606  rightCoord = maxCoord;
607 
608  // Calculate rightIndex/rightCoord to have wanted size
609  findBinary
610  (
611  sizeTol,
612  sortedWeightedSizes,
613  sortedCoord,
614  leftIndex,
615  leftCoord,
616  maxCoord,
617  globalCurrentSize/n_[compI], // wanted size
618 
619  rightIndex,
620  rightCoord
621  );
622  localSize = rightIndex - leftIndex;
623  }
624 
625  if (debug)
626  {
627  Pout<< "For component " << compI << ", bin " << bin
628  << " copying" << endl
629  << "from " << leftCoord << " at local index "
630  << leftIndex << endl
631  << "to " << rightCoord << " localSize:"
632  << localSize << endl
633  << endl;
634  }
635 
636 
637  // Copy localSize elements starting from leftIndex.
638  labelList slice(localSize);
639 
640  forAll(slice, i)
641  {
642  label pointi = current[sortedCoord.indices()[leftIndex+i]];
643 
644  // Mark point into correct bin
645  finalDecomp[pointi] += bin*mult;
646 
647  // And collect for next sorting action
648  slice[i] = pointi;
649  }
650 
651  // Sort slice in next component
652  if (componentIndex < 2)
653  {
654  string oldPrefix;
655  if (debug)
656  {
657  oldPrefix = Pout.prefix();
658  Pout.prefix() = " " + oldPrefix;
659  }
660 
661  sortComponent
662  (
663  sizeTol,
664  weights,
665  points,
666  slice,
667  componentIndex+1,
668  mult*n_[compI], // Multiplier to apply to decomposition.
669  finalDecomp
670  );
671 
672  if (debug)
673  {
674  Pout.prefix() = oldPrefix;
675  }
676  }
677 
678  // Step to next bin.
679  leftIndex += localSize;
680  leftCoord = rightCoord;
681  }
682 }
683 
684 
686 (
687  const pointField& points
688 )
689 {
690  // construct a list for the final result
691  labelList finalDecomp(points.size(), 0);
692 
693  // Start off with every point sorted onto itself.
694  labelList slice(points.size());
695  forAll(slice, i)
696  {
697  slice[i] = i;
698  }
699 
700  pointField rotatedPoints(rotDelta_ & points);
701 
702  // Calculate tolerance of cell distribution. For large cases finding
703  // distribution to the cell exact would cause too many iterations so allow
704  // some slack.
705  label allSize = points.size();
706  reduce(allSize, sumOp<label>());
707 
708  const label sizeTol = max(1, label(1e-3*allSize/nProcessors_));
709 
710  // Sort recursive
711  sortComponent
712  (
713  sizeTol,
714  rotatedPoints,
715  slice,
716  0, // Sort first component in decompOrder.
717  1, // Offset for different x bins.
718  finalDecomp
719  );
720 
721  return finalDecomp;
722 }
723 
724 
726 (
727  const pointField& points,
728  const scalarField& weights
729 )
730 {
731  // construct a list for the final result
732  labelList finalDecomp(points.size(), 0);
733 
734  // Start off with every point sorted onto itself.
735  labelList slice(points.size());
736  forAll(slice, i)
737  {
738  slice[i] = i;
739  }
740 
741  pointField rotatedPoints(rotDelta_ & points);
742 
743  // Calculate tolerance of cell distribution. For large cases finding
744  // distribution to the cell exact would cause too many iterations so allow
745  // some slack.
746  label allSize = points.size();
747  reduce(allSize, sumOp<label>());
748 
749  const label sizeTol = max(1, label(1e-3*allSize/nProcessors_));
750 
751  // Sort recursive
752  sortComponent
753  (
754  sizeTol,
755  weights,
756  rotatedPoints,
757  slice,
758  0, // Sort first component in decompOrder.
759  1, // Offset for different x bins.
760  finalDecomp
761  );
762 
763  return finalDecomp;
764 }
765 
766 
767 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
768 
770 (
771  const dictionary& decompositionDict
772 )
773 :
774  geomDecomp(decompositionDict, typeName),
775  decompOrder_()
776 {
777  setDecompOrder();
778 }
779 
780 
781 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Inter-processor communication reduction functions.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
scalar minValue
uint8_t direction
Definition: direction.H:45
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Holds information (coordinate and normal) regarding nearest wall point.
Macros for easy insertion into run-time selection tables.
const string & prefix() const
Return the prefix of the stream.
scalar maxValue
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Geometrical domain decomposition.
Definition: geomDecomp.H:47
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
hierarchGeomDecomp(const dictionary &decompositionDict)
Construct given the decomposition dictionary.
List< label > labelList
A List of labels.
Definition: labelList.H:56
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
virtual labelList decompose(const pointField &, const scalarField &weights)
Return for every coordinate the wanted processor number.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:318
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findLower(const ListType &, typename ListType::const_reference, const label stary, const BinaryOp &bop)
Find last element < given value in sorted list and return index,.
dimensioned< scalar > mag(const dimensioned< Type > &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
Namespace for OpenFOAM.
IOerror FatalIOError