93#include "CGAL/version.h"
94#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1050211000)
95#define BOOST_BIND_GLOBAL_PLACEHOLDERS
97#pragma clang diagnostic ignored "-Wbitwise-instead-of-logical"
98#pragma clang diagnostic ignored "-Wdeprecated-builtins"
99#pragma clang diagnostic ignored "-Wdeprecated-declarations"
101#include <CGAL/AABB_tree.h>
102#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1060011000)
103#include <CGAL/AABB_traits.h>
105#include <CGAL/AABB_traits_3.h>
107#include <CGAL/AABB_face_graph_triangle_primitive.h>
110typedef CGAL::AABB_face_graph_triangle_primitive
114#if defined(CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1060011000)
115typedef CGAL::AABB_traits<K, Primitive> Traits;
117typedef CGAL::AABB_traits_3<K, Primitive> Traits;
119typedef CGAL::AABB_tree<Traits> Tree;
122#if (CGAL_VERSION_NR < 1060000910)
123typedef boost::optional
128 Tree::Intersection_and_primitive_id<Segment>::Type
129> Segment_intersection;
138bool intersectSurfaces
144 bool hasMoved =
false;
146 for (label iter = 0; iter < 10; iter++)
148 Info<<
"Determining intersections of surface edges with itself" <<
endl;
197 Info<<
"Surface has been moved. Writing to " << newFile <<
endl;
207bool intersectSurfaces
215 bool hasMoved1 =
false;
216 bool hasMoved2 =
false;
218 for (label iter = 0; iter < 10; iter++)
220 Info<<
"Determining intersections of surf1 edges with surf2"
266 Info<<
"Determining intersections of surf2 edges with surf1"
309 if (nIters1 == 0 && nIters2 == 0)
311 Info<<
"** Resolved all intersections to be proper edge-face pierce"
320 Info<<
"Surface 1 has been moved. Writing to " << newFile
322 surf1.
write(newFile);
328 Info<<
"Surface 2 has been moved. Writing to " << newFile
330 surf2.
write(newFile);
333 return hasMoved1 || hasMoved2;
337label calcNormalDirection
340 const vector& otherNormal,
350 label nDir = ((
cross & fC0tofE0) > 0.0 ? 1 : -1);
352 nDir *= ((otherNormal & fC0tofE0) > 0.0 ? -1 : 1);
381 Info<<
"Determining intersections of surf1 edges with surf2 faces"
393 Info<<
"Determining intersections of surf2 edges with surf1 faces"
415 const label startEdgeI,
416 const label startFaceI,
420 const labelList& eFaces =
s.edgeFaces()[startEdgeI];
422 if (eFaces.
size() == 2)
425 if (eFaces[0] == startFaceI)
427 nextFaceI = eFaces[1];
429 else if (eFaces[1] == startFaceI)
431 nextFaceI = eFaces[0];
443 label index =
s.pointFaces()[pointI].find(nextFaceI);
445 if (pFacesZone[index] == -1)
448 pFacesZone[index] = zoneI;
451 const labelList& fEdges =
s.faceEdges()[nextFaceI];
453 label nextEdgeI = -1;
457 label edgeI = fEdges[i];
458 const edge&
e =
s.edges()[edgeI];
460 if (edgeI != startEdgeI && (
e[0] == pointI ||
e[1] == pointI))
471 <<
"Problem: cannot find edge out of " << fEdges
472 <<
"on face " << nextFaceI <<
" that uses point " << pointI
501 List<labelledTri> newFaces(
s);
502 label nNonManifold = 0;
513 for (; index < pFacesZone.
size(); index++)
515 if (pFacesZone[index] == -1)
517 label zoneI = nZones++;
518 pFacesZone[index] = zoneI;
520 label faceI =
pFaces[index];
526 label edgeI = fEdges[fEdgeI];
527 const edge&
e = edges[edgeI];
529 if (
e[0] == pointI ||
e[1] == pointI)
549 for (label zoneI = 1; zoneI < nZones; zoneI++)
551 label newPointI = newPoints.size();
552 newPointMap.append(
s.meshPoints()[pointI]);
553 newPoints.append(
s.points()[
s.meshPoints()[pointI]]);
557 if (pFacesZone[index] == zoneI)
559 label faceI =
pFaces[index];
563 if (localF[fp] == pointI)
565 newFaces[faceI][fp] = newPointI;
576 Info<<
"Duplicating " << nNonManifold <<
" points out of " <<
s.nPoints()
578 if (nNonManifold > 0)
580 triSurface dupSurf(newFaces,
s.patches(), newPoints,
true);
585 const labelList& dupMp = dupSurf.meshPoints();
590 label dupMeshPointI = dupMp[pointI];
591 label meshPointI = newPointMap[dupMeshPointI];
592 dupPointMap[pointI] = mpm[meshPointI];
596 forAll(dupPointMap, pointI)
598 const point& dupPt = dupSurf.points()[dupMp[pointI]];
599 label sLocalPointI = dupPointMap[pointI];
600 label sMeshPointI =
s.meshPoints()[sLocalPointI];
601 const point& sPt =
s.points()[sMeshPointI];
603 if (
mag(dupPt-sPt) > SMALL)
636 List<List<pointIndexHit>> intersections(edges.
size());
642 std::vector<Segment_intersection> segments;
645 const edge&
e = edges[edgeI];
650 K::Segment_3 segment_query
657 tree.all_intersections(segment_query, std::back_inserter(segments));
660 for (
const Segment_intersection& intersect : segments)
666 #
if (CGAL_VERSION_NR < 1060000910)
667 const Point* ptPtr = boost::get<Point>(&(intersect->first))
669 const Point* ptPtr = std::get_if<Point>(&(intersect->first))
675 CGAL::to_double(ptPtr->x()),
676 CGAL::to_double(ptPtr->y()),
677 CGAL::to_double(ptPtr->z())
680 #if defined (CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1041400000)
681 Polyhedron::Face_handle
f = (intersect->second);
684 Polyhedron::Face_handle
f = (intersect->second).first;
687 intersections[edgeI].append
697 classifications[edgeI].append(-1);
703 #
if (CGAL_VERSION_NR < 1060000910)
704 const Segment* sPtr = boost::get<Segment>(&(intersect->first))
706 const Segment* sPtr = std::get_if<Segment>(&(intersect->first))
710 #if defined (CGAL_VERSION_NR) && (CGAL_VERSION_NR < 1041400000)
711 Polyhedron::Face_handle
f = (intersect->second);
714 Polyhedron::Face_handle
f = (intersect->second).first;
726 CGAL::to_double(sPtr->source().x()),
727 CGAL::to_double(sPtr->source().y()),
728 CGAL::to_double(sPtr->source().z())
733 CGAL::to_double(sPtr->target().x()),
734 CGAL::to_double(sPtr->target().y()),
735 CGAL::to_double(sPtr->target().z())
739 intersections[edgeI].append
749 classifications[edgeI].append(2);
778 const Tree
tree(
p.facets_begin(),
p.facets_end(),
p);
786 for (label iter = 0; iter < 5; iter++)
789 Info<<
"Determining intersections of " << surf1.
nEdges()
790 <<
" edges with surface of " << label(
tree.size()) <<
" triangles"
792 cuts = edgeIntersectionsCGAL
799 Info<<
"Determined intersections:" <<
nl
800 <<
" edges : " << surf1.
nEdges() <<
nl
801 <<
" non-degenerate cuts : " << cuts.first() <<
nl
802 <<
" degenerate cuts : " << cuts.second() <<
nl
805 if (cuts.second() == 0)
810 Info<<
"Shuffling conflicting points" <<
endl;
815 const point p05(0.5, 0.5, 0.5);
824 const edge&
e = edges[edgeI];
828 surf1Points[
mp[
e[eI]]] += surf1PointTol[
e[eI]]*d;
862 forAll(subEdges, subEdgeI)
864 const edge& subE = subEdges[subEdgeI];
867 label start = pointMap[subE[0]];
868 label
end = pointMap[subE[1]];
869 const labelList& pEdges = pointEdges[start];
872 label edgeI = pEdges[pEdgeI];
873 const edge&
e = edges[edgeI];
875 if (
e.otherVertex(start) == end)
877 if (edgeMap[subEdgeI] == -1)
879 edgeMap[subEdgeI] = edgeI;
881 else if (edgeMap[subEdgeI] != edgeI)
884 << subE <<
" points:"
886 <<
" matches to " << edgeI
888 <<
" but also " << edgeMap[subEdgeI]
890 << edges[edgeMap[subEdgeI]].line(surf.
localPoints())
897 if (edgeMap[subEdgeI] == -1)
900 << subE <<
" at:" << subSurf.
localPoints()[subE[0]]
923 Info<<
"Intersect surface 1 edges with surface 2:" <<
nl;
924 Info<<
" constructing CGAL surface ..." <<
endl;
928 Info<<
" constructing CGAL tree ..." <<
endl;
929 const Tree
tree(
p.facets_begin(),
p.facets_end(),
p);
931 edgeIntersectionsCGAL
941 Info<<
"Intersect surface 2 edges with surface 1:" <<
nl;
942 Info<<
" constructing CGAL surface ..." <<
endl;
946 Info<<
" constructing CGAL tree ..." <<
endl;
947 const Tree
tree(
p.facets_begin(),
p.facets_end(),
p);
949 edgeIntersectionsCGAL
976 for (label iter = 0; iter < 10; iter++)
979 cuts1 = edgeIntersectionsAndShuffleCGAL
989 cuts2 = edgeIntersectionsAndShuffleCGAL
1007void calcEdgeCutsBitsCGAL
1023 Info<<
"Surface triangles " << surf1.
size()
1024 <<
" number of zones (orientation or non-manifold):"
1035 Info<<
"Surface triangles " << surf2.
size()
1036 <<
" number of zones (orientation or non-manifold):"
1041 if (nZones1 == 1 && nZones2 == 1)
1056 List<List<pointIndexHit>>(surf1.
nEdges()),
1061 List<List<pointIndexHit>>(surf2.
nEdges()),
1066 for (label zone1I = 0; zone1I < nZones1; zone1I++)
1073 if (zone1[faceI] == zone1I)
1075 includeMap1[faceI] =
true;
1094 dupNonManifoldPoints(subSurf1, pointMap1);
1099 subSurf1.meshPoints(),
1114 for (label zone2I = 0; zone2I < nZones2; zone2I++)
1121 if (zone2[faceI] == zone2I)
1123 includeMap2[faceI] =
true;
1143 subSurf2.meshPoints(),
1148 if (!subBb2.overlaps(subBb1))
1155 dupNonManifoldPoints(subSurf2, pointMap2);
1185 label subMeshPointI = subSurf2.meshPoints()[i];
1186 label meshPointI = surf2.
meshPoints()[pointMap2[i]];
1187 points2[meshPointI] = subSurf2.points()[subMeshPointI];
1193 edgeCuts1.
merge(subEdgeCuts1, edgeMap1, faceMap2);
1194 edgeCuts2.
merge(subEdgeCuts2, edgeMap2, faceMap1);
1203 label subMeshPointI = subSurf1.meshPoints()[i];
1204 label meshPointI = surf1.
meshPoints()[pointMap1[i]];
1205 points1[meshPointI] = subSurf1.points()[subMeshPointI];
1270 const bool surf1Baffle,
1271 const bool surf2Baffle,
1272 const bool invertedSpace,
1288 label nFeatEds = inter.cutEdges().size();
1296 List<DynamicList<label>> edgeNormals(nFeatEds);
1297 List<DynamicList<label>> normalDirections(nFeatEds);
1306 const label cutEdgeI = iter.val();
1308 const edge& fE = inter.cutEdges()[cutEdgeI];
1316 edgeDirections[cutEdgeI] = fE.
vec(inter.cutPoints());
1318 normals.append(norm1);
1319 eNormals.
append(normals.size() - 1);
1336 edgeDirections[cutEdgeI],
1338 inter.cutPoints()[fE.
start()]
1343 normals.append(norm2);
1344 eNormals.
append(normals.size() - 1);
1362 edgeDirections[cutEdgeI],
1364 inter.cutPoints()[fE.
start()]
1372 normals.append(norm2);
1390 edgeDirections[cutEdgeI],
1392 inter.cutPoints()[fE.
start()]
1397 eNormals.
append(normals.size() - 1);
1402 normals.append(norm1);
1420 edgeDirections[cutEdgeI],
1422 inter.cutPoints()[fE.
start()]
1427 eNormals.
append(normals.size() - 1);
1432 label internalStart = -1;
1433 label nIntOrExt = 0;
1440 label nEdNorms = edgeNormals[eI].size();
1446 else if (nEdNorms == 2)
1448 const vector& n0(normals[edgeNormals[eI][0]]);
1449 const vector& n1(normals[edgeNormals[eI][1]]);
1476 internalStart = nIntOrExt;
1484 internalStart = nIntOrExt;
1495 internalStart = nIntOrExt;
1500 <<
"Unsupported booleanSurface:booleanOpType and space "
1501 << action <<
" " << invertedSpace
1510 List<extendedFeatureEdgeMesh::sideVolumeType> normalVolumeTypesTmp
1515 forAll(edgeNormalsTmp, i)
1517 edgeNormalsTmp[i] = edgeNormals[i];
1520 forAll(normalDirectionsTmp, i)
1522 normalDirectionsTmp[i] = normalDirections[i];
1540 nIntOrExt + nFlat + nOpen,
1543 normalVolumeTypesTmp,
1545 normalDirectionsTmp,
1554int main(
int argc,
char *argv[])
1558 "Generates a extendedFeatureEdgeMesh for the interface created by"
1559 " a boolean operation on two surfaces."
1561 " [Compiled with CGAL]"
1563 " [Compiled without CGAL]"
1571 "One of (intersection | union | difference)"
1580 "Geometry scaling factor (both surfaces)"
1585 "Mark surface 1 as a baffle"
1591 "Mark surface 2 as a baffle"
1597 "Perturb surface points to escape degenerate intersections"
1604 "Do not use CGAL algorithms"
1606 "Ignored, compiled without CGAL"
1613 "Do the surfaces have inverted space orientation, "
1614 "i.e. a point at infinity is considered inside. "
1615 "This is only sensible for union and intersection."
1621 "((surface1 volumeType) .. (surfaceN volumeType))",
1622 "Trim resulting intersection with additional surfaces;"
1623 " volumeType is 'inside' (keep (parts of) edges that are inside)"
1624 ", 'outside' (keep (parts of) edges that are outside) or"
1625 " 'mixed' (keep all)"
1640 if (!validActions.
found(action))
1643 <<
"Unsupported action " << action <<
endl
1644 <<
"Supported actions:" << validActions <<
nl
1649 List<Pair<word>> surfaceAndSide;
1650 if (
args.readIfPresent(
"trim", surfaceAndSide))
1652 Info<<
"Trimming edges with " << surfaceAndSide <<
endl;
1657 const scalar scaleFactor =
args.getOrDefault<scalar>(
"scale", -1);
1660 Info<<
"Reading surface " << surf1Name <<
endl;
1671 if (scaleFactor > 0)
1673 Info<<
"Scaling : " << scaleFactor <<
nl;
1677 Info<< surf1Name <<
" statistics:" <<
endl;
1682 Info<<
"Reading surface " << surf2Name <<
endl;
1693 if (scaleFactor > 0)
1695 Info<<
"Scaling : " << scaleFactor <<
nl;
1699 Info<< surf2Name <<
" statistics:" <<
endl;
1703 const bool surf1Baffle =
args.found(
"surf1Baffle");
1704 const bool surf2Baffle =
args.found(
"surf2Baffle");
1709 const bool invertedSpace =
args.found(
"invertedSpace");
1714 <<
"Inverted space only makes sense for union or intersection."
1721 if (!
args.found(
"no-cgal"))
1723 calcEdgeCutsBitsCGAL
1727 args.found(
"perturb"),
1739 args.found(
"perturb"),
1760 sFeatFileName +
".extendedFeatureEdgeMesh",
1762 "extendedFeatureEdgeMesh",
1780 forAll(surfaceAndSide, trimI)
1782 const word& trimName = surfaceAndSide[trimI].
first();
1788 Info<<
"Reading trim surface " << trimName <<
endl;
1802 Info<< trimName <<
" statistics:" <<
endl;
1803 trimSurf.writeStats(
Info);
1821 Info<<
nl <<
"Writing extendedFeatureEdgeMesh to "
1833 sFeatFileName +
".eMesh",
1845 Info<<
nl <<
"Writing featureEdgeMesh to "
1846 << bfeMesh.objectRelPath() <<
endl;
1848 bfeMesh.regIOobject::write();
CGAL data structures used for triSurface handling.
CGAL::Segment_3< K > Segment
CGAL::Polyhedron_3< K, My_items > Polyhedron
Dynamically sized Field. Similar to DynamicList, but inheriting from a Field instead of a List.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
void append(const T &val)
Copy append an element to the end of this list.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
bool found(const word &key) const
Same as contains().
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
fileName objectRelPath() const
The object path relative to the case.
fileName path() const
The complete path for the object (with instance, local,...).
A HashTable to objects of type <T> with a label key.
const T & first() const noexcept
Access the first element.
const T & second() const noexcept
Access the second element.
label nEdges() const
Number of edges in patch.
label nPoints() const
Number of points supporting patch faces.
const labelListList & pointEdges() const
Return point-edge addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
const Field< point_type > & faceNormals() const
Return face unit normals for patch.
const Field< point_type > & points() const noexcept
Return reference to global points.
T & first()
Access first element of the list, position [0].
void size(const label n)
Older name for setAddressableSize.
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
static void addArgument(const string &argName, const string &usage="")
Append a (mandatory) argument to validArgs.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
static void noParallel()
Remove the parallel options.
static void addOption(const word &optName, const string ¶m="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
static void addNote(const string ¬e)
Add extra notes for the usage information.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
static const Enum< booleanOpType > booleanOpTypeNames
booleanOpType
Enumeration listing the possible volume operator types.
A bounding box defined in terms of min/max extrema points.
Holder of intersections of edges of a surface with another surface. Optionally shuffles around points...
const labelListList & classification() const
For every intersection the classification status.
void merge(const edgeIntersections &, const labelList &edgeMap, const labelList &faceMap, const bool merge=true)
Merge (or override) edge intersection for a subset.
static scalarField minEdgeLength(const triSurface &surf)
Calculate min edge length for every surface point.
label removeDegenerates(const label nIters, const triSurface &surf1, const triSurfaceSearch &query2, const scalarField &surf1PointTol, pointField &points1)
Resolve ties. Shuffles points so all edge - face intersections.
const pointField & points() const noexcept
Return points.
const edgeList & edges() const noexcept
Return edges.
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
vector vec(const UList< point > &pts) const
Return the vector (from first to second).
linePointRef line(const UList< point > &pts) const
Return edge line.
label start() const noexcept
The start (first) vertex label.
void writeObj(const fileName &prefix) const
Write all components of the extendedEdgeMesh as obj files.
virtual void writeStats(Ostream &os) const
Dump some information.
static scalar cosNormalAngleTol_
Angular closeness tolerance for treating normals as the same.
virtual bool write(const bool writeOnProc=true) const
Give precedence to the regIOobject write.
A class for handling file names.
static std::string stem(const std::string &str)
Return the basename, without extension.
A triFace with additional (region) index.
Basic surface-surface intersection description. Constructed from two surfaces it creates a descriptio...
IOoject and searching on triSurface.
static word meshSubDir
Return the mesh sub-directory name (usually "triSurface").
Helper class to search on triSurface.
Triangulated surface description with patch information.
triSurface subsetMesh(const UList< bool > &include, labelList &pointMap, labelList &faceMap) const
Return a new surface subsetted on the selected faces.
void write(Ostream &os) const
Write to Ostream in simple OpenFOAM format.
void writeStats(Ostream &os) const
Write some statistics.
virtual void movePoints(const pointField &pts)
Move points.
virtual void scalePoints(const scalar scaleFactor)
Scale points. A non-positive factor is ignored.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
An enumeration wrapper for classification of a location as being inside/outside of a volume.
static const Enum< volumeType::type > names
Names for the classification enumeration.
A class for handling words, derived from Foam::string.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
const dimensionedScalar mp
Proton mass.
List< edge > edgeList
List of edge.
Pair< label > labelPair
A pair of labels.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
List< labelList > labelListList
List of labelList.
List< label > labelList
A List of labels.
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
UIndirectList< label > labelUIndList
UIndirectList of labels.
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
List< bool > boolList
A List of bools.
static constexpr const zero Zero
Global zero (0).
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
void cross(FieldField< Field1, typename crossProduct< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
vectorField pointField
pointField is a vectorField.
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
constexpr char nl
The newline '\n' character (0x0a).
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} pointMap[start]=pointMap[end]=Foam::min(start, end);} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};constexpr label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< DynamicList< face > > pFaces[nBCs]
Tree tree(triangles.begin(), triangles.end())
Foam::argList args(argc, argv)
#define forAll(list, i)
Loop across all elements in list.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.