cellShapeControlMesh.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) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2022 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 "cellShapeControlMesh.H"
31 #include "pointIOField.H"
32 #include "scalarIOField.H"
33 #include "triadIOField.H"
34 #include "tetrahedron.H"
35 #include "plane.H"
36 #include "transform.H"
37 #include "meshTools.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 defineTypeNameAndDebug(cellShapeControlMesh, 0);
44 }
45 
46 Foam::word Foam::cellShapeControlMesh::meshSubDir = "cellShapeControlMesh";
47 
48 
49 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
50 
51 //Foam::tensor Foam::cellShapeControlMesh::requiredAlignment
52 //(
53 // const Foam::point& pt,
54 // const searchableSurfaces& allGeometry,
55 // const conformationSurfaces& geometryToConformTo
56 //) const
57 //{
58 // pointIndexHit surfHit;
59 // label hitSurface;
60 //
61 // geometryToConformTo.findSurfaceNearest
62 // (
63 // pt,
64 // sqr(GREAT),
65 // surfHit,
66 // hitSurface
67 // );
68 //
69 // if (!surfHit.hit())
70 // {
71 // FatalErrorInFunction
72 // << "findSurfaceNearest did not find a hit across the surfaces."
73 // << exit(FatalError) << endl;
74 // }
75 //
76 // // Primary alignment
77 //
78 // vectorField norm(1);
79 //
80 // allGeometry[hitSurface].getNormal
81 // (
82 // List<pointIndexHit>(1, surfHit),
83 // norm
84 // );
85 //
86 // const vector np = norm[0];
87 //
88 // // Generate equally spaced 'spokes' in a circle normal to the
89 // // direction from the vertex to the closest point on the surface
90 // // and look for a secondary intersection.
91 //
92 // const vector d = surfHit.hitPoint() - pt;
93 //
94 // const tensor Rp = rotationTensor(vector(0,0,1), np);
95 //
96 // const label s = 36;//foamyHexMeshControls().alignmentSearchSpokes();
97 //
98 // scalar closestSpokeHitDistance = GREAT;
99 //
100 // pointIndexHit closestSpokeHit;
101 //
102 // label closestSpokeSurface = -1;
103 //
104 // const scalar spanMag = geometryToConformTo.globalBounds().mag();
105 //
106 // for (label i = 0; i < s; i++)
107 // {
108 // vector spoke
109 // (
110 // Foam::cos(i*constant::mathematical::twoPi/s),
111 // Foam::sin(i*constant::mathematical::twoPi/s),
112 // 0
113 // );
114 //
115 // spoke *= spanMag;
116 //
117 // spoke = Rp & spoke;
118 //
119 // pointIndexHit spokeHit;
120 //
121 // label spokeSurface = -1;
122 //
123 // // internal spoke
124 //
125 // geometryToConformTo.findSurfaceNearestIntersection
126 // (
127 // pt,
128 // pt + spoke,
129 // spokeHit,
130 // spokeSurface
131 // );
132 //
133 // if (spokeHit.hit())
134 // {
135 // scalar spokeHitDistance = spokeHit.point().dist(pt);
136 //
137 // if (spokeHitDistance < closestSpokeHitDistance)
138 // {
139 // closestSpokeHit = spokeHit;
140 // closestSpokeSurface = spokeSurface;
141 // closestSpokeHitDistance = spokeHitDistance;
142 // }
143 // }
144 //
145 // //external spoke
146 //
147 // Foam::point mirrorPt = pt + 2*d;
148 //
149 // geometryToConformTo.findSurfaceNearestIntersection
150 // (
151 // mirrorPt,
152 // mirrorPt + spoke,
153 // spokeHit,
154 // spokeSurface
155 // );
156 //
157 // if (spokeHit.hit())
158 // {
159 // scalar spokeHitDistance = spokeHit.point().dist(mirrorPt);
160 //
161 // if (spokeHitDistance < closestSpokeHitDistance)
162 // {
163 // closestSpokeHit = spokeHit;
164 // closestSpokeSurface = spokeSurface;
165 // closestSpokeHitDistance = spokeHitDistance;
166 // }
167 // }
168 // }
169 //
170 // if (closestSpokeSurface == -1)
171 // {
177 //
178 // return I;
179 // }
180 //
181 // // Auxiliary alignment generated by spoke intersection normal.
182 //
183 // allGeometry[closestSpokeSurface].getNormal
184 // (
185 // List<pointIndexHit>(1, closestSpokeHit),
186 // norm
187 // );
188 //
189 // const vector& na = norm[0];
190 //
191 // // Secondary alignment
192 // vector ns = np ^ na;
193 //
194 // if (mag(ns) < SMALL)
195 // {
196 // FatalErrorInFunction
197 // << "Parallel normals detected in spoke search." << nl
198 // << "point: " << pt << nl
199 // << "closest surface point: " << surfHit.point() << nl
200 // << "closest spoke hit: " << closestSpokeHit.point() << nl
201 // << "np: " << surfHit.point() + np << nl
202 // << "ns: " << closestSpokeHit.point() + na << nl
203 // << exit(FatalError);
204 // }
205 //
206 // ns /= mag(ns);
207 //
208 // tensor Rs = rotationTensor((Rp & vector(0,1,0)), ns);
209 //
210 // return (Rs & Rp);
211 //}
212 
213 
215 {
216  label nRemoved = 0;
217  for
218  (
219  CellSizeDelaunay::Finite_vertices_iterator vit =
220  finite_vertices_begin();
221  vit != finite_vertices_end();
222  ++vit
223  )
224  {
225  std::list<Vertex_handle> verts;
226  adjacent_vertices(vit, std::back_inserter(verts));
227 
228  bool removePt = true;
229  for
230  (
231  std::list<Vertex_handle>::iterator aVit = verts.begin();
232  aVit != verts.end();
233  ++aVit
234  )
235  {
236  Vertex_handle avh = *aVit;
237 
238  scalar diff =
239  mag(avh->targetCellSize() - vit->targetCellSize())
240  /max(vit->targetCellSize(), 1e-6);
241 
242  if (diff > 0.05)
243  {
244  removePt = false;
245  }
246  }
247 
248  if (removePt)
249  {
250  remove(vit);
251  nRemoved++;
252  }
253  }
254 
255  return nRemoved;
256 }
257 
258 
260 {
261  tmp<pointField> tcellCentres(new pointField(number_of_finite_cells()));
262  pointField& cellCentres = tcellCentres.ref();
263 
264  label count = 0;
265  for
266  (
267  CellSizeDelaunay::Finite_cells_iterator c = finite_cells_begin();
268  c != finite_cells_end();
269  ++c
270  )
271  {
272  if (c->hasFarPoint())
273  {
274  continue;
275  }
276 
277  scalarList bary;
279 
280  const Foam::point centre = topoint
281  (
282  CGAL::centroid<baseK>
283  (
284  c->vertex(0)->point(),
285  c->vertex(1)->point(),
286  c->vertex(2)->point(),
287  c->vertex(3)->point()
288  )
289  );
290 
291  cellCentres[count++] = centre;
292  }
293 
294  cellCentres.resize(count);
295 
296  return tcellCentres;
297 }
298 
299 
301 {
302  OFstream str
303  (
304  "refinementTriangulation_"
306  + ".obj"
307  );
308 
309  label count = 0;
310 
311  Info<< "Write refinementTriangulation" << endl;
312 
313  for
314  (
315  CellSizeDelaunay::Finite_edges_iterator e = finite_edges_begin();
316  e != finite_edges_end();
317  ++e
318  )
319  {
320  Cell_handle c = e->first;
321  Vertex_handle vA = c->vertex(e->second);
322  Vertex_handle vB = c->vertex(e->third);
323 
324  // Don't write far edges
325  if (vA->farPoint() || vB->farPoint())
326  {
327  continue;
328  }
329 
330  // Don't write unowned edges
331  if (vA->referred() && vB->referred())
332  {
333  continue;
334  }
335 
336  pointFromPoint p1 = topoint(vA->point());
337  pointFromPoint p2 = topoint(vB->point());
338 
339  meshTools::writeOBJ(str, p1, p2, count);
340  }
341 
342  if (is_valid())
343  {
344  Info<< " Triangulation is valid" << endl;
345  }
346  else
347  {
349  << "Triangulation is not valid"
350  << abort(FatalError);
351  }
352 }
353 
354 
355 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
356 
357 Foam::cellShapeControlMesh::cellShapeControlMesh(const Time& runTime)
358 :
359  DistributedDelaunayMesh<CellSizeDelaunay>
360  (
361  runTime,
362  meshSubDir
363  ),
364  runTime_(runTime)
365 {
366  if (this->vertexCount())
367  {
368  fvMesh mesh
369  (
370  IOobject
371  (
372  meshSubDir,
373  runTime.timeName(),
374  runTime,
377  )
378  );
379 
380  if (mesh.nPoints() == this->vertexCount())
381  {
382  IOobject io
383  (
384  "sizes",
385  runTime.timeName(),
386  meshSubDir,
387  runTime,
391  );
392 
393  if (io.typeHeaderOk<pointScalarField>(true))
394  {
396 
397  triadIOField alignments
398  (
399  IOobject
400  (
401  "alignments",
402  mesh.time().timeName(),
403  meshSubDir,
404  mesh.time(),
408  )
409  );
410 
411  if (alignments.size() == this->vertexCount())
412  {
413  for
414  (
415  Finite_vertices_iterator vit = finite_vertices_begin();
416  vit != finite_vertices_end();
417  ++vit
418  )
419  {
420  vit->targetCellSize() = sizes[vit->index()];
421  vit->alignment() = alignments[vit->index()];
422  }
423  }
424  else
425  {
427  << "Cell alignments point field " << alignments.size()
428  << " is not the same size as the number of vertices"
429  << " in the mesh " << this->vertexCount()
430  << abort(FatalError);
431  }
432  }
433  }
434  }
435 }
436 
437 
438 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
439 
441 (
442  const Foam::point& pt,
443  barycentric& bary,
444  Cell_handle& ch
445 ) const
446 {
447  // Use the previous cell handle as a hint on where to start searching
448  // Giving a hint causes strange errors...
449  ch = locate(toPoint(pt));
450 
451  if (dimension() > 2 && !is_infinite(ch))
452  {
453  oldCellHandle_ = ch;
454 
455  tetPointRef tet
456  (
457  topoint(ch->vertex(0)->point()),
458  topoint(ch->vertex(1)->point()),
459  topoint(ch->vertex(2)->point()),
460  topoint(ch->vertex(3)->point())
461  );
462 
463  bary = tet.pointToBarycentric(pt);
464  }
465 }
466 
467 
469 {
470  DynamicList<Foam::point> pts(number_of_vertices());
471 
472  for
473  (
474  Finite_vertices_iterator vit = finite_vertices_begin();
475  vit != finite_vertices_end();
476  ++vit
477  )
478  {
479  if (vit->real())
480  {
481  pts.append(topoint(vit->point()));
482  }
483  }
484 
485  boundBox bb(pts);
486 
487  return bb;
488 }
489 
490 
492 (
493  const backgroundMeshDecomposition& decomposition
494 )
495 {
496  DynamicList<Foam::point> points(number_of_vertices());
497  DynamicList<scalar> sizes(number_of_vertices());
498  DynamicList<tensor> alignments(number_of_vertices());
499 
500  DynamicList<Vb> farPts(8);
501 
502  for
503  (
504  Finite_vertices_iterator vit = finite_vertices_begin();
505  vit != finite_vertices_end();
506  ++vit
507  )
508  {
509  if (vit->real())
510  {
511  points.append(topoint(vit->point()));
512  sizes.append(vit->targetCellSize());
513  alignments.append(vit->alignment());
514  }
515  else if (vit->farPoint())
516  {
517  farPts.append
518  (
519  Vb
520  (
521  vit->point(),
522  -1,
523  Vb::vtFar,
525  )
526  );
527 
528  farPts.last().targetCellSize() = vit->targetCellSize();
529  farPts.last().alignment() = vit->alignment();
530  }
531  }
532 
533  autoPtr<mapDistribute> mapDist =
535  (
536  decomposition,
537  points
538  );
539 
540  mapDist().distribute(sizes);
541  mapDist().distribute(alignments);
542 
543  // Reset the entire tessellation
545 
546 
547  // Internal points have to be inserted first
548  DynamicList<Vb> verticesToInsert(points.size());
549 
550 
551  forAll(farPts, ptI)
552  {
553  verticesToInsert.append(farPts[ptI]);
554  }
555 
556 
557  forAll(points, pI)
558  {
559  verticesToInsert.append
560  (
561  Vb
562  (
563  toPoint(points[pI]),
564  -1,
567  )
568  );
569 
570  verticesToInsert.last().targetCellSize() = sizes[pI];
571  verticesToInsert.last().alignment() = alignments[pI];
572  }
573 
574  Info<< nl << " Inserting distributed background tessellation..." << endl;
575 
576  this->rangeInsertWithInfo
577  (
578  verticesToInsert.begin(),
579  verticesToInsert.end(),
580  true
581  );
582 
583  sync(decomposition.procBounds());
584 
585  Info<< " Total number of vertices after redistribution "
586  << returnReduce(label(number_of_vertices()), sumOp<label>()) << endl;
587 }
588 
589 
591 {
592  tensorField alignmentsTmp(number_of_vertices(), Zero);
593 
594  label count = 0;
595  for
596  (
597  Finite_vertices_iterator vit = finite_vertices_begin();
598  vit != finite_vertices_end();
599  ++vit
600  )
601  {
602  alignmentsTmp[count++] = vit->alignment();
603  }
604 
605  return alignmentsTmp;
606 }
607 
608 
610 {
611  Info<< "Writing " << meshSubDir << endl;
612 
613  // Reindex the cells
614  label cellCount = 0;
615  for
616  (
617  Finite_cells_iterator cit = finite_cells_begin();
618  cit != finite_cells_end();
619  ++cit
620  )
621  {
622  if (!cit->hasFarPoint() && !is_infinite(cit))
623  {
624  cit->cellIndex() = cellCount++;
625  }
626  }
627 
628  labelPairLookup vertexMap;
629  labelList cellMap;
630 
632  (
633  meshSubDir,
634  vertexMap,
635  cellMap
636  );
637  const polyMesh& mesh = meshPtr();
638 
639  pointScalarField sizes
640  (
641  IOobject
642  (
643  "sizes",
644  mesh.time().timeName(),
645  meshSubDir,
646  mesh.time(),
649  ),
652  );
653 
654  triadIOField alignments
655  (
656  IOobject
657  (
658  "alignments",
659  mesh.time().timeName(),
660  meshSubDir,
661  mesh.time(),
664  ),
665  sizes.size()
666  );
667 
668  // Write alignments
669 // OFstream str(runTime_.path()/"alignments.obj");
670 
671  for
672  (
673  Finite_vertices_iterator vit = finite_vertices_begin();
674  vit != finite_vertices_end();
675  ++vit
676  )
677  {
678  if (!vit->farPoint())
679  {
680  // Populate sizes
681  sizes[vertexMap[labelPair(vit->index(), vit->procIndex())]] =
682  vit->targetCellSize();
683 
684  alignments[vertexMap[labelPair(vit->index(), vit->procIndex())]] =
685  vit->alignment();
686 
687 // // Write alignments
688 // const tensor& alignment = vit->alignment();
689 // pointFromPoint pt = topoint(vit->point());
690 //
691 // if
692 // (
693 // alignment.x() == triad::unset[0]
694 // || alignment.y() == triad::unset[0]
695 // || alignment.z() == triad::unset[0]
696 // )
697 // {
698 // Info<< "Bad alignment = " << vit->info();
699 //
700 // vit->alignment() = tensor::I;
701 //
702 // Info<< "New alignment = " << vit->info();
703 //
704 // continue;
705 // }
706 //
707 // meshTools::writeOBJ(str, pt, alignment.x() + pt);
708 // meshTools::writeOBJ(str, pt, alignment.y() + pt);
709 // meshTools::writeOBJ(str, pt, alignment.z() + pt);
710  }
711  }
712 
713  mesh.write();
714  sizes.write();
715  alignments.write();
716 }
717 
718 
720 (
721  const autoPtr<backgroundMeshDecomposition>& decomposition
722 ) const
723 {
724  // Loop over all the tets and estimate the cell count in each one
725 
726  scalar cellCount = 0;
727 
728  for
729  (
730  Finite_cells_iterator cit = finite_cells_begin();
731  cit != finite_cells_end();
732  ++cit
733  )
734  {
735  if (!cit->hasFarPoint() && !is_infinite(cit))
736  {
737  // TODO: Check if tet centre is on the processor..
738  CGAL::Tetrahedron_3<baseK> tet
739  (
740  cit->vertex(0)->point(),
741  cit->vertex(1)->point(),
742  cit->vertex(2)->point(),
743  cit->vertex(3)->point()
744  );
745 
746  pointFromPoint centre = topoint(CGAL::centroid(tet));
747 
748  if
749  (
751  && !decomposition().positionOnThisProcessor(centre)
752  )
753  {
754  continue;
755  }
756 
757  scalar volume = CGAL::to_double(tet.volume());
758 
759  scalar averagedPointCellSize = 0;
760  //scalar averagedPointCellSize = 1;
761 
762  // Get an average volume by averaging the cell size of the vertices
763  for (label vI = 0; vI < 4; ++vI)
764  {
765  averagedPointCellSize += cit->vertex(vI)->targetCellSize();
766  //averagedPointCellSize *= cit->vertex(vI)->targetCellSize();
767  }
768 
769  averagedPointCellSize /= 4;
770  //averagedPointCellSize = ::sqrt(averagedPointCellSize);
771 
772 // if (averagedPointCellSize < SMALL)
773 // {
774 // Pout<< "Volume = " << volume << endl;
775 //
776 // for (label vI = 0; vI < 4; ++vI)
777 // {
778 // Pout<< "Point " << vI
779 // << ", point = " << topoint(cit->vertex(vI)->point())
780 // << ", size = " << cit->vertex(vI)->targetCellSize()
781 // << endl;
782 // }
783 // }
784 
785  cellCount += volume/pow(averagedPointCellSize, 3);
786  }
787  }
788 
789  return cellCount;
790 }
791 
792 
793 // ************************************************************************* //
List< scalar > scalarList
List of scalar.
Definition: scalarList.H:32
scalar diff(const triad &A, const triad &B)
Return a quantity of the difference between two triads.
Definition: triad.C:373
autoPtr< polyMesh > createMesh(const fileName &name, labelPairLookup &vertexMap, labelList &cellMap, const bool writeDelaunayData=true) const
Create an fvMesh from the triangulation.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
pointFromPoint topoint(const Point &P)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:500
tmp< pointField > cellCentres() const
Get the centres of all the tets.
tetrahedron< point, const point & > tetPointRef
A tetrahedron using referred points.
Definition: tetrahedron.H:72
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1004
bool distribute(const boundBox &bb)
void barycentricCoords(const Foam::point &pt, barycentric &bary, Cell_handle &ch) const
Calculate and return the barycentric coordinates for.
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
Barycentric< scalar > barycentric
A scalar version of the templated Barycentric.
Definition: barycentric.H:41
Ignore writing from objectRegistry::writeObject()
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1029
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:362
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
tensorField dumpAlignments() const
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
Definition: indexedVertex.H:59
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
const pointField & points
Foam::autoPtr< Foam::dynamicFvMesh > meshPtr
A class for handling words, derived from Foam::string.
Definition: word.H:63
3D tensor transformation operations.
Reading is optional [identical to LAZY_READ].
void reset()
Clear the entire triangulation.
PointFrompoint toPoint(const Foam::point &p)
void distribute(const backgroundMeshDecomposition &decomposition)
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
errorManip< error > abort(error &err)
Definition: errorManip.H:139
virtual bool write(const bool writeOnProc=true) const
Write mesh using IO settings from time.
Definition: fvMesh.C:1069
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
IOField< triad > triadIOField
IO for a Field of triad.
Definition: triadIOField.H:32
CellSizeDelaunay::Vertex_handle Vertex_handle
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:770
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
defineTypeNameAndDebug(combustionModel, 0)
CGAL::Delaunay_triangulation_3< K, Tds, FastLocator > CellSizeDelaunay
static word meshSubDir
Return the mesh sub-directory name (usually "cellShapeControlMesh")
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
CellSizeDelaunay::Cell_handle Cell_handle
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionedScalar c
Speed of light in a vacuum.
Nothing to be read.
Automatically write from objectRegistry::writeObject()
HashTable< label, labelPair, Foam::Hash< labelPair > > labelPairLookup
This is a Map of a labelPair to a label. Used for e.g. for face1, face2 to shared edge...
boundBox bounds() const
messageStream Info
Information stream (stdout output on master, null elsewhere)
label estimateCellCount(const autoPtr< backgroundMeshDecomposition > &decomposition) const
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
Do not request registration (bool: false)
Namespace for OpenFOAM.
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133