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