tetIndices Class Reference

Storage and named access for the indices of a tet which is part of the decomposition of a cell. More...

Public Member Functions

constexpr tetIndices () noexcept
 Default construct, with invalid labels (-1) More...
 
constexpr tetIndices (label celli, label facei, label tetPointi) noexcept
 Construct from components. More...
 
 ~tetIndices ()=default
 Destructor. More...
 
label cell () const noexcept
 Return the cell index. More...
 
label & cell () noexcept
 Non-const access to the cell index. More...
 
label face () const noexcept
 Return the face index. More...
 
label & face () noexcept
 Non-const access to the face index. More...
 
label tetPt () const noexcept
 Return the characterising tet point index. More...
 
label & tetPt () noexcept
 Non-const access to the characterising tet point index. More...
 
triFace faceTriIs (const polyMesh &mesh, const bool warn=true) const
 Return the indices corresponding to the tri on the face for. More...
 
triFace triIs (const polyMesh &mesh, const bool warn=true) const
 
tetPointRef tet (const polyMesh &mesh) const
 Return the geometry corresponding to this tet. More...
 
triPointRef faceTri (const polyMesh &mesh) const
 
triPointRef oldFaceTri (const polyMesh &mesh) const
 

Static Public Member Functions

static label compare (const tetIndices &a, const tetIndices &b) noexcept
 

Friends

Istreamoperator>> (Istream &, tetIndices &)
 
Ostreamoperator<< (Ostream &, const tetIndices &)
 

Detailed Description

Storage and named access for the indices of a tet which is part of the decomposition of a cell.

Tets are designated by

  • cell (of course)
  • face on cell
  • three points on face (faceBasePt, facePtA, facePtB) When constructing from a mesh and index in the face (tetPtI):
    • faceBasePt is the mesh.tetBasePtIs() base point
    • facePtA is tetPtI away from faceBasePt
    • facePtB is next one after/before facePtA e.g.:

      +—+ |2 /| | / | |/ 1| <- tetPt (so 1 for first triangle, 2 for second) +—+ ^ faceBasePt

Source files

Definition at line 83 of file tetIndices.H.

Constructor & Destructor Documentation

◆ tetIndices() [1/2]

constexpr tetIndices ( )
inlineconstexprnoexcept

Default construct, with invalid labels (-1)

Definition at line 31 of file tetIndicesI.H.

◆ tetIndices() [2/2]

constexpr tetIndices ( label  celli,
label  facei,
label  tetPointi 
)
inlineconstexprnoexcept

Construct from components.

Definition at line 39 of file tetIndicesI.H.

◆ ~tetIndices()

~tetIndices ( )
default

Destructor.

Member Function Documentation

◆ cell() [1/2]

◆ cell() [2/2]

label & cell ( )
inlinenoexcept

Non-const access to the cell index.

Definition at line 136 of file tetIndices.H.

◆ face() [1/2]

label face ( ) const
inlinenoexcept

Return the face index.

Definition at line 139 of file tetIndices.H.

Referenced by polyMesh::findTetFacePt(), interpolationCellPoint< Type >::interpolate(), interpolationCellPointWallModified< Type >::interpolate(), Foam::operator<<(), Foam::operator>>(), wallBoundedStreamLine::track(), and wallBoundedParticle::trackToEdge().

Here is the caller graph for this function:

◆ face() [2/2]

label & face ( )
inlinenoexcept

Non-const access to the face index.

Definition at line 142 of file tetIndices.H.

◆ tetPt() [1/2]

label tetPt ( ) const
inlinenoexcept

Return the characterising tet point index.

Definition at line 145 of file tetIndices.H.

Referenced by polyMesh::findTetFacePt(), Foam::operator<<(), Foam::operator>>(), and wallBoundedStreamLine::track().

Here is the caller graph for this function:

◆ tetPt() [2/2]

label & tetPt ( )
inlinenoexcept

Non-const access to the characterising tet point index.

Definition at line 148 of file tetIndices.H.

◆ faceTriIs()

Foam::triFace faceTriIs ( const polyMesh mesh,
const bool  warn = true 
) const
inline

Return the indices corresponding to the tri on the face for.

this tet. The normal of the tri points out of the cell

Definition at line 54 of file tetIndicesI.H.

References Foam::endl(), f(), polyMesh::faceOwner(), polyMesh::faces(), UList< T >::fcIndex(), mesh, UList< T >::size(), polyMesh::tetBasePtIs(), Foam::Warning, and WarningInFunction.

Referenced by Dual< Type >::add(), Moment< Type >::add(), Dual< Type >::Dual(), cellPointWeight::findTetrahedron(), cellPointWeight::findTriangle(), Dual< Type >::interpolate(), Moment< Type >::interpolate(), Dual< Type >::interpolateGrad(), and Moment< Type >::Moment().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ triIs()

Foam::triFace triIs ( const polyMesh mesh,
const bool  warn = true 
) const
inline

Return the local indices corresponding to the tri on the face for this tet. The normal of the tri points out of the cell

Definition at line 94 of file tetIndicesI.H.

References Foam::endl(), f(), polyMesh::faceOwner(), polyMesh::faces(), UList< T >::fcIndex(), mesh, UList< T >::size(), polyMesh::tetBasePtIs(), Foam::Warning, and WarningInFunction.

Referenced by wallBoundedParticle::trackToEdge().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ tet()

Foam::tetPointRef tet ( const polyMesh mesh) const
inline

Return the geometry corresponding to this tet.

Definition at line 134 of file tetIndicesI.H.

References primitiveMesh::cellCentres(), mesh, polyMesh::points(), and tetIndices::tet().

Referenced by nearWallFields::calcAddressing(), Dual< Type >::Dual(), polyMeshTetDecomposition::findTet(), cellPointWeight::findTetrahedron(), interpolation< Type >::interpolate(), Moment< Type >::Moment(), tetIndices::tet(), and wallBoundedParticle::trackToEdge().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ faceTri()

Foam::triPointRef faceTri ( const polyMesh mesh) const
inline

Return the geometry corresponding to the tri on the face for this tet. The normal of the tri points out of the cell

Definition at line 149 of file tetIndicesI.H.

References mesh, and polyMesh::points().

Referenced by cellPointWeight::findTriangle(), FreeStream< CloudType >::inflow(), polyMesh::pointInCell(), wallBoundedStreamLine::track(), and wallBoundedParticle::trackToEdge().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ oldFaceTri()

Foam::triPointRef oldFaceTri ( const polyMesh mesh) const
inline

Return the geometry corresponding to the tri on the face for this tet using the old positions

Definition at line 163 of file tetIndicesI.H.

References mesh, and polyMesh::oldPoints().

Here is the call graph for this function:

◆ compare()

Foam::label compare ( const tetIndices a,
const tetIndices b 
)
inlinestaticnoexcept

Compare tetIndices for equality. Compares cell, face, tetPt elements in order, stopping at the first inequality.

Returns
negative/zero/positive from the last element compared

Definition at line 180 of file tetIndicesI.H.

References b, and Foam::diff().

Referenced by Foam::operator<(), Foam::operator<=(), Foam::operator>(), and Foam::operator>=().

Here is the call graph for this function:
Here is the caller graph for this function:

Friends And Related Function Documentation

◆ operator>>

Istream & operator>> ( Istream ,
tetIndices  
)
friend

◆ operator<<

Ostream & operator<< ( Ostream ,
const tetIndices  
)
friend

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