Go to the documentation of this file.
37 #include "triSurface.H"
60 if (elems[elemI] == elem)
82 boolList cutFace(mesh_.nFaces(),
false);
89 forAll(mesh_.edges(), edgeI)
91 if (debug && (edgeI % 10000 == 0))
93 Pout<<
"Intersecting mesh edge " << edgeI <<
" with surface"
97 const edge&
e = mesh_.edges()[edgeI];
99 const point& p0 = mesh_.points()[
e.start()];
100 const point& p1 = mesh_.points()[
e.end()];
106 const labelList& myFaces = mesh_.edgeFaces()[edgeI];
110 label faceI = myFaces[myFaceI];
114 cutFace[faceI] =
true;
124 Pout<<
"Intersected edges of mesh with surface in = "
132 labelList allFaces(mesh_.nFaces() - nCutFaces);
140 allFaces[allFaceI++] = faceI;
146 Pout<<
"Testing " << allFaceI <<
" faces for piercing by surface"
153 scalar tol = 1
e-6 * allBb.
avgDim();
182 if (debug && (edgeI % 10000 == 0))
184 Pout<<
"Intersecting surface edge " << edgeI
185 <<
" with mesh faces" <<
endl;
187 const edge&
e = edges[edgeI];
189 const point& start = localPoints[
e.start()];
190 const point& end = localPoints[
e.end()];
192 vector edgeNormal(end - start);
193 const scalar edgeMag =
mag(edgeNormal);
194 const vector smallVec = 1
e-9*edgeNormal;
196 edgeNormal /= edgeMag+VSMALL;
215 cutFace[faceI] =
true;
223 if (((pt-start) & edgeNormal) >= edgeMag)
233 Pout<<
"Detected an additional " << nAddFaces <<
" faces cut"
236 Pout<<
"Intersected edges of surface with mesh faces in = "
261 forAll(piercedFace, faceI)
263 if (piercedFace[faceI])
265 cellInfoList[mesh_.faceOwner()[faceI]] =
268 if (mesh_.isInternalFace(faceI))
270 cellInfoList[mesh_.faceNeighbour()[faceI]] =
281 labelHashSet outsideFacesMap(outsidePts.size() * 6 * 2);
283 forAll(outsidePts, outsidePtI)
286 label cellI = queryMesh.
findCell(outsidePts[outsidePtI], -1,
false);
291 <<
"outsidePoint " << outsidePts[outsidePtI]
292 <<
" is not inside any cell"
293 <<
nl <<
"It might be on a face or outside the geometry"
299 cellInfoList[cellI] =
cellInfo(cellClassification::OUTSIDE);
302 const labelList& myFaces = mesh_.cells()[cellI];
305 outsideFacesMap.
insert(myFaces[myFaceI]);
320 cellInfo(cellClassification::OUTSIDE)
329 mesh_.globalData().nTotalCells()+1
337 label t = allInfo[cellI].type();
339 if (t == cellClassification::NOTSET)
341 t = cellClassification::INSIDE;
343 operator[](cellI) = t;
350 const label meshType,
355 pointSide.
setSize(mesh_.nPoints());
357 forAll(mesh_.pointCells(), pointI)
359 const labelList& pCells = mesh_.pointCells()[pointI];
361 pointSide[pointI] = UNSET;
367 if (
type == meshType)
369 if (pointSide[pointI] == UNSET)
371 pointSide[pointI] =
MESH;
373 else if (pointSide[pointI] == NONMESH)
375 pointSide[pointI] = MIXED;
382 if (pointSide[pointI] == UNSET)
384 pointSide[pointI] = NONMESH;
386 else if (pointSide[pointI] ==
MESH)
388 pointSide[pointI] = MIXED;
404 const faceList& faces = mesh_.faces();
406 const cell& cFaces = mesh_.cells()[cellI];
410 const face&
f = faces[cFaces[cFaceI]];
414 if (pointSide[
f[fp]] != MIXED)
428 const label meshType,
433 const labelList& own = mesh_.faceOwner();
434 const labelList& nbr = mesh_.faceNeighbour();
435 const faceList& faces = mesh_.faces();
437 outsideFaces.
setSize(mesh_.nFaces());
438 outsideOwner.
setSize(mesh_.nFaces());
443 for (
label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
445 label ownType = operator[](own[faceI]);
446 label nbrType = operator[](nbr[faceI]);
448 if (ownType == meshType && nbrType != meshType)
450 outsideFaces[outsideI] = faces[faceI];
451 outsideOwner[outsideI] = own[faceI];
454 else if (ownType != meshType && nbrType == meshType)
456 outsideFaces[outsideI] = faces[faceI];
457 outsideOwner[outsideI] = nbr[faceI];
464 for (
label faceI = mesh_.nInternalFaces(); faceI < mesh_.nFaces(); faceI++)
466 if (
operator[](own[faceI]) == meshType)
468 outsideFaces[outsideI] = faces[faceI];
469 outsideOwner[outsideI] = own[faceI];
473 outsideFaces.
setSize(outsideI);
474 outsideOwner.
setSize(outsideI);
495 markFaces(surfQuery),
511 if (mesh_.nCells() != size())
514 <<
"Number of elements of cellType argument is not equal to the"
536 const label meshType,
562 for (
label iter = 0; iter < nLayers; iter++)
568 classifyPoints(meshType, newCellType, pointSide);
573 if (pointSide[pointI] == MIXED)
576 const labelList& pCells = mesh_.pointCells()[pointI];
586 newCellType[pCells[i]] = meshType;
600 forAll(newCellType, cellI)
604 if (newCellType[cellI] != meshType)
608 operator[](cellI) = fillType;
621 const label meshType,
625 boolList hasMeshType(mesh_.nPoints(),
false);
629 forAll(mesh_.pointCells(), pointI)
631 const labelList& myCells = mesh_.pointCells()[pointI];
636 label type = operator[](myCells[myCellI]);
638 if (
type == meshType)
640 hasMeshType[pointI] =
true;
651 forAll(hasMeshType, pointI)
653 if (hasMeshType[pointI])
655 const labelList& myCells = mesh_.pointCells()[pointI];
659 if (
operator[](myCells[myCellI]) != meshType)
661 operator[](myCells[myCellI]) = fillType;
678 const label meshType,
679 const label fillType,
683 label nTotChanged = 0;
685 for (
label iter = 0; iter < maxIter; iter++)
691 classifyPoints(meshType, *
this, pointSide);
698 if (pointSide[pointI] == MIXED)
700 const labelList& pCells = mesh_.pointCells()[pointI];
704 label cellI = pCells[i];
706 if (
operator[](cellI) == meshType)
708 if (usesMixedPointsOnly(pointSide, cellI))
710 operator[](cellI) = fillType;
718 nTotChanged += nChanged;
720 Pout<<
"removeHangingCells : changed " << nChanged
721 <<
" hanging cells" <<
endl;
735 const label meshType,
736 const label fillType,
740 label nTotChanged = 0;
742 for (
label iter = 0; iter < maxIter; iter++)
749 getMeshOutside(meshType, outsideFaces, outsideOwner);
762 const labelList& eFaces = edgeFaces[edgeI];
764 if (eFaces.
size() > 2)
771 label patchFaceI = eFaces[i];
773 label ownerCell = outsideOwner[patchFaceI];
775 if (
operator[](ownerCell) == meshType)
777 operator[](ownerCell) = fillType;
787 nTotChanged += nChanged;
789 Pout<<
"fillRegionEdges : changed " << nChanged
790 <<
" cells using multiply connected edges" <<
endl;
804 const label meshType,
805 const label fillType,
809 label nTotChanged = 0;
811 for (
label iter = 0; iter < maxIter; iter++)
818 getMeshOutside(meshType, outsideFaces, outsideOwner);
835 const label patchPointI = meshPointMap[iter.key()];
845 const label ownerCell = outsideOwner[patchFaceI];
847 if (
operator[](ownerCell) == meshType)
849 operator[](ownerCell) = fillType;
857 nTotChanged += nChanged;
859 Pout<<
"fillRegionPoints : changed " << nChanged
860 <<
" cells using multiply connected points" <<
endl;
874 os <<
"Cells:" << size() <<
endl
875 <<
" notset : " << count(*
this, NOTSET) <<
endl
876 <<
" cut : " << count(*
this, CUT) <<
endl
877 <<
" inside : " << count(*
this, INSIDE) <<
endl
878 <<
" outside : " << count(*
this, OUTSIDE) <<
endl;
Starts timing CPU usage and return elapsed time from start.
const labelListList & pointFaces() const
Return point-face addressing.
label index() const
Return index.
const labelListList & edgeFaces() const
Return edge-face addressing.
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
List< label > labelList
A List of labels.
const point & max() const
Maximum describing the bounding box.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
List< Key > toc() const
Return the table of contents.
#define forAll(list, i)
Loop across all elements in list.
bool usesMixedPointsOnly(const List< pointStatus > &, const label cellI) const
Return true if cell uses only points with status=mixed.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Standard boundBox + extra functionality for use in octree.
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
bool hit() const
Is there a hit.
const Type & shapes() const
Reference to shape.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
boolList markFaces(const triSurfaceSearch &) const
Mark all faces intersected by or intersecting surface.
void operator=(const UList< T > &)
Assignment from UList operator. Takes linear time.
Helper class to search on triSurface.
Mesh consisting of general polyhedral cells.
label trimCutCells(const label nLayers, const label meshType, const label fillType)
label fillHangingCells(const label meshType, const label fillType, const label maxIter)
Find hanging cells (cells with all points on outside) and set their.
static label count(const labelList &elems, const label elem)
Count number of occurrences of elem in list.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("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]=cellShape(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];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } 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};const 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< SLList< face > > pFaces[nBCs]
void classifyPoints(const label meshType, const labelList &cellType, List< pointStatus > &pointSide) const
Use cell status to classify points as being internal to meshType,.
This class describes the interaction of (usually) a face and a point. It carries the info of a succes...
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
cellClassification(const polyMesh &mesh, const meshSearch &meshQuery, const triSurfaceSearch &surfQuery, const pointField &outsidePoints)
Construct from mesh and surface and point(s) on outside.
Triangulated surface description with patch information.
pointIndexHit findLine(const bool findAny, const point &treeStart, const point &treeEnd, const label startNodeI, const direction startOctantI, const FindIntersectOp &fiOp, const bool verbose=false) const
Find any or nearest intersection.
void operator=(const cellClassification &)
const Point & hitPoint() const
Return hit point.
const indexedOctree< treeDataTriSurface > & tree() const
Demand driven construction of the octree.
Non-pointer based hierarchical recursive searching.
label fillRegionEdges(const label meshType, const label fillType, const label maxIter)
Find regionEdges and fill one neighbour. Iterate until nothing.
const point & min() const
Minimum describing the bounding box.
void markCells(const meshSearch &queryMesh, const boolList &piercedFace, const pointField &outsidePts)
Divide cells into cut/inside/outside by using MeshWave from cut.
errorManip< error > abort(error &err)
'Cuts' a mesh with a surface.
void writeStats(Ostream &os) const
Write statistics on cell types to Ostream.
errorManipArg< error, int > exit(error &err, const int errNo=1)
void setSize(const label)
Reset size of List.
Holds information regarding type of cell. Used in inside/outside determination in cellClassification.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
prefixOSstream Pout(cout, "Pout")
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
label findCell(const point &location, const label seedCellI=-1, const bool useTreeSearch=true) const
Find cell containing location.
Encapsulation of data needed to search for faces.
label fillRegionPoints(const label meshType, const label fillType, const label maxIter)
Find regionPoints and fill all neighbours. Iterate until nothing.
Implements a timeout mechanism via sigalarm.
const triSurface & surface() const
Return reference to the surface.
const dimensionedScalar e
Elementary charge.
scalar avgDim() const
Average length/height/width dimension.
bool insert(const Key &key)
Insert a new entry.
A face is a list of labels corresponding to mesh vertices.
Various functions to operate on Lists.
void size(const label)
Override size to be inconsistent with allocated storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
label growSurface(const label meshType, const label fillType)
Sets vertex neighbours of meshType cells to fillType.
const Map< label > & meshPointMap() const
Mesh point map. Given the global point index find its.
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
bool checkPointManifold(const bool report=false, labelHashSet *setPtr=NULL) const
Checks primitivePatch for faces sharing point but not edge.
defineTypeNameAndDebug(combustionModel, 0)
A cell is defined as a list of faces with extra functionality.
const List< Type > & allCellInfo() const
Get allCellInfo.
A list of faces which address into the list of points.
void getMeshOutside(const label meshType, faceList &, labelList &) const
Get faces (and its 'owner') inbetween cells of differing type.