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