67 void Foam::conformalVoronoiMesh::cellSizeMeshOverlapsBackground()
const 69 const cellShapeControlMesh& cellSizeMesh =
70 cellShapeControl_.shapeControlMesh();
72 DynamicList<Foam::point> pts(number_of_vertices());
76 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
77 vit != finite_vertices_end();
81 if (vit->internalOrBoundaryPoint() && !vit->referred())
83 pts.append(
topoint(vit->point()));
89 boundBox cellSizeMeshBb = cellSizeMesh.bounds();
91 bool fullyContained =
true;
93 if (!cellSizeMeshBb.contains(bb))
95 Pout<<
"Triangulation not fully contained in cell size mesh." 98 Pout<<
"Cell Size Mesh Bounds = " << cellSizeMesh.bounds() <<
endl;
99 Pout<<
"foamyHexMesh Bounds = " << bb <<
endl;
101 fullyContained =
false;
104 reduce(fullyContained, andOp<unsigned int>());
106 Info<<
"Triangulation is " 107 << (fullyContained ?
"fully" :
"not fully")
108 <<
" contained in the cell size mesh" 113 void Foam::conformalVoronoiMesh::insertInternalPoints
119 label nPoints = points.size();
123 reduce(nPoints, sumOp<label>());
126 Info<<
" " << nPoints <<
" points to insert..." <<
endl;
130 List<Foam::point> transferPoints(points.size());
134 transferPoints[pI] =
topoint(points[pI]);
141 decomposition_().distributePoints(transferPoints)
144 transferPoints.clear();
146 map().distribute(points);
149 label nVert = number_of_vertices();
151 insert(points.begin(), points.end());
153 label nInserted(number_of_vertices() - nVert);
157 reduce(nInserted, sumOp<label>());
160 Info<<
" " << nInserted <<
" points inserted" 161 <<
", failed to insert " << nPoints - nInserted
163 << 100.0*(nPoints - nInserted)/(nInserted + small)
168 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
169 vit != finite_vertices_end();
175 vit->index() = getNewVertexIndex();
191 autoPtr<mapDistribute> mapDist =
192 decomposition_().distributePoints(vertices);
206 label preReinsertionSize(number_of_vertices());
208 Map<label> oldToNewIndices =
213 label(number_of_vertices()) - preReinsertionSize,
217 Info<<
" Reinserted " << nReinserted <<
" vertices out of " 221 return oldToNewIndices;
225 void Foam::conformalVoronoiMesh::insertSurfacePointPairs
227 const pointIndexHitAndFeatureList& surfaceHits,
228 const fileName fName,
237 const label featureIndex = surfaceHits[i].second();
239 allGeometry_[featureIndex].getNormal
241 List<pointIndexHit>(1, surfaceHit),
245 const vector& normal = norm[0];
247 const Foam::point& surfacePt(surfaceHit.hitPoint());
250 geometryToConformTo_.meshableSide(featureIndex, surfaceHit);
254 createBafflePointPair
256 pointPairDistance(surfacePt),
267 pointPairDistance(surfacePt),
278 pointPairDistance(surfacePt),
288 << meshableSide <<
", bad" 293 if (foamyHexMeshControls().objOutput() && fName !=
fileName::null)
300 void Foam::conformalVoronoiMesh::insertEdgePointGroups
302 const pointIndexHitAndFeatureList& edgeHits,
303 const fileName fName,
309 if (edgeHits[i].first().hit())
311 const extendedFeatureEdgeMesh& feMesh
313 geometryToConformTo_.features()[edgeHits[i].second()]
328 if (foamyHexMeshControls().objOutput() && fName !=
fileName::null)
335 bool Foam::conformalVoronoiMesh::nearFeaturePt(
const Foam::point& pt)
const 337 scalar exclusionRangeSqr = featurePointExclusionDistanceSqr(pt);
342 geometryToConformTo_.findFeaturePointNearest
354 bool Foam::conformalVoronoiMesh::surfacePtNearFeatureEdge
359 scalar exclusionRangeSqr = surfacePtExclusionDistanceSqr(pt);
364 geometryToConformTo_.findEdgeNearest
376 void Foam::conformalVoronoiMesh::insertInitialPoints()
378 Info<<
nl <<
"Inserting initial points" <<
endl;
380 timeCheck(
"Before initial points call");
382 List<Point> initPts = initialPointsMethod_->initialPoints();
384 timeCheck(
"After initial points call");
388 insertInternalPoints(initPts);
390 if (initialPointsMethod_->fixInitialPoints())
394 Finite_vertices_iterator vit = finite_vertices_begin();
395 vit != finite_vertices_end();
403 if (foamyHexMeshControls().objOutput())
407 time().
path()/
"initialPoints.obj",
415 void Foam::conformalVoronoiMesh::distribute()
422 DynamicList<Foam::point>
points(number_of_vertices());
423 DynamicList<Foam::indexedVertexEnum::vertexType> types
427 DynamicList<scalar> sizes(number_of_vertices());
428 DynamicList<tensor> alignments(number_of_vertices());
432 Finite_vertices_iterator vit = finite_vertices_begin();
433 vit != finite_vertices_end();
439 points.append(
topoint(vit->point()));
440 types.append(vit->type());
441 sizes.append(vit->targetCellSize());
442 alignments.append(vit->alignment());
446 autoPtr<mapDistribute> mapDist =
449 mapDist().distribute(types);
450 mapDist().distribute(sizes);
451 mapDist().distribute(alignments);
456 Info<<
nl <<
" Inserting distributed tessellation" <<
endl;
460 DynamicList<Vb> verticesToInsert(points.size());
464 verticesToInsert.append
475 verticesToInsert.last().targetCellSize() = sizes[pI];
476 verticesToInsert.last().alignment() = alignments[pI];
479 this->rangeInsertWithInfo
481 verticesToInsert.begin(),
482 verticesToInsert.end(),
486 Info<<
" Total number of vertices after redistribution " 489 label(number_of_vertices()), sumOp<label>()
495 void Foam::conformalVoronoiMesh::buildCellSizeAndAlignmentMesh()
497 controlMeshRefinement meshRefinement
502 smoothAlignmentSolver meshAlignmentSmoother
504 cellShapeControl_.shapeControlMesh()
507 meshRefinement.initialMeshPopulation(decomposition_);
509 cellShapeControlMesh& cellSizeMesh = cellShapeControl_.shapeControlMesh();
513 if (!distributeBackground(cellSizeMesh))
516 cellSizeMesh.distribute(decomposition_);
520 const dictionary& motionControlDict
521 = foamyHexMeshControls().foamyHexMeshDict().subDict(
"motionControl");
523 const label nMaxIter =
524 motionControlDict.lookup<
label>(
"maxRefinementIterations");
526 Info<<
"Maximum number of refinement iterations : " << nMaxIter <<
endl;
528 for (
label i = 0; i < nMaxIter; ++i)
530 label nAdded = meshRefinement.refineMesh(decomposition_);
532 reduce(nAdded, sumOp<label>());
536 cellSizeMesh.distribute(decomposition_);
539 Info<<
" Iteration " << i
540 <<
" Added = " << nAdded <<
" points" 552 if (!distributeBackground(cellSizeMesh))
554 cellSizeMesh.distribute(decomposition_);
558 const label maxSmoothingIterations =
559 motionControlDict.lookup<
label>(
"maxSmoothingIterations");
561 meshAlignmentSmoother.smoothAlignments(maxSmoothingIterations);
563 Info<<
"Background cell size and alignment mesh:" <<
endl;
564 cellSizeMesh.printInfo(
Info);
566 Info<<
"Triangulation is " 567 << (cellSizeMesh.is_valid() ?
"valid" :
"not valid!" )
570 if (foamyHexMeshControls().writeCellShapeControlMesh())
573 cellSizeMesh.write();
576 if (foamyHexMeshControls().printVertexInfo())
578 cellSizeMesh.printVertexInfo(
Info);
591 void Foam::conformalVoronoiMesh::setVertexSizeAndAlignment()
593 Info<<
nl <<
"Calculating target cell alignment and size" <<
endl;
597 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
598 vit != finite_vertices_end();
602 if (vit->internalOrBoundaryPoint())
606 cellShapeControls().cellSizeAndAlignment
609 vit->targetCellSize(),
617 Foam::face Foam::conformalVoronoiMesh::buildDualFace
619 const Delaunay::Finite_edges_iterator& eit
622 Cell_circulator ccStart = incident_cells(*eit);
623 Cell_circulator cc1 = ccStart;
624 Cell_circulator cc2 = cc1;
630 DynamicList<label> verticesOnFace;
632 label nUniqueVertices = 0;
638 cc1->hasFarPoint() || cc2->hasFarPoint()
639 || is_infinite(cc1) || is_infinite(cc2)
642 Cell_handle c = eit->first;
643 Vertex_handle vA = c->vertex(eit->second);
644 Vertex_handle vB = c->vertex(eit->third);
650 <<
"Dual face uses circumcenter defined by a " 651 <<
"Delaunay tetrahedron with no internal " 652 <<
"or boundary points. Defining Delaunay edge ends: " 659 label cc1I = cc1->cellIndex();
660 label cc2I = cc2->cellIndex();
664 if (
findIndex(verticesOnFace, cc1I) == -1)
669 verticesOnFace.append(cc1I);
676 }
while (cc1 != ccStart);
678 verticesOnFace.shrink();
680 if (verticesOnFace.size() >= 3 && nUniqueVertices < 3)
695 verticesOnFace.setSize(nUniqueVertices);
698 return face(verticesOnFace);
702 Foam::label Foam::conformalVoronoiMesh::maxFilterCount
704 const Delaunay::Finite_edges_iterator& eit
707 Cell_circulator ccStart = incident_cells(*eit);
708 Cell_circulator cc = ccStart;
714 if (cc->hasFarPoint())
716 Cell_handle c = eit->first;
717 Vertex_handle vA = c->vertex(eit->second);
718 Vertex_handle vB = c->vertex(eit->third);
721 <<
"Dual face uses circumcenter defined by a " 722 <<
"Delaunay tetrahedron with no internal " 723 <<
"or boundary points. Defining Delaunay edge ends: " 729 if (cc->filterCount() > maxFC)
731 maxFC = cc->filterCount();
736 }
while (cc != ccStart);
742 bool Foam::conformalVoronoiMesh::ownerAndNeighbour
750 bool reverse =
false;
755 label dualCellIndexA = vA->index();
757 if (!vA->internalOrBoundaryPoint() || vA->referred())
759 if (!vA->constrained())
765 label dualCellIndexB = vB->index();
767 if (!vB->internalOrBoundaryPoint() || vB->referred())
769 if (!vB->constrained())
775 if (dualCellIndexA == -1 && dualCellIndexB == -1)
778 <<
"Attempting to create a face joining " 779 <<
"two unindexed dual cells " 782 else if (dualCellIndexA == -1 || dualCellIndexB == -1)
786 if (dualCellIndexA == -1)
788 owner = dualCellIndexB;
794 owner = dualCellIndexA;
801 if (dualCellIndexB > dualCellIndexA)
803 owner = dualCellIndexA;
804 neighbour = dualCellIndexB;
808 owner = dualCellIndexB;
809 neighbour = dualCellIndexA;
825 const dictionary& foamyHexMeshDict
828 DistributedDelaunayMesh<
Delaunay>(runTime),
830 rndGen_(64293*Pstream::myProcNo()),
831 foamyHexMeshControls_(foamyHexMeshDict),
836 "cvSearchableSurfaces",
838 searchableSurface::geometryDir(runTime),
843 foamyHexMeshDict.subDict(
"geometry"),
844 foamyHexMeshDict.lookupOrDefault(
"singleRegionName", true)
851 foamyHexMeshDict.subDict(
"surfaceConformation")
856 ? new backgroundMeshDecomposition
860 geometryToConformTo_,
861 foamyHexMeshControls().foamyHexMeshDict().subDict
863 "backgroundMeshDecomposition" 871 foamyHexMeshControls_,
877 ftPtConformer_(*this),
878 edgeLocationTreePtr_(),
879 surfacePtLocationTreePtr_(),
880 surfaceConformationVertices_(),
883 initialPointsMethod::
New 885 foamyHexMeshDict.subDict(
"initialPoints"),
888 geometryToConformTo_,
897 foamyHexMeshDict.subDict(
"motionControl"),
903 faceAreaWeightModel::
New 905 foamyHexMeshDict.subDict(
"motionControl")
926 buildCellSizeAndAlignmentMesh();
928 insertInitialPoints();
930 insertFeaturePoints(
true);
932 setVertexSizeAndAlignment();
934 cellSizeMeshOverlapsBackground();
939 distributeBackground(*
this);
941 buildSurfaceConformation();
945 distributeBackground(*
this);
949 sync(decomposition_().procBounds());
954 storeSurfaceConformation();
960 cellSizeMeshOverlapsBackground();
986 new backgroundMeshDecomposition
990 geometryToConformTo_,
993 "backgroundMeshDecomposition" 999 insertInitialPoints();
1001 insertFeaturePoints();
1006 distributeBackground(*
this);
1008 buildSurfaceConformation();
1012 distributeBackground(*
this);
1016 sync(decomposition_().procBounds());
1019 cellSizeMeshOverlapsBackground();
1032 scalar relaxation = relaxationModel_->relaxation();
1034 Info<<
nl <<
"Relaxation = " << relaxation <<
endl;
1036 pointField dualVertices(number_of_finite_cells());
1043 Delaunay::Finite_cells_iterator cit = finite_cells_begin();
1044 cit != finite_cells_end();
1050 if (cit->anyInternalOrBoundaryDualVertex())
1054 dualVertices[cit->cellIndex()] = cit->dual();
1057 if (cit->hasFarPoint())
1065 setVertexSizeAndAlignment();
1067 timeCheck(
"Determined sizes and alignments");
1069 Info<<
nl <<
"Determining vertex displacements" <<
endl;
1073 cartesianDirections[0] =
vector(1, 0, 0);
1074 cartesianDirections[1] =
vector(0, 1, 0);
1075 cartesianDirections[2] =
vector(0, 0, 1);
1079 number_of_vertices(),
1083 PackedBoolList pointToBeRetained
1085 number_of_vertices(),
1089 DynamicList<Point> pointsToInsert(number_of_vertices());
1093 Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1094 eit != finite_edges_end();
1105 vA->internalPoint() && !vA->referred()
1106 && vB->internalOrBoundaryPoint()
1109 vB->internalPoint() && !vB->referred()
1110 && vA->internalOrBoundaryPoint()
1117 Field<vector> alignmentDirsA
1119 vA->alignment().T() & cartesianDirections
1121 Field<vector> alignmentDirsB
1123 vB->alignment().T() & cartesianDirections
1126 Field<vector> alignmentDirs(alignmentDirsA);
1128 forAll(alignmentDirsA, aA)
1130 const vector& a = alignmentDirsA[aA];
1132 scalar maxDotProduct = 0.0;
1134 forAll(alignmentDirsB, aB)
1136 const vector& b = alignmentDirsB[aB];
1138 const scalar dotProduct = a &
b;
1140 if (
mag(dotProduct) > maxDotProduct)
1142 maxDotProduct =
mag(dotProduct);
1144 alignmentDirs[aA] = a +
sign(dotProduct)*
b;
1146 alignmentDirs[aA] /=
mag(alignmentDirs[aA]);
1153 scalar rABMag =
mag(rAB);
1161 vA->internalPoint() && !vA->referred() && !vA->fixed()
1162 && vB->internalPoint() && !vB->referred() && !vB->fixed()
1172 pointToBeRetained[vA->index()] ==
true 1173 && pointToBeRetained[vB->index()] == true
1178 if (internalPointIsInside(pt))
1180 pointsToInsert.append(
toPoint(pt));
1185 if (vA->internalPoint() && !vA->referred() && !vA->fixed())
1187 pointToBeRetained[vA->index()] =
false;
1190 if (vB->internalPoint() && !vB->referred() && !vB->fixed())
1192 pointToBeRetained[vB->index()] =
false;
1200 forAll(alignmentDirs, aD)
1202 vector& alignmentDir = alignmentDirs[aD];
1204 scalar dotProd = rAB & alignmentDir;
1214 const scalar alignmentDotProd = dotProd/rABMag;
1222 scalar targetCellSize =
1225 scalar targetFaceArea =
sqr(targetCellSize);
1227 const vector originalAlignmentDir = alignmentDir;
1240 vector delta = alignmentDir - 0.5*rAB;
1242 face dualFace = buildDualFace(eit);
1248 const scalar faceArea = dualFace.mag(dualVertices);
1253 originalAlignmentDir,
1276 delta *= faceAreaWeightModel_->faceAreaWeight
1278 faceArea/targetFaceArea
1284 (vA->internalPoint() && vB->internalPoint())
1285 && (!vA->referred() || !vB->referred())
1315 if (internalPointIsInside(newPt))
1327 pointsToInsert.append(
toPoint(newPt));
1332 pointsToInsert.append(
toPoint(newPt));
1340 (vA->internalPoint() && !vA->referred())
1341 || (vB->internalPoint() && !vB->referred())
1370 pointToBeRetained[vA->index()] ==
true 1371 && pointToBeRetained[vB->index()] == true
1376 if (internalPointIsInside(pt))
1378 pointsToInsert.append(
toPoint(pt));
1390 pointToBeRetained[vA->index()] =
false;
1400 pointToBeRetained[vB->index()] =
false;
1414 displacementAccumulator[vA->index()] += 2*
delta;
1418 displacementAccumulator[vA->index()] +=
delta;
1431 displacementAccumulator[vB->index()] -= 2*
delta;
1435 displacementAccumulator[vB->index()] -=
delta;
1444 Info<<
"Limit displacements" <<
endl;
1449 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1450 vit != finite_vertices_end();
1454 if (vit->internalPoint() && !vit->referred() && !vit->fixed())
1456 if (pointToBeRetained[vit->index()] ==
true)
1461 displacementAccumulator[vit->index()]
1467 vector totalDisp =
gSum(displacementAccumulator);
1468 scalar totalDist =
gSum(
mag(displacementAccumulator));
1470 displacementAccumulator *= relaxation;
1472 Info<<
"Sum displacements" <<
endl;
1474 label nPointsToRetain = 0;
1475 label nPointsToRemove = 0;
1479 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1480 vit != finite_vertices_end();
1484 if (vit->internalPoint() && !vit->referred() && !vit->fixed())
1486 if (pointToBeRetained[vit->index()] ==
true)
1500 + displacementAccumulator[vit->index()]
1503 if (internalPointIsInside(pt))
1505 pointsToInsert.append(
toPoint(pt));
1514 pointsToInsert.shrink();
1518 nPointsToRetain - nPointsToRemove,
1521 <<
" internal points are outside the domain. " 1522 <<
"They will not be inserted." <<
endl;
1527 Info<<
"Writing point displacement vectors to file." <<
endl;
1535 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1536 vit != finite_vertices_end();
1540 if (vit->internalPoint() && !vit->referred())
1542 if (pointToBeRetained[vit->index()] ==
true)
1547 << displacementAccumulator[vit->index()][0] <<
" " 1548 << displacementAccumulator[vit->index()][1] <<
" " 1549 << displacementAccumulator[vit->index()][2] <<
" " 1561 Info<<
nl<<
"Inserting displaced tessellation" <<
endl;
1563 insertFeaturePoints(
true);
1565 insertInternalPoints(pointsToInsert,
true);
1595 Delaunay::Finite_vertices_iterator vit = finite_vertices_begin();
1596 vit != finite_vertices_end();
1600 if (!vit->referred() && !usedIndices.insert(vit->index()))
1603 <<
"Index already used! Could not insert: " <<
nl 1621 if (reconformToSurface())
1647 OBJstream multipleIntersections
1649 "multipleIntersections_" 1656 Delaunay::Finite_edges_iterator eit = finite_edges_begin();
1657 eit != finite_edges_end();
1668 List<pointIndexHit> surfHits;
1681 surfHits.size() >= 2
1683 surfHits.size() == 0
1685 (vA->externalBoundaryPoint()
1686 && vB->internalBoundaryPoint())
1687 || (vB->externalBoundaryPoint()
1688 && vA->internalBoundaryPoint())
1706 if (
time().writeTime())
1712 <<
"Total displacement = " << totalDisp <<
nl 1713 <<
"Total distance = " << totalDist <<
nl 1718 void Foam::conformalVoronoiMesh::checkCoPlanarCells()
const 1720 typedef CGAL::Exact_predicates_exact_constructions_kernel Kexact;
1721 typedef CGAL::Point_3<Kexact> PointExact;
1725 Pout<<
"Triangulation is invalid!" <<
endl;
1728 OFstream str(
"badCells.obj");
1734 Delaunay::Finite_cells_iterator cit = finite_cells_begin();
1735 cit != finite_cells_end();
1744 <<
" quality = " << quality <<
nl 1745 <<
" dual = " <<
topoint(cit->dual()) << endl;
1749 FixedList<PointExact, 4> cellVerticesExact(PointExact(0,0,0));
1750 forAll(cellVerticesExact, vI)
1752 cellVerticesExact[vI] = PointExact
1754 cit->vertex(vI)->point().x(),
1755 cit->vertex(vI)->point().y(),
1756 cit->vertex(vI)->point().z()
1760 PointExact synchronisedDual = CGAL::circumcenter<Kexact>
1762 cellVerticesExact[0],
1763 cellVerticesExact[1],
1764 cellVerticesExact[2],
1765 cellVerticesExact[3]
1770 CGAL::to_double(synchronisedDual.x()),
1771 CGAL::to_double(synchronisedDual.y()),
1772 CGAL::to_double(synchronisedDual.z())
1775 Info<<
"inexact = " << cit->dual() <<
nl 1776 <<
"exact = " << exactPt <<
endl;
1780 Pout<<
"There are " << badCells <<
" bad cells out of " 1781 << number_of_finite_cells() <<
endl;
1784 label nNonGabriel = 0;
1787 Delaunay::Finite_facets_iterator fit = finite_facets_begin();
1788 fit != finite_facets_end();
1792 if (!is_Gabriel(*fit))
1798 Pout<<
"There are " << nNonGabriel <<
" non-Gabriel faces out of " 1799 << number_of_finite_facets() <<
endl;
dimensionedScalar sign(const dimensionedScalar &ds)
CGAL::Delaunay_triangulation_3< K, Tds, CompactLocator > Delaunay
label cellCount() const
Return the cell count (the next unique cell index)
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
pointFromPoint topoint(const Point &P)
errorManipArg< error, int > exit(error &err, const int errNo=1)
A face is a list of labels corresponding to mesh vertices.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
static const fileName null
An empty fileName.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
const dimensionedScalar b
Wien displacement law constant: default SI units: [m K].
Ostream & endl(Ostream &os)
Add newline and flush stream.
bool distribute(const boundBox &bb)
scalar faceAreaRatioCoeff() const
Return the faceAreaRatioCoeff.
PointIndexHit< point > pointIndexHit
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Vector< scalar > vector
A scalar version of the templated Vector.
Initialise the NamedEnum HashTable from the static list of names.
scalar removalDistCoeff() const
Return removalDistCoeff.
void printVertexInfo(Ostream &os) const
Write vertex statistics in the form of a table to stream.
void insert(const scalar, DynamicList< floatScalar > &)
Append scalar to given DynamicList.
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
An indexed form of CGAL::Triangulation_vertex_base_3<K> used to keep track of the Delaunay vertices i...
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
vectorField pointField
pointField is a vectorField.
bool uninitialised(const VertexType &v)
Type gSum(const FieldField< Field, Type > &f)
line< point, const point & > linePointRef
Line using referred points.
void updateCellSizeAndFaceArea(vector &alignmentDir, scalar &targetFaceArea, scalar &targetCellSize) const
const cellAspectRatioControl & aspectRatio() const
sideVolumeType
Normals point to the outside.
void reset()
Clear the entire triangulation.
PointFrompoint toPoint(const Foam::point &p)
List< label > labelList
A List of labels.
errorManip< error > abort(error &err)
Istream and Ostream manipulators taking arguments.
void reverse(UList< T > &, const label n)
defineTypeNameAndDebug(combustionModel, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
scalar cosInsertionAcceptanceAngle() const
Return the cosInsertionAcceptanceAngle.
void resetCellCount()
Set the cell count to zero.
static bool & parRun()
Is this a parallel run?
Foam::scalar averageCellSize(const VertexType &vA, const VertexType &vB)
Return the target cell size from that stored on a pair of Delaunay vertices,.
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
void updateDeltaVector(const vector &alignmentDir, const scalar targetCellSize, const scalar rABMag, vector &delta) const
label getNewCellIndex() const
Create a new unique cell index and return.
dimensioned< scalar > mag(const dimensioned< Type > &)
Field< vector > vectorField
Specialisation of Field<T> for vector.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const doubleScalar e
Elementary charge.
scalar coplanarTet(Cell &c, const scalar tol=1e-12)
Map< label > insertPoints(const List< Vb > &vertices, const bool reIndex)
Insert the list of vertices (calls rangeInsertWithInfo)
scalar insertionDistCoeff() const
Return the insertionDistCoeff.
void sync(const boundBox &bb)
Refer vertices so that the processor interfaces are consistent.
InfoProxy< IOstream > info() const
Return info proxy.
fileName path(UMean.rootPath()/UMean.caseName()/functionObjects::writeFile::outputPrefix/"graphs"/UMean.instance())
A HashTable to objects of type <T> with a label key.