Cell-face mesh analysis engine. More...
#include <primitiveMesh.H>

Public Member Functions | |
| ClassName ("primitiveMesh") | |
| primitiveMesh (const label nPoints, const label nInternalFaces, const label nFaces, const label nCells) | |
| Construct from components. | |
| virtual | ~primitiveMesh () |
| Destructor. | |
| void | reset (const label nPoints, const label nInternalFaces, const label nFaces, const label nCells) |
| Reset this primitiveMesh given the primitive array sizes. | |
| void | reset (const label nPoints, const label nInternalFaces, const label nFaces, const label nCells, cellList &cells) |
| Reset this primitiveMesh given the primitive array sizes and cells. | |
| void | resetGeometry (pointField &&faceCentres, pointField &&faceAreas, pointField &&cellCentres, scalarField &&cellVolumes) |
| Reset the local geometry. | |
| virtual bool | init (const bool doInit) |
| Initialise all non-demand-driven data. | |
| label | nPoints () const noexcept |
| Number of mesh points. | |
| label | nEdges () const |
| Number of mesh edges. | |
| label | nFaces () const noexcept |
| Number of mesh faces. | |
| label | nCells () const noexcept |
| Number of mesh cells. | |
| label | nInternalFaces () const noexcept |
| Number of internal faces. | |
| label | nBoundaryFaces () const noexcept |
| Number of boundary faces (== nFaces - nInternalFaces). | |
| label | nInternalPoints () const noexcept |
| Points not on boundary. | |
| label | nInternal0Edges () const |
| Internal edges (i.e. not on boundary face) using no boundary point. | |
| label | nInternal1Edges () const |
| Internal edges using 0 or 1 boundary point. | |
| label | nInternalEdges () const |
| Internal edges using 0,1 or 2 boundary points. | |
| virtual const pointField & | points () const =0 |
| Return mesh points. | |
| virtual const faceList & | faces () const =0 |
| Return faces. | |
| virtual const labelList & | faceOwner () const =0 |
| Face face-owner addressing. | |
| virtual const labelList & | faceNeighbour () const =0 |
| Face face-neighbour addressing. | |
| virtual const pointField & | oldPoints () const =0 |
| Return old points for mesh motion. | |
| const cellShapeList & | cellShapes () const |
| Return cell shapes. | |
| const edgeList & | edges () const |
| Return mesh edges. Uses calcEdges. | |
| const labelListList & | cellCells () const |
| const labelListList & | edgeCells () const |
| const labelListList & | pointCells () const |
| const cellList & | cells () const |
| const labelListList & | edgeFaces () const |
| const labelListList & | pointFaces () const |
| const labelListList & | cellEdges () const |
| const labelListList & | faceEdges () const |
| const labelListList & | pointEdges () const |
| const labelListList & | pointPoints () const |
| const labelListList & | cellPoints () const |
| const vectorField & | cellCentres () const |
| const vectorField & | faceCentres () const |
| const scalarField & | cellVolumes () const |
| const vectorField & | faceAreas () const |
| void | movePoints (const pointField &p, const pointField &oldP) |
| Move points. | |
| bool | isInternalFace (const label faceIndex) const noexcept |
| Return true if given face label is internal to the mesh. | |
| virtual bool | checkUpperTriangular (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check face ordering. | |
| virtual bool | checkCellsZipUp (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check cell zip-up. | |
| virtual bool | checkFaceVertices (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check uniqueness of face vertices. | |
| virtual bool | checkPoints (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check for unused points. | |
| virtual bool | checkFaceFaces (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check face-face connectivity. | |
| virtual bool | checkClosedBoundary (const bool report=false) const |
| Check boundary for closedness. | |
| virtual bool | checkClosedCells (const bool report=false, labelHashSet *setPtr=nullptr, labelHashSet *highAspectSetPtr=nullptr, const Vector< label > &solutionD=Vector< label >::one) const |
| Check cells for closedness. | |
| virtual bool | checkFaceAreas (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check for negative face areas. | |
| virtual bool | checkCellVolumes (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check for negative cell volumes. | |
| virtual bool | checkFaceOrthogonality (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check for non-orthogonality. | |
| virtual bool | checkFacePyramids (const bool report=false, const scalar minPyrVol=-SMALL, labelHashSet *setPtr=nullptr) const |
| Check face pyramid volume. | |
| virtual bool | checkFaceSkewness (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check face skewness. | |
| virtual bool | checkFaceAngles (const bool report=false, const scalar maxSin=10, labelHashSet *setPtr=nullptr) const |
| Check face angles. | |
| virtual bool | checkFaceFlatness (const bool report, const scalar warnFlatness, labelHashSet *setPtr) const |
| Check face warpage: decompose face and check ratio between. | |
| virtual bool | checkPointNearness (const bool report, const scalar reportDistSqr, labelHashSet *setPtr=nullptr) const |
| Check for point-point-nearness,. | |
| virtual bool | checkEdgeLength (const bool report, const scalar minLenSqr, labelHashSet *setPtr=nullptr) const |
| Check edge length. | |
| virtual bool | checkConcaveCells (const bool report=false, labelHashSet *setPtr=nullptr) const |
| Check for concave cells by the planes of faces. | |
| virtual bool | checkTopology (const bool report=false) const |
| Check mesh topology for correctness. | |
| virtual bool | checkGeometry (const bool report=false) const |
| Check mesh geometry (& implicitly topology) for correctness. | |
| virtual bool | checkMesh (const bool report=false) const |
| Check mesh for correctness. Returns false for no error. | |
| boundBox | cellBb (const label celli) const |
| The bounding box for given cell index. | |
| bool | pointInCellBB (const point &p, label celli, scalar inflationFraction=0) const |
| Return true if the point in the cell bounding box. | |
| bool | pointInCell (const point &p, label celli) const |
| Return true if the point is in the cell. | |
| label | findNearestCell (const point &location) const |
| Find the cell with the nearest cell centre to location. | |
| label | findCell (const point &location) const |
| Find cell enclosing this location (-1 if not in mesh). | |
| void | printAllocated () const |
| Print a list of all the currently allocated mesh data. | |
| bool | hasCellShapes () const noexcept |
| bool | hasEdges () const noexcept |
| bool | hasCellCells () const noexcept |
| bool | hasEdgeCells () const noexcept |
| bool | hasPointCells () const noexcept |
| bool | hasCells () const noexcept |
| bool | hasEdgeFaces () const noexcept |
| bool | hasPointFaces () const noexcept |
| bool | hasCellEdges () const noexcept |
| bool | hasFaceEdges () const noexcept |
| bool | hasPointEdges () const noexcept |
| bool | hasPointPoints () const noexcept |
| bool | hasCellPoints () const noexcept |
| bool | hasCellCentres () const noexcept |
| bool | hasCellVolumes () const noexcept |
| bool | hasFaceCentres () const noexcept |
| bool | hasFaceAreas () const noexcept |
| const labelList & | cellCells (const label celli, DynamicList< label > &) const |
| cellCells using cells. | |
| const labelList & | cellCells (const label celli) const |
| const labelList & | cellPoints (const label celli, labelHashSet &, DynamicList< label > &) const |
| cellPoints using cells | |
| const labelList & | cellPoints (const label celli) const |
| const labelList & | pointCells (const label pointi, DynamicList< label > &) const |
| pointCells using pointFaces | |
| const labelList & | pointCells (const label pointi) const |
| const labelList & | pointPoints (const label pointi, DynamicList< label > &) const |
| pointPoints using edges, pointEdges | |
| const labelList & | pointPoints (const label pointi) const |
| const labelList & | faceEdges (const label facei, DynamicList< label > &) const |
| faceEdges using pointFaces, edges, pointEdges | |
| const labelList & | faceEdges (const label facei) const |
| const labelList & | edgeFaces (const label edgeI, DynamicList< label > &) const |
| edgeFaces using pointFaces, edges, pointEdges | |
| const labelList & | edgeFaces (const label edgeI) const |
| const labelList & | edgeCells (const label edgeI, DynamicList< label > &) const |
| edgeCells using pointFaces, edges, pointEdges | |
| const labelList & | edgeCells (const label edgeI) const |
| const labelList & | cellEdges (const label celli, labelHashSet &, DynamicList< label > &) const |
| cellEdges using cells, pointFaces, edges, pointEdges | |
| const labelList & | cellEdges (const label celli) const |
| virtual void | updateGeom () |
| Update all geometric data. | |
| void | clearGeom () |
| Clear geometry. | |
| void | clearCellGeom () |
| Clear cell-based geometry only. | |
| void | clearAddressing () |
| Clear topological data. | |
| void | clearOut () |
| Clear all geometry and addressing unnecessary for CFD. | |
Static Public Member Functions | |
| static void | calcCells (cellList &, const labelUList &own, const labelUList &nei, const label nCells=-1) |
| Helper function to calculate cell-face addressing from. | |
| static bool | calcPointOrder (label &nInternalPoints, labelList &pointMap, const faceList &, const label nInternalFaces, const label nPoints) |
| Helper function to calculate point ordering. Returns true. | |
| static scalar | setClosedThreshold (const scalar) |
| Set the closedness ratio warning threshold. | |
| static scalar | setAspectThreshold (const scalar) |
| Set the aspect ratio warning threshold. | |
| static scalar | setNonOrthThreshold (const scalar) |
| Set the non-orthogonality warning threshold in degrees. | |
| static scalar | setSkewThreshold (const scalar) |
| Set the skewness warning threshold as percentage. | |
Static Public Attributes | |
| static const unsigned | cellsPerEdge_ = 4 |
| Estimated number of cells per edge. | |
| static const unsigned | cellsPerPoint_ = 8 |
| Estimated number of cells per point. | |
| static const unsigned | facesPerCell_ = 6 |
| Estimated number of faces per cell. | |
| static const unsigned | facesPerEdge_ = 4 |
| Estimated number of faces per edge. | |
| static const unsigned | facesPerPoint_ = 12 |
| Estimated number of faces per point. | |
| static const unsigned | edgesPerCell_ = 12 |
| Estimated number of edges per cell. | |
| static const unsigned | edgesPerFace_ = 4 |
| Estimated number of edges per cell. | |
| static const unsigned | edgesPerPoint_ = 6 |
| Estimated number of edges per point. | |
| static const unsigned | pointsPerCell_ = 8 |
| Estimated number of points per cell. | |
| static const unsigned | pointsPerFace_ = 4 |
| Estimated number of points per face. | |
Protected Member Functions | |
| void | calcFaceCentresAndAreas () const |
| Calculate face centres and areas. | |
| void | calcCellCentresAndVols () const |
| Calculate cell centres and volumes. | |
| void | calcEdgeVectors () const |
| Calculate edge vectors. | |
| bool | checkDuplicateFaces (const label, const Map< label > &, label &nBaffleFaces, labelHashSet *) const |
| Check if all points on face are shared with another face. | |
| bool | checkCommonOrder (const label, const Map< label > &, labelHashSet *) const |
| Check that shared points are in consecutive order. | |
| bool | checkClosedBoundary (const vectorField &areas, const bool report, const bitSet &internalOrCoupledFaces) const |
| Check boundary for closedness. | |
| bool | checkClosedCells (const vectorField &faceAreas, const scalarField &cellVolumes, const bool report, labelHashSet *setPtr, labelHashSet *aspectSetPtr, const Vector< label > &meshD) const |
| Check cells for closedness. | |
| bool | checkFaceAreas (const vectorField &faceAreas, const bool report, const bool detailedReport, labelHashSet *setPtr) const |
| Check for negative face areas. | |
| bool | checkCellVolumes (const scalarField &vols, const bool report, const bool detailedReport, labelHashSet *setPtr) const |
| Check for negative cell volumes. | |
| bool | checkFaceOrthogonality (const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const |
| Check for non-orthogonality. | |
| bool | checkFacePyramids (const pointField &points, const vectorField &ctrs, const bool report, const bool detailedReport, const scalar minPyrVol, labelHashSet *setPtr) const |
| Check face pyramid volume. | |
| bool | checkFaceSkewness (const pointField &points, const vectorField &fCtrs, const vectorField &fAreas, const vectorField &cellCtrs, const bool report, labelHashSet *setPtr) const |
| Check face skewness. | |
| bool | checkFaceAngles (const pointField &points, const vectorField &faceAreas, const bool report, const scalar maxDeg, labelHashSet *setPtr) const |
| Check face angles. | |
| bool | checkFaceFlatness (const pointField &points, const vectorField &faceCentres, const vectorField &faceAreas, const bool report, const scalar warnFlatness, labelHashSet *setPtr) const |
| Check face warpage. | |
| bool | checkConcaveCells (const vectorField &fAreas, const pointField &fCentres, const bool report, labelHashSet *setPtr) const |
| Check for concave cells by the planes of faces. | |
| primitiveMesh () | |
| Construct null. | |
Static Protected Attributes | |
| static scalar | closedThreshold_ = 1.0e-6 |
| Static data to control mesh checking. | |
| static scalar | aspectThreshold_ = 1000 |
| Aspect ratio warning threshold. | |
| static scalar | nonOrthThreshold_ = 70 |
| Non-orthogonality warning threshold in deg. | |
| static scalar | skewThreshold_ = 4 |
| Skewness warning threshold. | |
| static scalar | planarCosAngle_ = 1.0e-6 |
| Threshold where faces are considered coplanar. | |
Cell-face mesh analysis engine.
Definition at line 75 of file primitiveMesh.H.
|
protected |
Construct null.
Definition at line 34 of file primitiveMesh.C.
| primitiveMesh | ( | const label | nPoints, |
| const label | nInternalFaces, | ||
| const label | nFaces, | ||
| const label | nCells ) |
Construct from components.
Definition at line 48 of file primitiveMesh.C.
References nCells(), nFaces(), nInternalFaces(), and nPoints().

|
virtual |
Destructor.
Definition at line 67 of file primitiveMesh.C.
References clearOut().

|
protected |
Calculate face centres and areas.
Definition at line 35 of file primitiveMeshFaceCentresAndAreas.C.
References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, primitiveMeshTools::makeFaceCentresAndAreas(), nFaces(), points(), and Foam::Pout.
Referenced by updateGeom().


|
protected |
Calculate cell centres and volumes.
Definition at line 32 of file primitiveMeshCellCentresAndVols.C.
References Foam::abort(), cellVols, Foam::endl(), faceAreas(), faceCentres(), Foam::FatalError, FatalErrorInFunction, primitiveMeshTools::makeCellCentresAndVols(), nCells(), Foam::Pout, and Foam::Zero.
Referenced by updateGeom().


|
protected |
Calculate edge vectors.
References cellVolumes(), faceAreas(), faceCentres(), and points().

|
protected |
Check if all points on face are shared with another face.
Definition at line 1293 of file primitiveMeshCheck.C.
References faces(), forAllConstIters, HashSet< Key, Hash >::insert(), and UList< T >::size().
Referenced by checkFaceFaces().


|
protected |
Check that shared points are in consecutive order.
Definition at line 1334 of file primitiveMeshCheck.C.
References faces(), UList< T >::fcIndex(), face::find(), forAll, forAllConstIters, HashSet< Key, Hash >::insert(), Foam::labelMax, UList< T >::rcIndex(), and UList< T >::size().
Referenced by checkFaceFaces().


|
protected |
Check boundary for closedness.
Definition at line 41 of file primitiveMeshCheck.C.
References closedThreshold_, Foam::cmptMag(), Foam::cmptMax(), DebugInFunction, Foam::endl(), Foam::Info, Foam::mag(), nInternalFaces(), Foam::reduce(), PackedList< Width >::size(), UList< T >::size(), and Foam::Zero.
Referenced by checkClosedBoundary(), checkGeometry(), and oldPoints().


|
protected |
Check cells for closedness.
Definition at line 92 of file primitiveMeshCheck.C.
References aspectThreshold_, primitiveMeshTools::cellClosedness(), cells, cellVolumes(), closedThreshold_, DebugInFunction, Foam::endl(), faceAreas(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::max(), Foam::min(), nFaces(), Foam::nl, and Foam::reduce().
Referenced by checkClosedCells(), checkGeometry(), and oldPoints().


|
protected |
Check for negative face areas.
Definition at line 220 of file primitiveMeshCheck.C.
References DebugInFunction, Foam::endl(), faceAreas(), faceNeighbour(), faceOwner(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), isInternalFace(), Foam::mag(), Foam::max(), Foam::min(), Foam::Pout, and Foam::reduce().
Referenced by checkFaceAreas(), checkGeometry(), polyMesh::checkMeshMotion(), and oldPoints().


|
protected |
Check for negative cell volumes.
Definition at line 293 of file primitiveMeshCheck.C.
References DebugInFunction, Foam::endl(), forAll, Foam::gSum(), Foam::Info, HashSet< Key, Hash >::insert(), Foam::max(), Foam::min(), Foam::Pout, and Foam::reduce().
Referenced by checkCellVolumes(), checkGeometry(), polyMesh::checkMeshMotion(), and oldPoints().


|
protected |
Check for non-orthogonality.
Definition at line 358 of file primitiveMeshCheck.C.
References Foam::acos(), Foam::cos(), DebugInFunction, Foam::degToRad(), Foam::endl(), primitiveMeshTools::faceOrthogonality(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::min(), nonOrthThreshold_, Foam::radToDeg(), Foam::reduce(), Foam::returnReduce(), UList< T >::size(), and Foam::sum().
Referenced by checkFaceOrthogonality(), checkGeometry(), and oldPoints().


|
protected |
Check face pyramid volume.
Definition at line 461 of file primitiveMeshCheck.C.
References cells, DebugInFunction, Foam::endl(), f(), faceNeighbour(), faceOwner(), primitiveMeshTools::facePyramidVolume(), faces(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), isInternalFace(), Foam::nl, points(), Foam::Pout, and Foam::reduce().
Referenced by checkFacePyramids(), checkGeometry(), polyMesh::checkMeshMotion(), and oldPoints().


|
protected |
Check face skewness.
Definition at line 558 of file primitiveMeshCheck.C.
References DebugInFunction, Foam::endl(), primitiveMeshTools::faceSkewness(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::max(), points(), Foam::reduce(), and skewThreshold_.
Referenced by checkFaceSkewness(), checkGeometry(), and oldPoints().


|
protected |
Check face angles.
Allows a slight non-convexity. E.g. maxDeg = 10 allows for angles < 190 (or 10 degrees concavity) (if truly concave and points not visible from face centre the face-pyramid check in checkMesh will fail)
Definition at line 626 of file primitiveMeshCheck.C.
References Foam::asin(), DebugInFunction, Foam::degToRad(), Foam::endl(), Foam::exit(), faceAreas(), primitiveMeshTools::faceConcavity(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::max(), Foam::min(), points(), Foam::radToDeg(), Foam::reduce(), and Foam::sin().
Referenced by checkFaceAngles(), and oldPoints().


|
protected |
Check face warpage.
Definition at line 701 of file primitiveMeshCheck.C.
References DebugInFunction, Foam::endl(), Foam::exit(), faceAreas(), faceCentres(), primitiveMeshTools::faceFlatness(), faces(), Foam::FatalError, FatalErrorInFunction, forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::mag(), Foam::min(), points(), and Foam::reduce().
Referenced by checkFaceFlatness(), and oldPoints().


|
protected |
Check for concave cells by the planes of faces.
Definition at line 801 of file primitiveMeshCheck.C.
References cells, DebugInFunction, Foam::endl(), faceOwner(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::mag(), Foam::max(), planarCosAngle_, and Foam::reduce().
Referenced by checkConcaveCells(), and oldPoints().


| ClassName | ( | "primitiveMesh" | ) |
| void reset | ( | const label | nPoints, |
| const label | nInternalFaces, | ||
| const label | nFaces, | ||
| const label | nCells ) |
Reset this primitiveMesh given the primitive array sizes.
Definition at line 150 of file primitiveMesh.C.
References calcPointOrder(), clearOut(), Foam::endl(), faces(), nCells(), nFaces(), nInternalFaces(), nInternalPoints(), nPoints(), and Foam::Pout.
Referenced by distanceSurface::createGeometry(), and hexRef8Data::hexRef8Data().


| void reset | ( | const label | nPoints, |
| const label | nInternalFaces, | ||
| const label | nFaces, | ||
| const label | nCells, | ||
| cellList & | cells ) |
Reset this primitiveMesh given the primitive array sizes and cells.
Definition at line 206 of file primitiveMesh.C.
References nCells(), nFaces(), nInternalFaces(), nPoints(), and reset().

| void resetGeometry | ( | pointField && | faceCentres, |
| pointField && | faceAreas, | ||
| pointField && | cellCentres, | ||
| scalarField && | cellVolumes ) |
Reset the local geometry.
Definition at line 227 of file primitiveMesh.C.
References Foam::abort(), cellCentres(), cellVolumes(), clearGeom(), Foam::endl(), faceAreas(), faceCentres(), Foam::FatalError, FatalErrorInFunction, and Foam::Pout.
Referenced by averageNeighbourFvGeometryScheme::movePoints(), highAspectRatioFvGeometryScheme::movePoints(), solidBodyFvGeometryScheme::movePoints(), and stabilisedFvGeometryScheme::movePoints().


|
inlinevirtual |
Initialise all non-demand-driven data.
Reimplemented in dynamicFvMesh, dynamicMotionSolverFvMesh, dynamicMotionSolverFvMeshAMI, dynamicMotionSolverListFvMesh, dynamicMotionSolverTopoFvMesh, dynamicMultiMotionSolverFvMesh, dynamicRefineFvMesh, fvMesh, interfaceTrackingFvMesh, movingConeTopoFvMesh, and polyMesh.
Definition at line 617 of file primitiveMesh.H.
Referenced by polyMesh::init().

|
inlinenoexcept |
Number of mesh points.
Definition at line 30 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by polyMeshAdder::add(), polyMeshAdder::add(), calcPointOrder(), checkEdgeLength(), motionSmootherAlgo::checkMesh(), distanceSurface::createGeometry(), Foam::createReconstructMap(), displacementInterpolationMotionSolver::curPoints(), polyMeshFilter::filterEdges(), oldPoints(), pointFaces(), primitiveMesh(), reset(), reset(), singleCellFvMesh::singleCellFvMesh(), singleCellFvMesh::singleCellFvMesh(), rigidBodyMeshMotion::solve(), rigidBodyMeshMotionSolver::solve(), sixDoFRigidBodyMotionSolver::solve(), inverseDistance::update(), and trackingInverseDistance::update().

|
inline |
Number of mesh edges.
Definition at line 60 of file primitiveMeshI.H.
References edges().
Referenced by edgeCells(), edgeFaces(), nInternal0Edges(), nInternal1Edges(), nInternalEdges(), and faBoundaryMesh::whichPatch().


|
inlinenoexcept |
Number of mesh faces.
Definition at line 83 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by polyMeshAdder::add(), polyMeshAdder::add(), snappyLayerDriver::addLayers(), calcFaceCentresAndAreas(), dynamicRefineFvMesh::calculateProtectedCells(), checkClosedCells(), checkFaceFaces(), polyMesh::checkMeshMotion(), Foam::createReconstructMap(), dynamicRefineFvMesh::extendMarkedCells(), interfaceTrackingFvMesh::freeSurfacePressureJump(), interfaceTrackingFvMesh::freeSurfaceSnGradU(), interfaceTrackingFvMesh::freeSurfaceSnGradUn(), dynamicRefineFvMesh::init(), fvMesh::lduAddr(), dynamicRefineFvMesh::mapFields(), fvMesh::mapFields(), dynamicRefineFvMesh::mapNewInternalFaces(), mappedPatchFieldBase< Type >::mappedField(), polyMesh::polyMesh(), polyMesh::polyMesh(), primitiveMesh(), reset(), reset(), displacementPointSmoothingMotionSolver::setFacesToMove(), displacementSmartPointSmoothingMotionSolver::setFacesToMove(), singleCellFvMesh::singleCellFvMesh(), singleCellFvMesh::singleCellFvMesh(), rawTopoChangerFvMesh::update(), mappedFixedInternalValueFvPatchField< Type >::updateCoeffs(), mappedVelocityFluxFixedValueFvPatchField::updateCoeffs(), Sampled< Type >::value(), and polyBoundaryMesh::whichPatchFace().

|
inlinenoexcept |
Number of mesh cells.
Definition at line 89 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by polyMeshAdder::add(), polyMeshAdder::add(), calcCellCentresAndVols(), calcCells(), dynamicRefineFvMesh::checkEightAnchorPoints(), polyMesh::checkMeshMotion(), Foam::createReconstructMap(), polyMesh::findCell(), findCell(), dynamicRefineFvMesh::init(), fvMesh::mapFields(), mappedPatchFieldBase< Type >::mappedField(), meshToMesh::mapSrcToTgt(), meshToMesh::mapTgtToSrc(), inverseDistance::markDonors(), trackingInverseDistance::markDonors(), dynamicRefineFvMesh::maxPointField(), oldPoints(), polyDualMesh::polyDualMesh(), polyDualMesh::polyDualMesh(), polyMesh::polyMesh(), primitiveMesh(), dynamicRefineFvMesh::refine(), reset(), reset(), polyMesh::resetPrimitives(), dynamicRefineFvMesh::unrefine(), cellVolumeWeight::update(), dynamicRefineFvMesh::updateTopology(), Sampled< Type >::value(), and meshToMeshMethod::writeConnectivity().

|
inlinenoexcept |
Number of internal faces.
Definition at line 71 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by polyMeshAdder::add(), snappyLayerDriver::addLayers(), calcPointOrder(), dynamicRefineFvMesh::calculateProtectedCells(), cellCells(), pointSmoother::cellQuality(), checkClosedBoundary(), polyBoundaryMesh::checkDefinition(), checkUpperTriangular(), inverseFaceDistanceDiffusivity::correct(), inversePointDistanceDiffusivity::correct(), dynamicRefineFvMesh::extendMarkedCells(), faceCoupleInfo::faceCoupleInfo(), boundaryMesh::getNearest(), dynamicRefineFvMesh::init(), cellToCellStencil::insertFaceCells(), cellToFaceStencil::insertFaceCells(), fvMeshSubset::interpolate(), volPointInterpolationAdjoint::interpolateBoundaryField(), volPointInterpolationAdjoint::interpolateSensitivitiesField(), volPointInterpolationAdjoint::interpolateSensitivitiesField(), volPointInterpolationAdjoint::makeBoundaryWeights(), dynamicRefineFvMesh::mapFields(), dynamicRefineFvMesh::mapNewInternalFaces(), fvMeshAdder::MapSurfaceField(), oldPoints(), boundaryMesh::patchify(), pointCells(), primitiveMesh(), dynamicRefineFvMesh::refine(), reset(), reset(), sampledMeshedSurface::sampleOnFaces(), rawTopoChangerFvMesh::update(), cellToCellStencil::validBoundaryFaces(), cellToFaceStencil::validBoundaryFaces(), and polyBoundaryMesh::whichPatchFace().

|
inlinenoexcept |
Number of boundary faces (== nFaces - nInternalFaces).
Definition at line 77 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by dynamicRefineFvMesh::calculateProtectedCells(), boundaryMesh::getNearest(), volPointInterpolationAdjoint::interpolateSensitivitiesField(), volPointInterpolationAdjoint::interpolateSensitivitiesField(), boundaryMesh::patchify(), polyDualMesh::polyDualMesh(), and polyDualMesh::polyDualMesh().

|
inlinenoexcept |
Points not on boundary.
Definition at line 24 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by calcPointOrder(), oldPoints(), and reset().

|
inline |
Internal edges (i.e. not on boundary face) using no boundary point.
Definition at line 36 of file primitiveMeshI.H.
References nEdges().

|
inline |
Internal edges using 0 or 1 boundary point.
Definition at line 44 of file primitiveMeshI.H.
References nEdges().

|
inline |
Internal edges using 0,1 or 2 boundary points.
Definition at line 52 of file primitiveMeshI.H.
References nEdges().
Referenced by faBoundaryMesh::checkDefinition().


|
pure virtual |
Return mesh points.
Implemented in polyMesh.
References points().
Referenced by calcEdgeVectors(), calcFaceCentresAndAreas(), cellBb(), checkEdgeLength(), checkFaceAngles(), checkFaceFlatness(), checkFacePyramids(), checkFaceSkewness(), checkPointNearness(), boundaryMesh::getNearest(), treeDataFace::findIntersectOp::operator()(), and points().


|
pure virtual |
Return faces.
Implemented in polyMesh.
References faces().
Referenced by calcPointOrder(), cellBb(), cellPoints(), checkCellsZipUp(), checkCommonOrder(), checkDuplicateFaces(), checkEdgeLength(), checkFaceFaces(), checkFaceFlatness(), checkFacePyramids(), checkFaceVertices(), edgeFaces(), faceEdges(), faceEdges(), faces(), treeDataFace::findNearest(), boundaryMesh::getNearest(), treeDataFace::findIntersectOp::operator()(), pointFaces(), and reset().


|
pure virtual |
Face face-owner addressing.
Implemented in polyMesh.
References faceOwner().
Referenced by cellCells(), checkConcaveCells(), checkFaceAreas(), checkFacePyramids(), checkUpperTriangular(), edgeCells(), faceOwner(), pointCells(), and pointInCell().


|
pure virtual |
Face face-neighbour addressing.
Implemented in polyMesh.
References faceNeighbour().
Referenced by cellCells(), checkFaceAreas(), checkFacePyramids(), checkUpperTriangular(), edgeCells(), faceNeighbour(), and pointCells().


|
pure virtual |
Return old points for mesh motion.
Implemented in polyMesh.
References calcPointOrder(), cellBb(), cellCells(), cellCentres(), cellEdges(), cellPoints(), cells(), cellShapes(), cellVolumes(), checkCellsZipUp(), checkCellVolumes(), checkClosedBoundary(), checkClosedCells(), checkConcaveCells(), checkEdgeLength(), checkFaceAngles(), checkFaceAreas(), checkFaceFaces(), checkFaceFlatness(), checkFaceOrthogonality(), checkFacePyramids(), checkFaceSkewness(), checkFaceVertices(), checkGeometry(), checkMesh(), checkPointNearness(), checkPoints(), checkTopology(), checkUpperTriangular(), clearAddressing(), clearCellGeom(), clearGeom(), clearOut(), edgeCells(), edgeFaces(), edges(), faceAreas(), faceCentres(), faceEdges(), findCell(), findNearestCell(), hasCellCells(), hasCellCentres(), hasCellEdges(), hasCellPoints(), hasCells(), hasCellShapes(), hasCellVolumes(), hasEdgeCells(), hasEdgeFaces(), hasEdges(), hasFaceAreas(), hasFaceCentres(), hasFaceEdges(), hasPointCells(), hasPointEdges(), hasPointFaces(), hasPointPoints(), isInternalFace(), movePoints(), nCells(), nInternalFaces(), nInternalPoints(), Foam::noexcept, nPoints(), oldPoints(), p, pointCells(), pointEdges(), pointFaces(), pointInCell(), pointInCellBB(), pointPoints(), printAllocated(), setAspectThreshold(), setClosedThreshold(), setNonOrthThreshold(), setSkewThreshold(), and updateGeom().
Referenced by movePoints(), and oldPoints().

| const Foam::cellShapeList & cellShapes | ( | ) | const |
Return cell shapes.
Definition at line 279 of file primitiveMesh.C.
Referenced by oldPoints().

| const Foam::edgeList & edges | ( | ) | const |
Return mesh edges. Uses calcEdges.
Definition at line 509 of file primitiveMeshEdges.C.
Referenced by geomCellLooper::cut(), topoCellLooper::cut(), edgeFaces(), faceEdges(), cellLooper::getFirstVertEdge(), nEdges(), oldPoints(), pointPoints(), boundaryMesh::setFeatureEdges(), meshCutAndRemove::setRefinement(), meshCutter::setRefinement(), dynamicRefineFvMesh::unrefine(), and edgeVertex::writeCut().

|
static |
Helper function to calculate cell-face addressing from.
face-cell addressing. If nCells is not provided it will scan for the maximum.
Definition at line 26 of file primitiveMeshCells.C.
References forAll, Foam::max(), nCells(), List< T >::setSize(), UList< T >::size(), and Foam::Zero.

|
static |
Helper function to calculate point ordering. Returns true.
if points already ordered, false and fills pointMap (old to new). Map splits points into those not used by any boundary face and those that are.
Definition at line 75 of file primitiveMesh.C.
References f(), faces(), forAll, nInternalFaces(), nInternalPoints(), nPoints(), and List< T >::resize_nocopy().
Referenced by oldPoints(), and reset().


| const Foam::labelListList & cellCells | ( | ) | const |
Definition at line 94 of file primitiveMeshCellCells.C.
Referenced by cellCells(), cellCells(), and oldPoints().

| const Foam::labelListList & edgeCells | ( | ) | const |
Definition at line 27 of file primitiveMeshEdgeCells.C.
References Foam::abort(), cellEdges(), Foam::endl(), Foam::FatalError, FatalErrorInFunction, Foam::invertManyToMany(), nEdges(), and Foam::Pout.
Referenced by edgeCells(), edgeCells(), and oldPoints().


| const Foam::labelListList & pointCells | ( | ) | const |
Definition at line 185 of file primitiveMeshPointCells.C.
Referenced by checkPoints(), cellDistFuncs::correctBoundaryCells(), cellDistFuncs::correctBoundaryPointCells(), oldPoints(), pointCells(), and pointCells().

| const Foam::cellList & cells | ( | ) | const |
Definition at line 132 of file primitiveMeshCells.C.
Referenced by cellBb(), fvMesh::ClassName(), inversePointDistanceDiffusivity::correct(), wallBoundedParticle::crossEdgeConnectedFace(), geomCellLooper::cut(), topoCellLooper::cut(), dynamicFvMesh::declareRunTimeSelectionTable(), dynamicFvMesh::dynamicFvMesh(), sampledSet::findNearFace(), fvMesh::fvMesh(), fvMesh::fvMesh(), Foam::markGrowFaceCellFace(), oldPoints(), treeDataCell::findIntersectOp::operator()(), polyMesh::polyMesh(), polyMesh::polyMesh(), polyMesh::polyMesh(), polyMesh::readUpdate(), meshCutter::setRefinement(), and meshToMeshMethod::writeConnectivity().

| const Foam::labelListList & edgeFaces | ( | ) | const |
Definition at line 27 of file primitiveMeshEdgeFaces.C.
References Foam::abort(), Foam::endl(), faceEdges(), Foam::FatalError, FatalErrorInFunction, Foam::invertManyToMany(), nEdges(), and Foam::Pout.
Referenced by polyDualMesh::calcFeatures(), edgeCells(), edgeFaces(), edgeFaces(), oldPoints(), boundaryMesh::setFeatureEdges(), meshCutAndRemove::setRefinement(), meshCutter::setRefinement(), and interfaceTrackingFvMesh::update().


| const Foam::labelListList & pointFaces | ( | ) | const |
Definition at line 27 of file primitiveMeshPointFaces.C.
References Foam::endl(), faces(), Foam::invertManyToMany(), nPoints(), and Foam::Pout.
Referenced by checkFaceFaces(), motionSmootherAlgo::checkMesh(), checkPoints(), edgeFaces(), cellLooper::getVertFacesNonEdge(), oldPoints(), pointCells(), and dynamicRefineFvMesh::unrefine().


| const Foam::labelListList & cellEdges | ( | ) | const |
Definition at line 114 of file primitiveMeshCellEdges.C.
Referenced by cellEdges(), cellEdges(), checkCellsZipUp(), geomCellLooper::cut(), edgeCells(), cellLooper::getMisAlignedEdge(), oldPoints(), and meshCutter::setRefinement().

| const Foam::labelListList & faceEdges | ( | ) | const |
Definition at line 532 of file primitiveMeshEdges.C.
References edges(), Foam::endl(), f(), faceEdges(), faces(), forAll, pointEdges(), Foam::Pout, List< T >::setSize(), and UList< T >::size().
Referenced by cellEdges(), edgeFaces(), faceEdges(), faceEdges(), faceEdges(), cellLooper::getFirstVertEdge(), cellLooper::getVertEdgesNonFace(), and oldPoints().


| const Foam::labelListList & pointEdges | ( | ) | const |
Definition at line 520 of file primitiveMeshEdges.C.
Referenced by polyDualMesh::calcFeatures(), faceEdges(), faceEdges(), cellLooper::getVertEdgesNonFace(), oldPoints(), pointPoints(), and dynamicRefineFvMesh::unrefine().

| const Foam::labelListList & pointPoints | ( | ) | const |
Definition at line 87 of file primitiveMeshPointPoints.C.
Referenced by oldPoints(), pointPoints(), and pointPoints().

| const Foam::labelListList & cellPoints | ( | ) | const |
Definition at line 104 of file primitiveMeshCellPoints.C.
Referenced by cellBb(), cellPoints(), cellPoints(), oldPoints(), tetOverlapVolume::overlappingCells(), and meshCutAndRemove::setRefinement().

| const Foam::vectorField & cellCentres | ( | ) | const |
Definition at line 78 of file primitiveMeshCellCentresAndVols.C.
Referenced by polyMesh::checkFaceOrthogonality(), checkFaceOrthogonality(), checkFacePyramids(), polyMesh::checkFaceSkewness(), checkFaceSkewness(), polyMesh::checkFaceWeight(), cellDistFuncs::correctBoundaryCells(), cellDistFuncs::correctBoundaryFaceCells(), cellDistFuncs::correctBoundaryPointCells(), mapNearestMethod::findNearestCell(), findNearestCell(), isoSurfacePoint::isoSurfacePoint(), fvMesh::makeC(), inverseDistance::markDonors(), trackingInverseDistance::markDonors(), polyMesh::movePoints(), polyMesh::oldCellCentres(), oldPoints(), sampledSet::pushIn(), resetGeometry(), polyMesh::updateMesh(), reactingOneDim::updateqr(), and meshToMeshMethod::writeConnectivity().

| const Foam::vectorField & faceCentres | ( | ) | const |
Definition at line 71 of file primitiveMeshFaceCentresAndAreas.C.
Referenced by calcCellCentresAndVols(), calcEdgeVectors(), CentredFitSnGradData< Polynomial >::calcFit(), FitData< FitDataType, ExtendedStencil, Polynomial >::calcFit(), sampledSet::calcSign(), checkConcaveCells(), checkFaceFlatness(), checkFaceFlatness(), polyMesh::checkFaceSkewness(), checkFaceSkewness(), polyMesh::checkFaceWeight(), boundaryMesh::getNearest(), isoSurfacePoint::isoSurfacePoint(), volPointInterpolationAdjoint::makeBoundaryWeights(), fvMesh::makeC(), fvMesh::makeCf(), oldPoints(), treeDataFace::findIntersectOp::operator()(), polyMesh::pointInCell(), pointInCell(), sampledSet::pushIn(), dynamicRefineFvMesh::refine(), resetGeometry(), and wallBoundedParticle::trackToEdge().

| const Foam::scalarField & cellVolumes | ( | ) | const |
Definition at line 90 of file primitiveMeshCellCentresAndVols.C.
Referenced by calcEdgeVectors(), checkCellVolumes(), checkClosedCells(), checkClosedCells(), polyMesh::checkVolRatio(), InflationInjection< CloudType >::InflationInjection(), KinematicParcel< ParcelType >::massCell(), oldPoints(), resetGeometry(), and fvMesh::V().

| const Foam::vectorField & faceAreas | ( | ) | const |
Definition at line 83 of file primitiveMeshFaceCentresAndAreas.C.
Referenced by calcCellCentresAndVols(), calcEdgeVectors(), polyMesh::checkCellDeterminant(), checkClosedBoundary(), checkClosedCells(), checkClosedCells(), checkConcaveCells(), checkFaceAngles(), checkFaceAngles(), checkFaceAreas(), checkFaceAreas(), checkFaceFlatness(), checkFaceFlatness(), polyMesh::checkFaceOrthogonality(), checkFaceOrthogonality(), polyMesh::checkFaceSkewness(), checkFaceSkewness(), polyMesh::checkFaceWeight(), boundaryMesh::getNearest(), fvMesh::makeSf(), oldPoints(), pointInCell(), and resetGeometry().

| void movePoints | ( | const pointField & | p, |
| const pointField & | oldP ) |
Move points.
Definition at line 265 of file primitiveMesh.C.
References oldPoints().
Referenced by polyMesh::movePoints(), and oldPoints().


|
inlinenoexcept |
Return true if given face label is internal to the mesh.
Definition at line 95 of file primitiveMeshI.H.
Referenced by checkFaceAreas(), checkFacePyramids(), edgeCells(), dynamicRefineFvMesh::mapFields(), Foam::markGrowFaceCellFace(), oldPoints(), and dynamicRefineFvMesh::unrefine().

|
virtual |
Check face ordering.
Definition at line 907 of file primitiveMeshCheck.C.
References cells, DebugInFunction, Foam::endl(), faceNeighbour(), faceOwner(), forAll, SortableList< T >::indices(), Foam::Info, HashSet< Key, Hash >::insert(), Foam::labelMax, nInternalFaces(), Foam::reduce(), UPstream::reduceOr(), UList< T >::size(), and SortableList< T >::sort().
Referenced by checkTopology(), and oldPoints().


|
virtual |
Check cell zip-up.
Definition at line 1065 of file primitiveMeshCheck.C.
References cellEdges(), cells, DebugInFunction, Foam::endl(), f(), faces(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::reduce(), and Foam::Zero.
Referenced by checkTopology(), and oldPoints().


|
virtual |
Check uniqueness of face vertices.
Definition at line 1158 of file primitiveMeshCheck.C.
References DebugInFunction, Foam::endl(), f(), faces(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::max(), Foam::min(), nPoints, Foam::reduce(), and UList< T >::size().
Referenced by checkTopology(), and oldPoints().


|
virtual |
Check for unused points.
Definition at line 1226 of file primitiveMeshCheck.C.
References DebugInFunction, UList< T >::empty(), Foam::endl(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), pointCells(), pointFaces(), and Foam::reduce().
Referenced by checkTopology(), and oldPoints().


|
virtual |
Check face-face connectivity.
Definition at line 1495 of file primitiveMeshCheck.C.
References checkCommonOrder(), checkDuplicateFaces(), HashTable< T, Key, Hash >::clear(), DebugInFunction, Foam::endl(), faces(), forAll, Foam::Info, nFaces(), pointFaces(), and Foam::reduce().
Referenced by checkTopology(), and oldPoints().


|
virtual |
Check boundary for closedness.
Definition at line 1591 of file primitiveMeshCheck.C.
References checkClosedBoundary(), and faceAreas().

|
virtual |
Check cells for closedness.
Definition at line 1597 of file primitiveMeshCheck.C.
References cellVolumes(), checkClosedCells(), and faceAreas().

|
virtual |
Check for negative face areas.
Definition at line 1617 of file primitiveMeshCheck.C.
References checkFaceAreas(), and faceAreas().

|
virtual |
Check for negative cell volumes.
Definition at line 1633 of file primitiveMeshCheck.C.
References cellVolumes(), and checkCellVolumes().

|
virtual |
Check for non-orthogonality.
Reimplemented in polyMesh.
Definition at line 1649 of file primitiveMeshCheck.C.
References cellCentres(), checkFaceOrthogonality(), and faceAreas().

|
virtual |
Check face pyramid volume.
Definition at line 1665 of file primitiveMeshCheck.C.
References cellCentres(), checkFacePyramids(), and points.

|
virtual |
Check face skewness.
Reimplemented in polyMesh.
Definition at line 1684 of file primitiveMeshCheck.C.
References cellCentres(), checkFaceSkewness(), faceAreas(), faceCentres(), and points.

|
virtual |
Check face angles.
Definition at line 1702 of file primitiveMeshCheck.C.
References checkFaceAngles(), faceAreas(), and points.

|
virtual |
Check face warpage: decompose face and check ratio between.
magnitude of sum of triangle areas and sum of magnitude of triangle areas.
Definition at line 1720 of file primitiveMeshCheck.C.
References checkFaceFlatness(), faceAreas(), faceCentres(), and points.

|
virtual |
Check for point-point-nearness,.
e.g. colocated points which may be part of baffles.
Definition at line 27 of file primitiveMeshCheckPointNearness.C.
References Foam::endl(), SortableList< T >::indices(), Foam::Info, HashSet< Key, Hash >::insert(), Foam::magSqr(), points(), Foam::reduce(), UList< T >::size(), and Foam::sqrt().
Referenced by oldPoints().


|
virtual |
Check edge length.
Definition at line 26 of file primitiveMeshCheckEdgeLength.C.
References Foam::endl(), f(), faces(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::magSqr(), Foam::max(), Foam::min(), nPoints(), points(), Foam::reduce(), Foam::returnReduce(), HashTable< T, Key, Hash >::size(), Foam::sqr(), Foam::sqrt(), and HashTable< T, Key, Hash >::transfer().
Referenced by oldPoints().


|
virtual |
Check for concave cells by the planes of faces.
Definition at line 1739 of file primitiveMeshCheck.C.
References checkConcaveCells(), faceAreas(), and faceCentres().

|
virtual |
Check mesh topology for correctness.
Returns false for no error.
Definition at line 1755 of file primitiveMeshCheck.C.
References checkCellsZipUp(), checkFaceFaces(), checkFaceVertices(), checkPoints(), checkUpperTriangular(), Foam::endl(), and Foam::Info.
Referenced by oldPoints().


|
virtual |
Check mesh geometry (& implicitly topology) for correctness.
Returns false for no error.
Definition at line 1785 of file primitiveMeshCheck.C.
References checkCellVolumes(), checkClosedBoundary(), checkClosedCells(), checkFaceAreas(), checkFaceOrthogonality(), checkFacePyramids(), checkFaceSkewness(), Foam::endl(), and Foam::Info.
Referenced by oldPoints().


|
virtual |
Check mesh for correctness. Returns false for no error.
Definition at line 1817 of file primitiveMeshCheck.C.
References Foam::checkGeometry(), Foam::checkTopology(), DebugInFunction, Foam::endl(), and Foam::Info.
Referenced by oldPoints(), polyMesh::polyMesh(), and polyMesh::polyMesh().


|
static |
Set the closedness ratio warning threshold.
Definition at line 1843 of file primitiveMeshCheck.C.
References closedThreshold_.
Referenced by oldPoints().

|
static |
Set the aspect ratio warning threshold.
Definition at line 1852 of file primitiveMeshCheck.C.
References aspectThreshold_.
Referenced by oldPoints().

|
static |
Set the non-orthogonality warning threshold in degrees.
Definition at line 1861 of file primitiveMeshCheck.C.
References nonOrthThreshold_.
Referenced by oldPoints().

|
static |
Set the skewness warning threshold as percentage.
of the face area vector
Definition at line 1870 of file primitiveMeshCheck.C.
References skewThreshold_.
Referenced by oldPoints().

| Foam::boundBox cellBb | ( | const label | celli | ) | const |
The bounding box for given cell index.
Definition at line 28 of file primitiveMeshFindCell.C.
References cellPoints(), cells(), faces(), hasCellPoints(), and points().
Referenced by oldPoints(), and pointInCellBB().


| bool pointInCellBB | ( | const point & | p, |
| label | celli, | ||
| scalar | inflationFraction = 0 ) const |
Return true if the point in the cell bounding box.
The bounding box may be isotropically inflated by the fraction inflationFraction
Definition at line 40 of file primitiveMeshFindCell.C.
References cellBb(), boundBox::contains(), boundBox::inflate(), and p.
Referenced by oldPoints().


| bool pointInCell | ( | const point & | p, |
| label | celli ) const |
Return true if the point is in the cell.
Definition at line 58 of file primitiveMeshFindCell.C.
References cells, f(), faceAreas(), faceCentres(), faceOwner(), forAll, and p.
Referenced by findCell(), oldPoints(), and polyMesh::pointInCell().


| Foam::label findNearestCell | ( | const point & | location | ) | const |
Find the cell with the nearest cell centre to location.
Definition at line 85 of file primitiveMeshFindCell.C.
References cellCentres(), Foam::magSqr(), and UList< T >::size().
Referenced by polyMesh::findCell(), findCell(), and oldPoints().


| Foam::label findCell | ( | const point & | location | ) | const |
Find cell enclosing this location (-1 if not in mesh).
Definition at line 112 of file primitiveMeshFindCell.C.
References findNearestCell(), n, nCells(), and pointInCell().
Referenced by oldPoints().


| void printAllocated | ( | ) | const |
Print a list of all the currently allocated mesh data.
Definition at line 26 of file primitiveMeshClear.C.
References Foam::endl(), and Foam::Pout.
Referenced by oldPoints().


|
inlinenoexcept |
Definition at line 104 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 110 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 116 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by cellCells(), and oldPoints().

|
inlinenoexcept |
Definition at line 122 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by edgeCells(), and oldPoints().

|
inlinenoexcept |
Definition at line 128 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints(), and pointCells().

|
inlinenoexcept |
Definition at line 134 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 140 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by edgeFaces(), and oldPoints().

|
inlinenoexcept |
Definition at line 146 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 152 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by cellEdges(), and oldPoints().

|
inlinenoexcept |
Definition at line 158 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by faceEdges(), and oldPoints().

|
inlinenoexcept |
Definition at line 164 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 170 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints(), and pointPoints().

|
inlinenoexcept |
Definition at line 176 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by cellBb(), cellPoints(), and oldPoints().

|
inlinenoexcept |
Definition at line 182 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 188 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 194 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

|
inlinenoexcept |
Definition at line 200 of file primitiveMeshI.H.
References Foam::noexcept.
Referenced by oldPoints().

| const Foam::labelList & cellCells | ( | const label | celli, |
| DynamicList< label > & | storage ) const |
cellCells using cells.
Definition at line 105 of file primitiveMeshCellCells.C.
References cellCells(), cells, DynamicList< T, SizeMin >::clear(), faceNeighbour(), faceOwner(), hasCellCells(), nInternalFaces(), and DynamicList< T, SizeMin >::push_back().

| const Foam::labelList & cellCells | ( | const label | celli | ) | const |
Definition at line 143 of file primitiveMeshCellCells.C.
References cellCells().

| const Foam::labelList & cellPoints | ( | const label | celli, |
| labelHashSet & | set, | ||
| DynamicList< label > & | storage ) const |
cellPoints using cells
Definition at line 115 of file primitiveMeshCellPoints.C.
References DynamicList< T, SizeMin >::capacity(), cellPoints(), cells, DynamicList< T, SizeMin >::clear(), faces(), hasCellPoints(), DynamicList< T, SizeMin >::push_back(), and DynamicList< T, SizeMin >::setCapacity().

| const Foam::labelList & cellPoints | ( | const label | celli | ) | const |
Definition at line 152 of file primitiveMeshCellPoints.C.
References cellPoints().

| const Foam::labelList & pointCells | ( | const label | pointi, |
| DynamicList< label > & | storage ) const |
pointCells using pointFaces
Definition at line 196 of file primitiveMeshPointCells.C.
References UList< T >::begin(), DynamicList< T, SizeMin >::clear(), UList< T >::end(), faceNeighbour(), faceOwner(), hasPointCells(), nInternalFaces(), pFaces, pointCells(), pointFaces(), DynamicList< T, SizeMin >::push_back(), DynamicList< T, SizeMin >::resize(), and UList< T >::size().

| const Foam::labelList & pointCells | ( | const label | pointi | ) | const |
Definition at line 239 of file primitiveMeshPointCells.C.
References pointCells().

| const Foam::labelList & pointPoints | ( | const label | pointi, |
| DynamicList< label > & | storage ) const |
pointPoints using edges, pointEdges
Definition at line 98 of file primitiveMeshPointPoints.C.
References DynamicList< T, SizeMin >::capacity(), DynamicList< T, SizeMin >::clear(), edges(), hasPointPoints(), pointEdges(), pointPoints(), DynamicList< T, SizeMin >::push_back(), DynamicList< T, SizeMin >::setCapacity(), and UList< T >::size().

| const Foam::labelList & pointPoints | ( | const label | pointi | ) | const |
Definition at line 130 of file primitiveMeshPointPoints.C.
References pointPoints().

| const Foam::labelList & faceEdges | ( | const label | facei, |
| DynamicList< label > & | storage ) const |
faceEdges using pointFaces, edges, pointEdges
Definition at line 593 of file primitiveMeshEdges.C.
References DynamicList< T, SizeMin >::capacity(), DynamicList< T, SizeMin >::clear(), f(), faceEdges(), faces(), Foam::findFirstCommonElementFromSortedLists(), forAll, hasFaceEdges(), pointEdges(), DynamicList< T, SizeMin >::push_back(), and DynamicList< T, SizeMin >::setCapacity().

| const Foam::labelList & faceEdges | ( | const label | facei | ) | const |
Definition at line 629 of file primitiveMeshEdges.C.
References faceEdges().

| const Foam::labelList & edgeFaces | ( | const label | edgeI, |
| DynamicList< label > & | storage ) const |
edgeFaces using pointFaces, edges, pointEdges
Definition at line 53 of file primitiveMeshEdgeFaces.C.
References DynamicList< T, SizeMin >::append(), DynamicList< T, SizeMin >::clear(), e, edgeFaces(), edges(), f(), faces(), hasEdgeFaces(), pointFaces(), and UList< T >::size().

| const Foam::labelList & edgeFaces | ( | const label | edgeI | ) | const |
Definition at line 119 of file primitiveMeshEdgeFaces.C.
References edgeFaces().

| const Foam::labelList & edgeCells | ( | const label | edgeI, |
| DynamicList< label > & | storage ) const |
edgeCells using pointFaces, edges, pointEdges
Definition at line 52 of file primitiveMeshEdgeCells.C.
References DynamicList< T, SizeMin >::clear(), UList< T >::contains(), edgeCells(), edgeFaces(), faceNeighbour(), faceOwner(), hasEdgeCells(), isInternalFace(), and DynamicList< T, SizeMin >::push_back().

| const Foam::labelList & edgeCells | ( | const label | edgeI | ) | const |
Definition at line 99 of file primitiveMeshEdgeCells.C.
References edgeCells().

| const Foam::labelList & cellEdges | ( | const label | celli, |
| labelHashSet & | set, | ||
| DynamicList< label > & | storage ) const |
cellEdges using cells, pointFaces, edges, pointEdges
Definition at line 635 of file primitiveMeshEdges.C.
References DynamicList< T, SizeMin >::capacity(), cellEdges(), cells, DynamicList< T, SizeMin >::clear(), faceEdges(), hasCellEdges(), DynamicList< T, SizeMin >::push_back(), and DynamicList< T, SizeMin >::setCapacity().

| const Foam::labelList & cellEdges | ( | const label | celli | ) | const |
Definition at line 671 of file primitiveMeshEdges.C.
References cellEdges().

|
virtual |
Update all geometric data.
Reimplemented in fvMesh.
Definition at line 290 of file primitiveMesh.C.
References calcCellCentresAndVols(), and calcFaceCentresAndAreas().
Referenced by basicFvGeometryScheme::movePoints(), highAspectRatioFvGeometryScheme::movePoints(), polyMesh::movePoints(), solidBodyFvGeometryScheme::movePoints(), and oldPoints().


| void clearGeom | ( | ) |
Clear geometry.
Definition at line 119 of file primitiveMeshClear.C.
References Foam::endl(), and Foam::Pout.
Referenced by polyMesh::clearGeom(), clearOut(), fvGeometryScheme::movePoints(), solidBodyFvGeometryScheme::movePoints(), oldPoints(), resetGeometry(), and polyMesh::updateGeomPoints().


| void clearCellGeom | ( | ) |
Clear cell-based geometry only.
Use with care! currently used by cyclicACMI
Definition at line 135 of file primitiveMeshClear.C.
References Foam::endl(), and Foam::Pout.
Referenced by oldPoints().


| void clearAddressing | ( | ) |
Clear topological data.
Definition at line 149 of file primitiveMeshClear.C.
References Foam::endl(), and Foam::Pout.
Referenced by polyMesh::clearAddressing(), clearOut(), and oldPoints().


| void clearOut | ( | ) |
Clear all geometry and addressing unnecessary for CFD.
Definition at line 178 of file primitiveMeshClear.C.
References clearAddressing(), and clearGeom().
Referenced by fvMesh::clearOut(), oldPoints(), polyMesh::removeBoundary(), reset(), and ~primitiveMesh().


|
staticprotected |
Static data to control mesh checking.
Cell closedness warning threshold
set as the fraction of un-closed area to closed area
Definition at line 307 of file primitiveMesh.H.
Referenced by checkClosedBoundary(), checkClosedCells(), and setClosedThreshold().
|
staticprotected |
Aspect ratio warning threshold.
Definition at line 312 of file primitiveMesh.H.
Referenced by checkClosedCells(), and setAspectThreshold().
|
staticprotected |
Non-orthogonality warning threshold in deg.
Definition at line 317 of file primitiveMesh.H.
Referenced by checkFaceOrthogonality(), and setNonOrthThreshold().
|
staticprotected |
Skewness warning threshold.
Definition at line 322 of file primitiveMesh.H.
Referenced by checkFaceSkewness(), and setSkewThreshold().
|
staticprotected |
Threshold where faces are considered coplanar.
Definition at line 327 of file primitiveMesh.H.
Referenced by checkConcaveCells().
|
static |
Estimated number of cells per edge.
Definition at line 510 of file primitiveMesh.H.
|
static |
Estimated number of cells per point.
Definition at line 515 of file primitiveMesh.H.
|
static |
Estimated number of faces per cell.
Definition at line 520 of file primitiveMesh.H.
|
static |
Estimated number of faces per edge.
Definition at line 525 of file primitiveMesh.H.
|
static |
Estimated number of faces per point.
Definition at line 530 of file primitiveMesh.H.
|
static |
Estimated number of edges per cell.
Definition at line 535 of file primitiveMesh.H.
|
static |
Estimated number of edges per cell.
Definition at line 540 of file primitiveMesh.H.
|
static |
Estimated number of edges per point.
Definition at line 545 of file primitiveMesh.H.
|
static |
Estimated number of points per cell.
Definition at line 550 of file primitiveMesh.H.
|
static |
Estimated number of points per face.
Definition at line 555 of file primitiveMesh.H.