Storage and named access for the indices of a tet which is part of the decomposition of a cell. More...
Public Member Functions | |
tetIndices () | |
Construct null. More... | |
tetIndices (label cellI, label faceI, label faceBasePtI, label facePtAI, label facePtBI, label tetPtI) | |
Construct from components. More... | |
tetIndices (label cellI, label faceI, label tetPtI, const polyMesh &mesh) | |
Construct from cell, face, pt description of tet. More... | |
~tetIndices () | |
Destructor. More... | |
label | cell () const |
Return the cell. More... | |
label | face () const |
Return the face. More... | |
label | faceBasePt () const |
Return the face base point. More... | |
label | facePtA () const |
Return face point A. More... | |
label | facePtB () const |
Return face point B. More... | |
label | tetPt () const |
Return the characterising tetPtI. More... | |
tetPointRef | tet (const polyMesh &mesh) const |
Return the geometry corresponding to this tet from the. More... | |
tetPointRef | oldTet (const polyMesh &mesh) const |
Return the geometry corresponding to this tet from the. More... | |
triPointRef | faceTri (const polyMesh &mesh) const |
Return the geometry corresponding to the tri on the. More... | |
triFace | faceTriIs (const polyMesh &mesh) const |
Return the point indices corresponding to the tri on the mesh. More... | |
triPointRef | oldFaceTri (const polyMesh &mesh) const |
Return the geometry corresponding to the tri on the. More... | |
label & | cell () |
Return non-const access to the cell. More... | |
label & | face () |
Return non-const access to the face. More... | |
label & | faceBasePt () |
Return non-const access to the face base point. More... | |
label & | facePtA () |
Return non-const access to face point A. More... | |
label & | facePtB () |
Return non-const access to face point B. More... | |
label & | tetPt () |
Return non-const access to the characterising tetPtI. More... | |
bool | operator== (const tetIndices &) const |
bool | operator!= (const tetIndices &) const |
Private Attributes | |
label | cellI_ |
Cell that this is a decomposed tet of. More... | |
label | faceI_ |
Face that holds this decomposed tet. More... | |
label | faceBasePtI_ |
Base point on the face. More... | |
label | facePtAI_ |
Point on the face such that the right-hand circulation. More... | |
label | facePtBI_ |
Point on the face such that the right-hand circulation. More... | |
label | tetPtI_ |
Point on the face, *relative to the base point*, which. More... | |
Friends | |
Istream & | operator>> (Istream &, tetIndices &) |
Ostream & | operator<< (Ostream &, const tetIndices &) |
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
facePtB is next one after/before facePtA e.g.:
+—+ |2 /| | / | |/ 1| <- tetPt (so 1 for first triangle, 2 for second) +—+ ^ faceBasePt
Definition at line 73 of file tetIndices.H.
tetIndices | ( | ) |
Construct null.
Definition at line 30 of file tetIndices.C.
tetIndices | ( | label | cellI, |
label | faceI, | ||
label | faceBasePtI, | ||
label | facePtAI, | ||
label | facePtBI, | ||
label | tetPtI | ||
) |
Construct from components.
Definition at line 42 of file tetIndices.C.
tetIndices | ( | label | cellI, |
label | faceI, | ||
label | tetPtI, | ||
const polyMesh & | mesh | ||
) |
Construct from cell, face, pt description of tet.
Definition at line 61 of file tetIndices.C.
References f(), polyMesh::faceOwner(), polyMesh::faces(), mesh, pFaces, List::size(), and polyMesh::tetBasePtIs().
~tetIndices | ( | ) |
Destructor.
Definition at line 102 of file tetIndices.C.
|
inline |
Return the cell.
Definition at line 30 of file tetIndicesI.H.
References tetIndices::cellI_.
Referenced by Basic< Type >::add(), Moment< Type >::add(), Dual< Type >::add(), polyMeshTetDecomposition::faceTetIndices(), polyMeshTetDecomposition::findTet(), interpolationCellPointWallModified< Type >::interpolate(), interpolationCellPoint< Foam::Vector >::interpolate(), Basic< Type >::interpolate(), interpolation< Foam::Vector >::interpolate(), Moment< Type >::interpolate(), Dual< Type >::interpolate(), Basic< Type >::interpolateGrad(), Moment< Type >::interpolateGrad(), Dual< Type >::interpolateGrad(), tetIndices::operator==(), wallBoundedStreamLine::track(), particle< Type >::trackToFace(), and polyMeshTetDecomposition::triangleTetIndices().
|
inline |
Return the face.
Definition at line 36 of file tetIndicesI.H.
Referenced by Dual< Type >::Dual(), polyMeshTetDecomposition::faceTetIndices(), polyMesh::findTetFacePt(), cellPointWeight::findTetrahedron(), DSMCCloud< DSMCParcel< ParcelType > >::initialise(), interpolationCellPointWallModified< Type >::interpolate(), interpolationCellPoint< Foam::Vector >::interpolate(), Moment< Type >::Moment(), tetIndices::operator==(), Dual< Type >::tetGeometry(), wallBoundedStreamLine::track(), wallBoundedParticle::trackToEdge(), particle< Type >::trackToFace(), polyMeshTetDecomposition::triangleTetIndices(), and AveragingMethod< Foam::Vector >::write().
|
inline |
Return the face base point.
Definition at line 42 of file tetIndicesI.H.
Referenced by Dual< Type >::Dual(), polyMeshTetDecomposition::faceTetIndices(), cellPointWeight::findTetrahedron(), cellPointWeight::findTriangle(), particle< Type >::hitWallFaces(), interpolationCellPoint< Foam::Vector >::interpolate(), Moment< Type >::Moment(), tetIndices::operator==(), Dual< Type >::tetGeometry(), wallBoundedParticle::trackToEdge(), polyMeshTetDecomposition::triangleTetIndices(), and AveragingMethod< Foam::Vector >::write().
|
inline |
Return face point A.
Definition at line 48 of file tetIndicesI.H.
Referenced by Dual< Type >::Dual(), polyMeshTetDecomposition::faceTetIndices(), cellPointWeight::findTetrahedron(), cellPointWeight::findTriangle(), interpolationCellPoint< Foam::Vector >::interpolate(), Moment< Type >::Moment(), tetIndices::operator==(), Dual< Type >::tetGeometry(), wallBoundedParticle::trackToEdge(), polyMeshTetDecomposition::triangleTetIndices(), and AveragingMethod< Foam::Vector >::write().
|
inline |
Return face point B.
Definition at line 54 of file tetIndicesI.H.
Referenced by Dual< Type >::Dual(), polyMeshTetDecomposition::faceTetIndices(), cellPointWeight::findTetrahedron(), cellPointWeight::findTriangle(), interpolationCellPoint< Foam::Vector >::interpolate(), Moment< Type >::Moment(), tetIndices::operator==(), Dual< Type >::tetGeometry(), wallBoundedParticle::trackToEdge(), polyMeshTetDecomposition::triangleTetIndices(), and AveragingMethod< Foam::Vector >::write().
|
inline |
Return the characterising tetPtI.
Definition at line 60 of file tetIndicesI.H.
Referenced by polyMeshTetDecomposition::faceTetIndices(), polyMesh::findTetFacePt(), particle< Type >::hitWallFaces(), FreeStream< CloudType >::inflow(), DSMCCloud< DSMCParcel< ParcelType > >::initialise(), tetIndices::operator==(), wallBoundedStreamLine::track(), and polyMeshTetDecomposition::triangleTetIndices().
|
inline |
Return the geometry corresponding to this tet from the.
supplied mesh
Definition at line 66 of file tetIndicesI.H.
References primitiveMesh::cellCentres(), f(), polyMesh::faces(), mesh, pFaces, and polyMesh::points().
Referenced by Dual< Type >::Dual(), polyMeshTetDecomposition::findTet(), cellPointWeight::findTetrahedron(), particle< Type >::hitWallFaces(), DSMCCloud< DSMCParcel< ParcelType > >::initialise(), interpolationCellPoint< Foam::Vector >::interpolate(), Moment< Type >::Moment(), particle< Type >::movingTetLambda(), Dual< Type >::tetGeometry(), particle< Type >::trackToFace(), Implicit< CloudType >::velocityCorrection(), and AveragingMethod< Foam::Vector >::write().
|
inline |
Return the geometry corresponding to this tet from the.
supplied mesh using the old positions
Definition at line 84 of file tetIndicesI.H.
References primitiveMesh::cells(), f(), polyMesh::faces(), mesh, polyMesh::oldPoints(), and pFaces.
Referenced by particle< Type >::movingTetLambda().
|
inline |
Return the geometry corresponding to the tri on the.
mesh face for this tet from the supplied mesh. Normal of the tri points out of the cell.
Definition at line 109 of file tetIndicesI.H.
References f(), polyMesh::faces(), mesh, pFaces, and polyMesh::points().
Referenced by cellPointWeight::findTriangle(), particle< Type >::hitWallFaces(), solidParticle::hitWallPatch(), FreeStream< CloudType >::inflow(), KinematicCloud< CloudType >::patchData(), polyMesh::pointInCell(), wallBoundedStreamLine::track(), wallBoundedParticle::trackToEdge(), and particle< Type >::trackToFace().
|
inline |
Return the point indices corresponding to the tri on the mesh.
face for this tet from the supplied mesh. Normal of the tri points out of the cell.
Definition at line 125 of file tetIndicesI.H.
References f(), polyMesh::faces(), mesh, pFaces, and triFace().
|
inline |
Return the geometry corresponding to the tri on the.
mesh face for this tet from the supplied mesh using the old position
Definition at line 140 of file tetIndicesI.H.
References f(), polyMesh::faces(), mesh, polyMesh::oldPoints(), and pFaces.
Referenced by particle< Type >::hitWallFaces(), and KinematicCloud< CloudType >::patchData().
|
inline |
Return non-const access to the face base point.
|
inline |
Return non-const access to face point A.
|
inline |
Return non-const access to face point B.
|
inline |
Return non-const access to the characterising tetPtI.
|
inline |
Definition at line 194 of file tetIndicesI.H.
References tetIndices::cell(), tetIndices::face(), tetIndices::faceBasePt(), tetIndices::facePtA(), tetIndices::facePtB(), and tetIndices::tetPt().
|
inline |
Definition at line 208 of file tetIndicesI.H.
|
friend |
|
friend |
|
private |
Cell that this is a decomposed tet of.
Definition at line 78 of file tetIndices.H.
Referenced by tetIndices::cell().
|
private |
Face that holds this decomposed tet.
Definition at line 81 of file tetIndices.H.
|
private |
Base point on the face.
Definition at line 84 of file tetIndices.H.
|
private |
Point on the face such that the right-hand circulation.
{faceBasePtI_, facePtAI_, facePtBI_} forms a triangle that points out of the tet
Definition at line 89 of file tetIndices.H.
|
private |
Point on the face such that the right-hand circulation.
{faceBasePtI_, facePtAI_, facePtBI_} forms a triangle that points out of the tet
Definition at line 94 of file tetIndices.H.
|
private |
Point on the face, *relative to the base point*, which.
characterises this tet on the face.
Definition at line 98 of file tetIndices.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.