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