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