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