FreeFOAM The Cross-Platform CFD Toolkit
primitiveMesh Class Reference

Cell-face mesh analysis engine. More...

#include <OpenFOAM/primitiveMesh.H>


Detailed Description

+ Inheritance diagram for primitiveMesh:

List of all members.

Public Member Functions

 ClassName ("primitiveMesh")
 primitiveMesh (const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
 Construct from components.
virtual ~primitiveMesh ()
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 reset (const label nPoints, const label nInternalFaces, const label nFaces, const label nCells, const Xfer< cellList > &cells)
 Reset this primitiveMesh given the primitive array sizes and cells.
label nPoints () const
label nEdges () const
label nInternalFaces () const
label nFaces () const
label nCells () const
label nInternalPoints () const
 Points not on boundary.
label nInternal0Edges () const
 Internal edges (i.e. not on boundary face) using.
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 pointFieldpoints () const =0
 Return mesh points.
virtual const faceListfaces () const =0
 Return faces.
virtual const labelListfaceOwner () const =0
 Face face-owner addresing.
virtual const labelListfaceNeighbour () const =0
 Face face-neighbour addressing.
virtual const pointFieldoldPoints () const =0
 Return old points for mesh motion.
const cellShapeListcellShapes () const
 Return cell shapes.
const edgeListedges () const
 Return mesh edges. Uses calcEdges.
const labelListListcellCells () const
const labelListListedgeCells () const
const labelListListpointCells () const
const cellListcells () const
const labelListListedgeFaces () const
const labelListListpointFaces () const
const labelListListcellEdges () const
const labelListListfaceEdges () const
const labelListListpointEdges () const
const labelListListpointPoints () const
const labelListListcellPoints () const
const vectorFieldcellCentres () const
const vectorFieldfaceCentres () const
const scalarFieldcellVolumes () const
const vectorFieldfaceAreas () const
tmp< scalarFieldmovePoints (const pointField &p, const pointField &oldP)
 Move points, returns volumes swept by faces in motion.
bool isInternalFace (const label faceIndex) const
 Return true if given face label is internal to the mesh.
bool checkCellsZipUp (const bool report=false, labelHashSet *setPtr=NULL) const
 Check cell zip-up.
bool checkFaceVertices (const bool report=false, labelHashSet *setPtr=NULL) const
 Check uniqueness of face vertices.
bool checkFaceFaces (const bool report=false, labelHashSet *setPtr=NULL) const
 Check face-face connectivity.
bool checkUpperTriangular (const bool report=false, labelHashSet *setPtr=NULL) const
 Check face ordering.
bool checkClosedBoundary (const bool report=false) const
 Check boundary for closedness.
bool checkClosedCells (const bool report=false, labelHashSet *setPtr=NULL, labelHashSet *highAspectSetPtr=NULL, const Vector< label > &solutionD=Vector< label >::one) const
 Check cells for closedness.
bool checkFaceAreas (const bool report=false, labelHashSet *setPtr=NULL) const
 Check for negative face areas.
bool checkCellVolumes (const bool report=false, labelHashSet *setPtr=NULL) const
 Check for negative cell volumes.
bool checkFaceOrthogonality (const bool report=false, labelHashSet *setPtr=NULL) const
 Check for non-orthogonality.
bool checkFacePyramids (const bool report=false, const scalar minPyrVol=-SMALL, labelHashSet *setPtr=NULL) const
 Check face pyramid volume.
bool checkFaceSkewness (const bool report=false, labelHashSet *setPtr=NULL) const
 Check face skewness.
bool checkFaceAngles (const bool report=false, const scalar maxSin=10, labelHashSet *setPtr=NULL) const
 Check face angles.
bool checkFaceFlatness (const bool report, const scalar warnFlatness, labelHashSet *setPtr) const
 Check face warpage: decompose face and check ratio between.
bool checkEdgeAlignment (const bool report, const Vector< label > &directions, labelHashSet *setPtr=NULL) const
 Check edge alignment for 1D/2D cases.
bool checkPoints (const bool report=false, labelHashSet *setPtr=NULL) const
 Check for unused points.
bool checkPointNearness (const bool report, const scalar reportDistSqr, labelHashSet *setPtr=NULL) const
 Check for point-point-nearness,.
bool checkEdgeLength (const bool report, const scalar minLenSqr, labelHashSet *setPtr=NULL) const
 Check edge length.
bool checkCellDeterminant (const bool report=false, labelHashSet *setPtr=NULL, const Vector< label > &solutionD=Vector< label >::one) const
 Check cell determinant.
bool checkTopology (const bool report=false) const
 Check mesh topology for correctness.
bool checkGeometry (const bool report=false) const
 Check mesh geometry (& implicitly topology) for correctness.
bool checkMesh (const bool report=false) const
 Check mesh for correctness. Returns false for no error.
bool checkMeshMotion (const pointField &newPoints, const bool report=false) const
 Check mesh motion for correctness given motion points.
bool pointInCellBB (const point &p, label celli) const
 Is the point in the cell bounding box.
bool pointInCell (const point &p, label celli) const
 Is the point 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
bool hasEdges () const
bool hasCellCells () const
bool hasEdgeCells () const
bool hasPointCells () const
bool hasCells () const
bool hasEdgeFaces () const
bool hasPointFaces () const
bool hasCellEdges () const
bool hasFaceEdges () const
bool hasPointEdges () const
bool hasPointPoints () const
bool hasCellPoints () const
bool hasCellCentres () const
bool hasFaceCentres () const
bool hasCellVolumes () const
bool hasFaceAreas () const
const labelListcellCells (const label cellI, DynamicList< label > &) const
 cellCells using cells.
const labelListcellCells (const label cellI) const
const labelListcellPoints (const label cellI, DynamicList< label > &) const
 cellPoints using cells
const labelListcellPoints (const label cellI) const
const labelListpointCells (const label pointI, DynamicList< label > &) const
 pointCells using pointFaces
const labelListpointCells (const label pointI) const
const labelListpointPoints (const label pointI, DynamicList< label > &) const
 pointPoints using edges, pointEdges
const labelListpointPoints (const label pointI) const
const labelListfaceEdges (const label faceI, DynamicList< label > &) const
 faceEdges using pointFaces, edges, pointEdges
const labelListfaceEdges (const label faceI) const
const labelListedgeFaces (const label edgeI, DynamicList< label > &) const
 edgeFaces using pointFaces, edges, pointEdges
const labelListedgeFaces (const label edgeI) const
const labelListedgeCells (const label edgeI, DynamicList< label > &) const
 edgeCells using pointFaces, edges, pointEdges
const labelListedgeCells (const label edgeI) const
const labelListcellEdges (const label cellI, DynamicList< label > &) const
 cellEdges using cells, pointFaces, edges, pointEdges
const labelListcellEdges (const label cellI) const
void clearGeom ()
 Clear geometry.
void clearAddressing ()
 Clear topological data.
void clearOut ()
 Clear all geometry and addressing unnecessary for CFD.

Static Public Member Functions

static void calcCells (cellList &, const unallocLabelList &own, const unallocLabelList &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

 primitiveMesh ()
 Construct null.

Constructor & Destructor Documentation

primitiveMesh ( )
protected

Construct null.

Definition at line 38 of file primitiveMesh.C.

primitiveMesh ( const label  nPoints,
const label  nInternalFaces,
const label  nFaces,
const label  nCells 
)

Construct from components.

Definition at line 78 of file primitiveMesh.C.

~primitiveMesh ( )
virtual

Definition at line 119 of file primitiveMesh.C.


Member Function Documentation

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 207 of file primitiveMesh.C.

References Foam::endl(), nPoints, and Foam::Pout.

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 reset ( const label  nPoints,
const label  nInternalFaces,
const label  nFaces,
const label  nCells,
const Xfer< cellList > &  cells 
)

Reset this primitiveMesh given the primitive array sizes and cells.

label nInternalFaces ( ) const
inline

Definition at line 82 of file primitiveMeshI.H.

Referenced by polyMeshAdder::add(), polyTopoChange::addMesh(), scotchDecomp::calcCSR(), polyDualMesh::calcFeatures(), primitiveMesh::checkClosedBoundary(), meshRefinement::checkCoupledFaceZones(), polyBoundaryMesh::checkDefinition(), polyMeshGeometry::checkFaceDotProduct(), polyMeshGeometry::checkFaceSkewness(), polyMeshGeometry::checkFaceTwist(), polyMeshGeometry::checkFaceWeights(), faceZone::checkParallelSync(), polyMeshGeometry::checkVolRatio(), inversePointDistanceDiffusivity::correct(), inverseFaceDistanceDiffusivity::correct(), fvMeshDistribute::distribute(), dynamicRefineFvMesh::dynamicRefineFvMesh(), dynamicRefineFvMesh::extendMarkedCells(), faceCoupleInfo::faceCoupleInfo(), localPointRegion::findDuplicateFaces(), boundaryMesh::getNearest(), cellToCellStencil::insertFaceCells(), cellToFaceStencil::insertFaceCells(), fvMeshSubset::interpolate(), primitiveMesh::isInternalFace(), MeshedSurface< Face >::MeshedSurface(), octreeDataFace::octreeDataFace(), boundaryMesh::patchify(), fvMeshDistribute::printCoupleInfo(), fvMeshDistribute::printMeshInfo(), quadraticFitSnGradData::quadraticFitSnGradData(), boundaryMesh::read(), surfaceMesh::size(), fvSurfaceMapper::size(), syncTools::syncBoundaryFaceList(), syncTools::syncEdgeMap(), syncTools::syncFaceList(), triSurfaceTools::triangulate(), triSurfaceTools::triangulateFaceCentre(), cellToCellStencil::validBoundaryFaces(), cellToFaceStencil::validBoundaryFaces(), extendedCellToFaceStencil::weightedSum(), extendedUpwindCellToFaceStencil::weightedSum(), and polyBoundaryMesh::whichPatch().

label nFaces ( ) const
inline

Definition at line 88 of file primitiveMeshI.H.

Referenced by polyMeshAdder::add(), meshRefinement::addPatch(), polyDualMesh::calcFeatures(), polyTopoChange::changeMesh(), meshRefinement::checkCoupledFaceZones(), polyMeshGeometry::checkFaceDotProduct(), polyMeshGeometry::checkFaceSkewness(), polyMeshGeometry::checkFaceTwist(), polyMeshGeometry::checkFaceWeights(), motionSmoother::checkMesh(), faceZone::checkParallelSync(), polyMeshGeometry::checkVolRatio(), autoSnapDriver::createZoneBaffles(), fvMeshDistribute::distribute(), autoLayerDriver::doLayers(), autoSnapDriver::doSnap(), dynamicRefineFvMesh::dynamicRefineFvMesh(), evaluateError::evaluateError(), Foam::MULES::explicitSolve(), dynamicRefineFvMesh::extendMarkedCells(), faceMapper::faceMapper(), faceSet::faceSet(), faceZoneSet::faceZoneSet(), syncTools::getMasterFaces(), boundaryMesh::getNearest(), Foam::MULES::implicitSolve(), polyTopoChange::makeMesh(), faceSet::maxSize(), faceZoneSet::maxSize(), MeshedSurface< Face >::MeshedSurface(), octreeDataFace::octreeDataFace(), boundaryMesh::patchify(), autoSnapDriver::preSmoothPatch(), fvMeshDistribute::printMeshInfo(), boundaryMesh::read(), autoSnapDriver::repatchToSurface(), autoSnapDriver::scaleMesh(), syncTools::syncBoundaryFaceList(), syncTools::syncEdgeMap(), syncTools::syncFaceList(), triSurfaceTools::triangulate(), triSurfaceTools::triangulateFaceCentre(), directMappedVelocityFluxFixedValueFvPatchField::updateCoeffs(), directMappedFixedValueFvPatchField< Type >::updateCoeffs(), and polyBoundaryMesh::whichPatch().

label nInternalPoints ( ) const
inline

Points not on boundary.

Definition at line 35 of file primitiveMeshI.H.

label nInternal0Edges ( ) const
inline

Internal edges (i.e. not on boundary face) using.

no boundary point

Definition at line 47 of file primitiveMeshI.H.

References primitiveMesh::nEdges().

label nInternal1Edges ( ) const
inline

Internal edges using 0 or 1 boundary point.

Definition at line 55 of file primitiveMeshI.H.

References primitiveMesh::nEdges().

label nInternalEdges ( ) const
inline

Internal edges using 0,1 or 2 boundary points.

Definition at line 63 of file primitiveMeshI.H.

References primitiveMesh::nEdges().

virtual const pointField& oldPoints ( ) const
pure virtual

Return old points for mesh motion.

Implemented in polyMesh.

const Foam::cellShapeList & cellShapes ( ) const

Return cell shapes.

Definition at line 339 of file primitiveMesh.C.

static void calcCells ( cellList ,
const unallocLabelList own,
const unallocLabelList nei,
const label  nCells = -1 
)
static

Helper function to calculate cell-face addressing from.

face-cell addressing. If nCells is not provided it will scan for the maximum.

bool calcPointOrder ( label &  nInternalPoints,
labelList pointMap,
const faceList faces,
const label  nInternalFaces,
const label  nPoints 
)
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 128 of file primitiveMesh.C.

References f(), forAll, List< T >::setSize(), and List< T >::size().

const Foam::labelListList & cellCells ( ) const

Definition at line 100 of file primitiveMeshCellCells.C.

const Foam::labelListList & pointPoints ( ) const

Definition at line 91 of file primitiveMeshPointPoints.C.

Foam::tmp< Foam::scalarField > movePoints ( const pointField p,
const pointField oldP 
)

Move points, returns volumes swept by faces in motion.

Definition at line 305 of file primitiveMesh.C.

References Foam::abort(), f(), Foam::FatalError, FatalErrorIn, forAll, nPoints, and List< T >::size().

Referenced by polyMesh::movePoints().

bool checkCellsZipUp ( const bool  report = false,
labelHashSet setPtr = NULL 
) const

Check cell zip-up.

Definition at line 1402 of file primitiveMeshCheck.C.

References cells, Foam::endl(), f(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::reduce(), and List< T >::size().

bool checkFaceVertices ( const bool  report = false,
labelHashSet setPtr = NULL 
) const

Check uniqueness of face vertices.

Definition at line 1503 of file primitiveMeshCheck.C.

References Foam::endl(), f(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::max(), Foam::min(), nPoints, Foam::reduce(), and List< T >::size().

bool checkFaceFaces ( const bool  report = false,
labelHashSet setPtr = NULL 
) const
bool checkUpperTriangular ( const bool  report = false,
labelHashSet setPtr = NULL 
) const

Check face ordering.

Definition at line 1237 of file primitiveMeshCheck.C.

References cells, Foam::endl(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::reduce(), and List< T >::size().

bool checkClosedBoundary ( const bool  report = false) const
bool checkClosedCells ( const bool  report = false,
labelHashSet setPtr = NULL,
labelHashSet highAspectSetPtr = NULL,
const Vector< label > &  solutionD = Vector<label>::one 
) const
bool checkFaceAreas ( const bool  report = false,
labelHashSet setPtr = NULL 
) const

Check for negative face areas.

Definition at line 291 of file primitiveMeshCheck.C.

References Foam::endl(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::mag(), Foam::max(), Foam::min(), and Foam::reduce().

bool checkCellVolumes ( const bool  report = false,
labelHashSet setPtr = NULL 
) const

Check for negative cell volumes.

Definition at line 350 of file primitiveMeshCheck.C.

References Foam::endl(), forAll, Foam::gSum(), Foam::Info, HashSet< Key, Hash >::insert(), Foam::max(), Foam::min(), and Foam::reduce().

bool checkFaceOrthogonality ( const bool  report = false,
labelHashSet setPtr = NULL 
) const
bool checkFacePyramids ( const bool  report = false,
const scalar  minPyrVol = -SMALL,
labelHashSet setPtr = NULL 
) const

Check face pyramid volume.

Definition at line 537 of file primitiveMeshCheck.C.

References Foam::endl(), f(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), p, points, and Foam::reduce().

bool checkFaceSkewness ( const bool  report = false,
labelHashSet setPtr = NULL 
) const
bool checkFaceAngles ( const bool  report = false,
const scalar  maxSin = 10,
labelHashSet setPtr = NULL 
) const
bool checkFaceFlatness ( const bool  report,
const scalar  warnFlatness,
labelHashSet setPtr 
) const

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 972 of file primitiveMeshCheck.C.

References Foam::endl(), Foam::exit(), f(), Foam::FatalError, FatalErrorIn, forAll, Foam::Info, HashSet< Key, Hash >::insert(), Foam::mag(), Foam::min(), p, points, Foam::reduce(), and List< T >::size().

bool checkPoints ( const bool  report = false,
labelHashSet setPtr = NULL 
) const

Check for unused points.

Definition at line 766 of file primitiveMeshCheck.C.

References UList< T >::empty(), Foam::endl(), forAll, Foam::Info, HashSet< Key, Hash >::insert(), pc, and Foam::reduce().

bool checkPointNearness ( const bool  report,
const scalar  reportDistSqr,
labelHashSet setPtr = NULL 
) const

Check for point-point-nearness,.

e.g. colocated points which may be part of baffles.

Definition at line 32 of file primitiveMeshCheckPointNearness.C.

References Foam::endl(), SortableList< T >::indices(), Foam::Info, HashSet< Key, Hash >::insert(), magSqr(), points, reduce(), List< T >::size(), and Foam::sqrt().

bool checkEdgeLength ( const bool  report,
const scalar  minLenSqr,
labelHashSet setPtr = NULL 
) const
bool checkCellDeterminant ( const bool  report = false,
labelHashSet setPtr = NULL,
const Vector< label > &  solutionD = Vector<label>::one 
) const
bool checkTopology ( const bool  report = false) const

Check mesh topology for correctness.

Returns false for no error.

Definition at line 2059 of file primitiveMeshCheck.C.

References Foam::endl(), and Foam::Info.

bool checkGeometry ( const bool  report = false) const

Check mesh geometry (& implicitly topology) for correctness.

Returns false for no error.

Definition at line 2092 of file primitiveMeshCheck.C.

References Foam::endl(), and Foam::Info.

bool checkMesh ( const bool  report = false) const

Check mesh for correctness. Returns false for no error.

Definition at line 2125 of file primitiveMeshCheck.C.

References Foam::checkGeometry(), Foam::checkTopology(), Foam::endl(), and Foam::Info.

Referenced by attachPolyTopoChanger::attach().

bool checkMeshMotion ( const pointField newPoints,
const bool  report = false 
) const

Check mesh motion for correctness given motion points.

Definition at line 41 of file primitiveMeshCheckMotion.C.

References Foam::acos(), cells, cellVols, d, Foam::endl(), f(), forAll, Foam::mag(), Foam::min(), Foam::nl, Foam::mathematicalConstant::pi(), Foam::Pout, Foam::sum(), and WarningIn.

Referenced by polyMesh::movePoints().

Foam::scalar setClosedThreshold ( const scalar  val)
static

Set the closedness ratio warning threshold.

Definition at line 2157 of file primitiveMeshCheck.C.

Foam::scalar setAspectThreshold ( const scalar  val)
static

Set the aspect ratio warning threshold.

Definition at line 2166 of file primitiveMeshCheck.C.

Foam::scalar setNonOrthThreshold ( const scalar  val)
static

Set the non-orthogonality warning threshold in degrees.

Definition at line 2175 of file primitiveMeshCheck.C.

Foam::scalar setSkewThreshold ( const scalar  val)
static

Set the skewness warning threshold as percentage.

of the face area vector

Definition at line 2184 of file primitiveMeshCheck.C.

bool pointInCellBB ( const point p,
label  celli 
) const
bool pointInCell ( const point p,
label  celli 
) const

Is the point in the cell.

Definition at line 65 of file primitiveMeshFindCell.C.

References cells, f(), and forAll.

Foam::label findNearestCell ( const point location) const

Find the cell with the nearest cell centre to location.

Definition at line 91 of file primitiveMeshFindCell.C.

References Foam::magSqr(), and List< T >::size().

Foam::label findCell ( const point location) const
void printAllocated ( ) const

Print a list of all the currently allocated mesh data.

Definition at line 32 of file primitiveMeshClear.C.

References Foam::endl(), and Foam::Pout.

bool hasCellShapes ( ) const
inline

Definition at line 106 of file primitiveMeshI.H.

bool hasEdges ( ) const
inline

Definition at line 112 of file primitiveMeshI.H.

bool hasCellCells ( ) const
inline

Definition at line 118 of file primitiveMeshI.H.

bool hasEdgeCells ( ) const
inline

Definition at line 124 of file primitiveMeshI.H.

bool hasPointCells ( ) const
inline

Definition at line 130 of file primitiveMeshI.H.

bool hasCells ( ) const
inline

Definition at line 136 of file primitiveMeshI.H.

bool hasEdgeFaces ( ) const
inline

Definition at line 142 of file primitiveMeshI.H.

bool hasPointFaces ( ) const
inline

Definition at line 148 of file primitiveMeshI.H.

bool hasCellEdges ( ) const
inline

Definition at line 154 of file primitiveMeshI.H.

bool hasFaceEdges ( ) const
inline

Definition at line 160 of file primitiveMeshI.H.

bool hasPointEdges ( ) const
inline

Definition at line 166 of file primitiveMeshI.H.

bool hasPointPoints ( ) const
inline

Definition at line 172 of file primitiveMeshI.H.

bool hasCellPoints ( ) const
inline

Definition at line 178 of file primitiveMeshI.H.

bool hasCellCentres ( ) const
inline

Definition at line 184 of file primitiveMeshI.H.

bool hasFaceCentres ( ) const
inline

Definition at line 190 of file primitiveMeshI.H.

bool hasCellVolumes ( ) const
inline

Definition at line 196 of file primitiveMeshI.H.

bool hasFaceAreas ( ) const
inline

Definition at line 202 of file primitiveMeshI.H.

const Foam::labelList & cellCells ( const label  cellI,
DynamicList< label > &  storage 
) const
const Foam::labelList & cellCells ( const label  cellI) const

Definition at line 151 of file primitiveMeshCellCells.C.

const Foam::labelList & cellPoints ( const label  cellI) const

Definition at line 101 of file primitiveMeshCellPoints.C.

const Foam::labelList & pointCells ( const label  pointI) const

Definition at line 174 of file primitiveMeshPointCells.C.

const Foam::labelList & pointPoints ( const label  pointI) const

Definition at line 135 of file primitiveMeshPointPoints.C.

const Foam::labelList & faceEdges ( const label  faceI) const

Definition at line 621 of file primitiveMeshEdges.C.

const Foam::labelList & edgeFaces ( const label  edgeI,
DynamicList< label > &  storage 
) const
const Foam::labelList & edgeFaces ( const label  edgeI) const

Definition at line 103 of file primitiveMeshEdgeFaces.C.

const Foam::labelList & edgeCells ( const label  edgeI,
DynamicList< label > &  storage 
) const

edgeCells using pointFaces, edges, pointEdges

Definition at line 58 of file primitiveMeshEdgeCells.C.

References DynamicList< T, SizeInc, SizeMult, SizeDiv >::append(), DynamicList< T, SizeInc, SizeMult, SizeDiv >::clear(), and forAll.

const Foam::labelList & edgeCells ( const label  edgeI) const

Definition at line 127 of file primitiveMeshEdgeCells.C.

const Foam::labelList & cellEdges ( const label  cellI) const

Definition at line 670 of file primitiveMeshEdges.C.

void clearGeom ( )

Clear geometry.

Reimplemented in polyMesh.

Definition at line 126 of file primitiveMeshClear.C.

References Foam::deleteDemandDrivenData(), Foam::endl(), and Foam::Pout.

Referenced by polyMesh::clearGeom().

void clearAddressing ( )

Clear topological data.

Reimplemented in polyMesh.

Definition at line 142 of file primitiveMeshClear.C.

References Foam::deleteDemandDrivenData(), Foam::endl(), and Foam::Pout.

Referenced by polyMesh::clearAddressing().

void clearOut ( )

Clear all geometry and addressing unnecessary for CFD.

Reimplemented in polyMesh, and fvMesh.

Definition at line 171 of file primitiveMeshClear.C.


Member Data Documentation

const unsigned cellsPerEdge_ = 4
static

Estimated number of cells per edge.

Definition at line 308 of file primitiveMesh.H.

const unsigned cellsPerPoint_ = 8
static

Estimated number of cells per point.

Definition at line 311 of file primitiveMesh.H.

const unsigned facesPerCell_ = 6
static

Estimated number of faces per cell.

Definition at line 314 of file primitiveMesh.H.

const unsigned facesPerEdge_ = 4
static

Estimated number of faces per edge.

Definition at line 317 of file primitiveMesh.H.

const unsigned facesPerPoint_ = 12
static

Estimated number of faces per point.

Definition at line 320 of file primitiveMesh.H.

const unsigned edgesPerCell_ = 12
static

Estimated number of edges per cell.

Definition at line 323 of file primitiveMesh.H.

const unsigned edgesPerFace_ = 4
static

Estimated number of edges per cell.

Definition at line 326 of file primitiveMesh.H.

const unsigned edgesPerPoint_ = 6
static

Estimated number of edges per point.

Definition at line 329 of file primitiveMesh.H.

const unsigned pointsPerCell_ = 8
static

Estimated number of points per cell.

Definition at line 332 of file primitiveMesh.H.

const unsigned pointsPerFace_ = 4
static

Estimated number of points per face.

Definition at line 335 of file primitiveMesh.H.


The documentation for this class was generated from the following files: