autoDensity.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) 2012-2016 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 "autoDensity.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 defineTypeNameAndDebug(autoDensity, 0);
38 (
39  initialPointsMethod,
40  autoDensity,
41  dictionary
42 );
43 
44 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 
46 void Foam::autoDensity::writeOBJ
47 (
48  const treeBoundBox& bb,
49  fileName name
50 ) const
51 {
52  OFstream str(time().path()/name + ".obj");
53 
54  Pout<< "Writing " << str.name() << endl;
55 
56  pointField bbPoints(bb.points());
57 
58  forAll(bbPoints, i)
59  {
60  meshTools::writeOBJ(str, bbPoints[i]);
61  }
62 
64  {
65  const edge& e = treeBoundBox::edges[i];
66 
67  str << "l " << e[0] + 1 << ' ' << e[1] + 1 << nl;
68  }
69 }
70 
71 
72 bool Foam::autoDensity::combinedOverlaps(const treeBoundBox& box) const
73 {
74  if (Pstream::parRun())
75  {
76  return
77  decomposition().overlapsThisProcessor(box)
78  || geometryToConformTo().overlaps(box);
79  }
80 
81  return geometryToConformTo().overlaps(box);
82 }
83 
84 
85 bool Foam::autoDensity::combinedInside(const point& p) const
86 {
87  if (Pstream::parRun())
88  {
89  return
90  decomposition().positionOnThisProcessor(p)
91  && geometryToConformTo().inside(p);
92  }
93 
94  return geometryToConformTo().inside(p);
95 }
96 
97 
98 Foam::Field<bool> Foam::autoDensity::combinedWellInside
99 (
100  const pointField& pts,
101  const scalarField& sizes
102 ) const
103 {
104  if (!Pstream::parRun())
105  {
106  return geometryToConformTo().wellInside
107  (
108  pts,
109  minimumSurfaceDistanceCoeffSqr_*sqr(sizes)
110  );
111  }
112 
113  Field<bool> inside(pts.size(), true);
114 
115  // Perform AND operation between testing the surfaces and the previous
116  // field, i.e the parallel result, or in serial, with true.
117 
118  Field<bool> insideA
119  (
120  geometryToConformTo().wellInside
121  (
122  pts,
123  minimumSurfaceDistanceCoeffSqr_*sqr(sizes)
124  )
125  );
126 
127  Field<bool> insideB
128  (
129  decomposition().positionOnThisProcessor(pts)
130  );
131 
132  // inside = insideA && insideB;
133 
134  // Pout<< insideA << nl << insideB << endl;
135 
136  forAll(inside, i)
137  {
138  // if (inside[i] != (insideA[i] && insideB[i]))
139  // {
140  // Pout<< i << " not equal " << " "
141  // << pts[i] << " " << sizes[i] << " "
142  // << insideA[i] << " "
143  // << insideB[i] << " "
144  // << inside[i]
145  // << endl;
146  // }
147 
148  inside[i] = (insideA[i] && insideB[i]);
149  }
150 
151  return inside;
152 }
153 
154 
155 bool Foam::autoDensity::combinedWellInside
156 (
157  const point& p,
158  scalar size
159 ) const
160 {
161  bool inside = true;
162 
163  if (Pstream::parRun())
164  {
165  inside = decomposition().positionOnThisProcessor(p);
166  }
167 
168  // Perform AND operation between testing the surfaces and the previous
169  // result, i.e the parallel result, or in serial, with true.
170  inside =
171  inside
172  && geometryToConformTo().wellInside
173  (
174  p,
175  minimumSurfaceDistanceCoeffSqr_*sqr(size)
176  );
177 
178  return inside;
179 }
180 
181 
182 Foam::label Foam::autoDensity::recurseAndFill
183 (
184  DynamicList<Vb::Point>& initialPoints,
185  const treeBoundBox& bb,
186  label levelLimit,
187  word recursionName
188 ) const
189 {
190  label treeDepth = 0;
191 
192  for (direction i = 0; i < 8; i++)
193  {
194  treeBoundBox subBB = bb.subBbox(i);
195 
196  word newName = recursionName + "_" + Foam::name(i);
197 
198  conformalVoronoiMesh::timeCheck(time(), newName, debug);
199 
200  if (combinedOverlaps(subBB))
201  {
202  if (levelLimit > 0)
203  {
204  treeDepth =
205  max
206  (
207  treeDepth,
208  recurseAndFill
209  (
210  initialPoints,
211  subBB,
212  levelLimit - 1,
213  newName
214  )
215  );
216  }
217  else
218  {
219  if (debug)
220  {
221  writeOBJ
222  (
223  subBB,
224  word(newName + "_overlap")
225  );
226 
227  Pout<< newName + "_overlap " << subBB << endl;
228  }
229 
230  if (!fillBox(initialPoints, subBB, true))
231  {
232  treeDepth =
233  max
234  (
235  treeDepth,
236  recurseAndFill
237  (
238  initialPoints,
239  subBB,
240  levelLimit - 1,
241  newName
242  )
243  );
244  }
245  }
246  }
247  else if (combinedInside(subBB.midpoint()))
248  {
249  if (debug)
250  {
251  writeOBJ
252  (
253  subBB,
254  newName + "_inside"
255  );
256 
257  Pout<< newName + "_inside " << subBB << endl;
258  }
259 
260  if (!fillBox(initialPoints, subBB, false))
261  {
262  treeDepth =
263  max
264  (
265  treeDepth,
266  recurseAndFill
267  (
268  initialPoints,
269  subBB,
270  levelLimit - 1,
271  newName
272  )
273  );
274  }
275  }
276  else
277  {
278  if (debug)
279  {
280  writeOBJ
281  (
282  subBB,
283  newName + "_outside"
284  );
285  }
286  }
287  }
288 
289  return treeDepth + 1;
290 }
291 
292 
293 bool Foam::autoDensity::fillBox
294 (
295  DynamicList<Vb::Point>& initialPoints,
296  const treeBoundBox& bb,
297  bool overlapping
298 ) const
299 {
300  const conformationSurfaces& geometry = geometryToConformTo();
301 
302  label initialSize = initialPoints.size();
303 
304  scalar maxCellSize = -GREAT;
305 
306  scalar minCellSize = GREAT;
307 
308  scalar maxDensity = 1/pow3(minCellSize);
309 
310  scalar volumeAdded = 0.0;
311 
312  const point& min = bb.min();
313 
314  vector span = bb.span();
315 
316  scalar totalVolume = bb.volume();
317 
318  label trialPoints = 0;
319 
320  bool wellInside = false;
321 
322  if (!overlapping)
323  {
324  // Check the nearest point on the surface to the box, if it is far
325  // enough away, then the surface sampling of the box can be skipped.
326  // Checking if the nearest piece of surface is at least 1.5*bb.span away
327  // from the bb.midpoint.
328 
329  pointIndexHit surfHit;
330  label hitSurface;
331 
332  geometry.findSurfaceNearest
333  (
334  bb.midpoint(),
335  2.25*magSqr(span),
336  surfHit,
337  hitSurface
338  );
339 
340  if (!surfHit.hit())
341  {
342  if (debug)
343  {
344  Pout<< "box wellInside, no need to sample surface." << endl;
345  }
346 
347  wellInside = true;
348  }
349  }
350 
351  if (!overlapping && !wellInside)
352  {
353  // If this is an inside box then it is possible to fill points very
354  // close to the boundary, to prevent this, check the corners and sides
355  // of the box so ensure that they are "wellInside". If not, set as an
356  // overlapping box.
357 
358  pointField corners(bb.points());
359 
360  scalarField cornerSizes = cellShapeControls().cellSize(corners);
361 
362  Field<bool> insideCorners = combinedWellInside(corners, cornerSizes);
363 
364  // Pout<< corners << nl << cornerSizes << nl << insideCorners << endl;
365 
366  forAll(insideCorners, i)
367  {
368  // Use the sizes to improve the min/max cell size estimate
369  scalar s = cornerSizes[i];
370 
371  if (s > maxCellSize)
372  {
373  maxCellSize = s;
374  }
375 
376  if (s < minCellSize)
377  {
378  minCellSize = max(s, minCellSizeLimit_);
379  }
380 
381  if (maxCellSize/minCellSize > maxSizeRatio_)
382  {
383  if (debug)
384  {
385  Pout<< "Abort fill at corner sample stage,"
386  << " minCellSize " << minCellSize
387  << " maxCellSize " << maxCellSize
388  << " maxSizeRatio " << maxCellSize/minCellSize
389  << endl;
390  }
391 
392  return false;
393  }
394 
395  if (!insideCorners[i])
396  {
397  // If one or more corners is not "wellInside", then treat this
398  // as an overlapping box.
399 
400  if (debug)
401  {
402  Pout<< "Inside box found to have some non-wellInside "
403  << "corners, using overlapping fill."
404  << endl;
405  }
406 
407  overlapping = true;
408 
409  break;
410  }
411  }
412 
413  if (!overlapping)
414  {
415  vector delta = span/(surfRes_ - 1);
416 
417  label nLine = 6*(surfRes_ - 2);
418 
419  pointField linePoints(nLine, Zero);
420 
421  scalarField lineSizes(nLine, 0.0);
422 
423  for (label i = 0; i < surfRes_; i++)
424  {
425  label lPI = 0;
426 
427  for (label j = 1; j < surfRes_ - 1 ; j++)
428  {
429  linePoints[lPI++] =
430  min
431  + vector(0, delta.y()*i, delta.z()*j);
432 
433  linePoints[lPI++] =
434  min
435  + vector
436  (
437  delta.x()*(surfRes_ - 1),
438  delta.y()*i,
439  delta.z()*j
440  );
441 
442  linePoints[lPI++] =
443  min
444  + vector(delta.x()*j, 0, delta.z()*i);
445 
446  linePoints[lPI++] =
447  min
448  + vector
449  (
450  delta.x()*j,
451  delta.y()*(surfRes_ - 1),
452  delta.z()*i
453  );
454 
455  linePoints[lPI++] =
456  min
457  + vector(delta.x()*i, delta.y()*j, 0);
458 
459  linePoints[lPI++] =
460  min
461  + vector
462  (
463  delta.x()*i,
464  delta.y()*j,
465  delta.z()*(surfRes_ - 1)
466  );
467  }
468 
469  lineSizes = cellShapeControls().cellSize(linePoints);
470 
471  Field<bool> insideLines = combinedWellInside
472  (
473  linePoints,
474  lineSizes
475  );
476 
477  forAll(insideLines, i)
478  {
479  // Use the sizes to improve the min/max cell size estimate
480  scalar s = lineSizes[i];
481 
482  if (s > maxCellSize)
483  {
484  maxCellSize = s;
485  }
486 
487  if (s < minCellSize)
488  {
489  minCellSize = max(s, minCellSizeLimit_);
490  }
491 
492  if (maxCellSize/minCellSize > maxSizeRatio_)
493  {
494  if (debug)
495  {
496  Pout<< "Abort fill at surface sample stage, "
497  << " minCellSize " << minCellSize
498  << " maxCellSize " << maxCellSize
499  << " maxSizeRatio " << maxCellSize/minCellSize
500  << endl;
501  }
502 
503  return false;
504  }
505 
506  if (!insideLines[i])
507  {
508  // If one or more surface points is not "wellInside",
509  // then treat this as an overlapping box.
510  overlapping = true;
511 
512  if (debug)
513  {
514  Pout<< "Inside box found to have some non-"
515  << "wellInside surface points, using "
516  << "overlapping fill."
517  << endl;
518  }
519 
520  break;
521  }
522  }
523  }
524  }
525  }
526 
527  if (overlapping)
528  {
529  // Sample the box to find an estimate of the min size, and a volume
530  // estimate when overlapping == true.
531 
532  pointField samplePoints
533  (
534  volRes_*volRes_*volRes_,
535  Zero
536  );
537 
538  vector delta = span/volRes_;
539 
540  label pI = 0;
541 
542  for (label i = 0; i < volRes_; i++)
543  {
544  for (label j = 0; j < volRes_; j++)
545  {
546  for (label k = 0; k < volRes_; k++)
547  {
548  // Perturb the points to avoid creating degenerate positions
549  // in the Delaunay tessellation.
550 
551  samplePoints[pI++] =
552  min
553  + vector
554  (
555  delta.x()
556  *(i + 0.5 + 0.1*(rndGen().scalar01() - 0.5)),
557  delta.y()
558  *(j + 0.5 + 0.1*(rndGen().scalar01() - 0.5)),
559  delta.z()
560  *(k + 0.5 + 0.1*(rndGen().scalar01() - 0.5))
561  );
562  }
563  }
564  }
565 
566  // Randomise the order of the points to (potentially) improve the speed
567  // of assessing the density ratio, and prevent a box being filled from a
568  // corner when only some these points are required.
569  shuffle(samplePoints);
570 
571  scalarField sampleSizes = cellShapeControls().cellSize(samplePoints);
572 
573  Field<bool> insidePoints = combinedWellInside
574  (
575  samplePoints,
576  sampleSizes
577  );
578 
579  label nInside = 0;
580 
581  forAll(insidePoints, i)
582  {
583  if (insidePoints[i])
584  {
585  nInside++;
586 
587  scalar s = sampleSizes[i];
588 
589  if (s > maxCellSize)
590  {
591  maxCellSize = s;
592  }
593 
594  if (s < minCellSize)
595  {
596  minCellSize = max(s, minCellSizeLimit_);
597  }
598 
599  if (maxCellSize/minCellSize > maxSizeRatio_)
600  {
601  if (debug)
602  {
603  Pout<< "Abort fill at sample stage,"
604  << " minCellSize " << minCellSize
605  << " maxCellSize " << maxCellSize
606  << " maxSizeRatio " << maxCellSize/minCellSize
607  << endl;
608  }
609 
610  return false;
611  }
612  }
613  }
614 
615  if (nInside == 0)
616  {
617  if (debug)
618  {
619  Pout<< "No sample points found inside box" << endl;
620  }
621 
622  return true;
623  }
624 
625  if (debug)
626  {
627  Pout<< scalar(nInside)/scalar(samplePoints.size())
628  << " full overlapping box" << endl;
629  }
630 
631  totalVolume *= scalar(nInside)/scalar(samplePoints.size());
632 
633  if (debug)
634  {
635  Pout<< "Total volume to fill = " << totalVolume << endl;
636  }
637 
638  // Using the sampledPoints as the first test locations as they are
639  // randomly shuffled, but unfiormly sampling space and have wellInside
640  // and size data already
641 
642  maxDensity = 1/pow3(max(minCellSize, SMALL));
643 
644  forAll(insidePoints, i)
645  {
646  if (insidePoints[i])
647  {
648  trialPoints++;
649 
650  const point& p = samplePoints[i];
651 
652  scalar localSize = sampleSizes[i];
653 
654  scalar localDensity = 1/pow3(localSize);
655 
656  // No need to look at max/min cell size here, already handled
657  // by sampling
658 
659  // Accept possible placements proportional to the relative
660  // local density
661 
662  // TODO - is there a lot of cost in the 1/density calc? Could
663  // assess on
664  // (1/maxDensity)/(1/localDensity) = minVolume/localVolume
665  if (localDensity/maxDensity > rndGen().scalar01())
666  {
667  scalar localVolume = 1/localDensity;
668 
669  if (volumeAdded + localVolume > totalVolume)
670  {
671  // Add the final box with a probability of to the ratio
672  // of the remaining volume to the volume to be added,
673  // i.e. insert a box of volume 0.5 into a remaining
674  // volume of 0.1 20% of the time.
675  scalar addProbability =
676  (totalVolume - volumeAdded)/localVolume;
677 
678  scalar r = rndGen().scalar01();
679 
680  if (debug)
681  {
682  Pout<< "totalVolume " << totalVolume << nl
683  << "volumeAdded " << volumeAdded << nl
684  << "localVolume " << localVolume << nl
685  << "addProbability " << addProbability << nl
686  << "random " << r
687  << endl;
688  }
689 
690  if (addProbability > r)
691  {
692  // Place this volume before finishing filling this
693  // box
694 
695  // Pout<< "Final volume probability break accept"
696  // << endl;
697 
698  initialPoints.append
699  (
700  Vb::Point(p.x(), p.y(), p.z())
701  );
702 
703  volumeAdded += localVolume;
704  }
705 
706  break;
707  }
708 
709  initialPoints.append(Vb::Point(p.x(), p.y(), p.z()));
710 
711  volumeAdded += localVolume;
712  }
713  }
714  }
715  }
716 
717  if (volumeAdded < totalVolume)
718  {
719  if (debug)
720  {
721  Pout<< "Adding random points, remaining volume "
722  << totalVolume - volumeAdded
723  << endl;
724  }
725 
726  maxDensity = 1/pow3(max(minCellSize, SMALL));
727 
728  while (true)
729  {
730  trialPoints++;
731 
732  point p = min + cmptMultiply(span, rndGen().vector01());
733 
734  scalar localSize = cellShapeControls().cellSize(p);
735 
736  bool insidePoint = false;
737 
738  if (!overlapping)
739  {
740  insidePoint = true;
741  }
742  else
743  {
744  // Determine if the point is "wellInside" the domain
745  insidePoint = combinedWellInside(p, localSize);
746  }
747 
748  if (insidePoint)
749  {
750  if (localSize > maxCellSize)
751  {
752  maxCellSize = localSize;
753  }
754 
755  if (localSize < minCellSize)
756  {
757  minCellSize = max(localSize, minCellSizeLimit_);
758 
759  localSize = minCellSize;
760 
761  // 1/(minimum cell size)^3, gives the maximum permissible
762  // point density
763  maxDensity = 1/pow3(max(minCellSize, SMALL));
764  }
765 
766  if (maxCellSize/minCellSize > maxSizeRatio_)
767  {
768  if (debug)
769  {
770  Pout<< "Abort fill at random fill stage,"
771  << " minCellSize " << minCellSize
772  << " maxCellSize " << maxCellSize
773  << " maxSizeRatio " << maxCellSize/minCellSize
774  << endl;
775  }
776 
777  // Discard any points already filled into this box by
778  // setting size of initialPoints back to its starting value
779  initialPoints.resize(initialSize);
780 
781  return false;
782  }
783 
784  scalar localDensity = 1/pow3(max(localSize, SMALL));
785 
786  // Accept possible placements proportional to the relative local
787  // density
788  if (localDensity/maxDensity > rndGen().scalar01())
789  {
790  scalar localVolume = 1/localDensity;
791 
792  if (volumeAdded + localVolume > totalVolume)
793  {
794  // Add the final box with a probability of to the ratio
795  // of the remaining volume to the volume to be added,
796  // i.e. insert a box of volume 0.5 into a remaining
797  // volume of 0.1 20% of the time.
798  scalar addProbability =
799  (totalVolume - volumeAdded)/localVolume;
800 
801  scalar r = rndGen().scalar01();
802 
803  if (debug)
804  {
805  Pout<< "totalVolume " << totalVolume << nl
806  << "volumeAdded " << volumeAdded << nl
807  << "localVolume " << localVolume << nl
808  << "addProbability " << addProbability << nl
809  << "random " << r
810  << endl;
811  }
812 
813  if (addProbability > r)
814  {
815  // Place this volume before finishing filling this
816  // box
817 
818  // Pout<< "Final volume probability break accept"
819  // << endl;
820 
821  initialPoints.append
822  (
823  Vb::Point(p.x(), p.y(), p.z())
824  );
825 
826  volumeAdded += localVolume;
827  }
828 
829  break;
830  }
831 
832  initialPoints.append(Vb::Point(p.x(), p.y(), p.z()));
833 
834  volumeAdded += localVolume;
835  }
836  }
837  }
838  }
839 
840  globalTrialPoints_ += trialPoints;
841 
842  if (debug)
843  {
844  Pout<< trialPoints
845  << " locations queried, " << initialPoints.size() - initialSize
846  << " points placed, ("
847  << scalar(initialPoints.size() - initialSize)
848  /scalar(max(trialPoints, 1))
849  << " success rate)." << nl
850  << "minCellSize " << minCellSize
851  << ", maxCellSize " << maxCellSize
852  << ", ratio " << maxCellSize/minCellSize
853  << nl << endl;
854  }
855 
856  return true;
857 }
858 
859 
860 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
861 
863 (
864  const dictionary& initialPointsDict,
865  const Time& runTime,
866  Random& rndGen,
867  const conformationSurfaces& geometryToConformTo,
868  const cellShapeControl& cellShapeControls,
869  const autoPtr<backgroundMeshDecomposition>& decomposition
870 )
871 :
872  initialPointsMethod
873  (
874  typeName,
875  initialPointsDict,
876  runTime,
877  rndGen,
878  geometryToConformTo,
879  cellShapeControls,
880  decomposition
881  ),
882  globalTrialPoints_(0),
883  minCellSizeLimit_
884  (
885  detailsDict().lookupOrDefault<scalar>("minCellSizeLimit", 0.0)
886  ),
887  minLevels_(readLabel(detailsDict().lookup("minLevels"))),
888  maxSizeRatio_(readScalar(detailsDict().lookup("maxSizeRatio"))),
889  volRes_(readLabel(detailsDict().lookup("sampleResolution"))),
890  surfRes_
891  (
892  detailsDict().lookupOrDefault<label>("surfaceSampleResolution", volRes_)
893  )
894 {
895  if (maxSizeRatio_ <= 1.0)
896  {
897  maxSizeRatio_ = 2.0;
898 
900  << "The maxSizeRatio must be greater than one to be sensible, "
901  << "setting to " << maxSizeRatio_
902  << endl;
903  }
904 }
905 
906 
907 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
908 
909 List<Vb::Point> autoDensity::initialPoints() const
910 {
911  treeBoundBox hierBB;
912 
913  // Pick up the bounds of this processor, or the whole geometry, depending
914  // on whether this is a parallel run.
915  if (Pstream::parRun())
916  {
917  hierBB = decomposition().procBounds();
918  }
919  else
920  {
921  // Extend the global box to move it off large plane surfaces
922  hierBB = geometryToConformTo().globalBounds().extend
923  (
924  rndGen(),
925  1e-6
926  );
927  }
928 
929  DynamicList<Vb::Point> initialPoints;
930 
931  Info<< nl << " " << typeName << endl;
932 
933  if (debug)
934  {
935  Pout<< " Filling box " << hierBB << endl;
936  }
937 
938  label treeDepth = recurseAndFill
939  (
940  initialPoints,
941  hierBB,
942  minLevels_ - 1,
943  "recursionBox"
944  );
945 
946  initialPoints.shrink();
947 
948  label nInitialPoints = initialPoints.size();
949 
950  if (Pstream::parRun())
951  {
952  reduce(nInitialPoints, sumOp<label>());
953  reduce(globalTrialPoints_, sumOp<label>());
954  }
955 
957  << indent << nInitialPoints << " points placed" << nl
958  << indent << globalTrialPoints_ << " locations queried" << nl
959  << indent
960  << scalar(nInitialPoints)/scalar(max(globalTrialPoints_, 1))
961  << " success rate" << nl
962  << indent
963  << returnReduce(treeDepth, maxOp<label>())
964  << " levels of recursion (maximum)"
965  << decrIndent << decrIndent
966  << endl;
967 
968  return initialPoints;
969 }
970 
971 
972 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
973 
974 } // End namespace Foam
975 
976 // ************************************************************************* //
cachedRandom rndGen(label(0),-1)
void shuffle(UList< T > &)
Definition: UList.C:143
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
uint8_t direction
Definition: direction.H:46
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:223
const double e
Elementary charge.
Definition: doubleFloat.H:78
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static void timeCheck(const Time &runTime, const string &description=string::null, const bool check=true)
Write the elapsedCpuTime and memory usage, with an optional.
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
label k
Boltzmann constant.
static const Form min
Definition: VectorSpace.H:113
Macros for easy insertion into run-time selection tables.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
stressControl lookup("compactNormalStress") >> compactNormalStress
autoDensity(const dictionary &initialPointsDict, const Time &runTime, Random &rndGen, const conformationSurfaces &geometryToConformTo, const cellShapeControl &cellShapeControls, const autoPtr< backgroundMeshDecomposition > &decomposition)
Construct from components.
Pre-declare SubField and related Field type.
Definition: Field.H:57
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
virtual List< Vb::Point > initialPoints() const
Return the initial points for the conformalVoronoiMesh.
static const edgeList edges
Edge to point addressing.
Definition: treeBoundBox.H:169
static const zero Zero
Definition: zero.H:91
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
Definition: doubleScalar.H:63
label readLabel(Istream &is)
Definition: label.H:64
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
static const char nl
Definition: Ostream.H:262
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:237
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar pow3(const dimensionedScalar &ds)
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:393
vector point
Point is a vector.
Definition: point.H:41
#define WarningInFunction
Report a warning using Foam::Warning.
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:230
Namespace for OpenFOAM.