Loading...
Searching...
No Matches
tetIndices Class Reference

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

#include <tetIndices.H>

Public Member Functions

constexpr tetIndices () noexcept
 Default construct, with invalid labels (-1).
constexpr tetIndices (label celli, label facei, label tetPointi) noexcept
 Construct from components.
 ~tetIndices ()=default
 Destructor.
label cell () const noexcept
 Return the cell index.
label & cell () noexcept
 Non-const access to the cell index.
label face () const noexcept
 Return the face index.
label & face () noexcept
 Non-const access to the face index.
label tetPt () const noexcept
 Return the characterising tet point index.
label & tetPt () noexcept
 Non-const access to the characterising tet point index.
triFace faceTriIs (const polyMesh &mesh, const bool warn=true) const
 The indices corresponding to the tri on the face for this tet. The normal of the tri points out of the cell.
triFace triIs (const polyMesh &mesh, const bool warn=true) const
 Local indices corresponding to the tri on the face for this tet. The normal of the tri points out of the cell.
tetPointRef tet (const polyMesh &mesh) const
 The tet geometry for this tet, where point0 is the cell centre.
tetPointRef oldTet (const polyMesh &mesh) const
 The tet geometry for this tet (using old positions), where point0 is the cell centre.
triPointRef faceTri (const polyMesh &mesh) const
 The triangle geometry for the face for this tet. The normal of the tri points out of the cell.
triPointRef oldFaceTri (const polyMesh &mesh) const
 The triangle geometry for the face for this tet (using old positions).
point barycentricToPoint (const polyMesh &mesh, const barycentric &bary) const
 The x/y/z position for given barycentric coordinates (where point0 is the cell centre).

Static Public Member Functions

static label compare (const tetIndices &a, const tetIndices &b) noexcept
 Compare tetIndices for equality. Compares cell, face, tetPt elements in order, stopping at the first inequality.

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 78 of file tetIndices.H.

Constructor & Destructor Documentation

◆ tetIndices() [1/2]

tetIndices ( )
inlineconstexprnoexcept

Default construct, with invalid labels (-1).

Definition at line 24 of file tetIndicesI.H.

References Foam::noexcept.

Referenced by compare(), operator<<, and operator>>.

Here is the caller graph for this function:

◆ tetIndices() [2/2]

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

Construct from components.

Definition at line 32 of file tetIndicesI.H.

References Foam::noexcept.

◆ ~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 151 of file tetIndices.H.

References Foam::noexcept.

◆ face() [1/2]

label face ( ) const
inlinenoexcept

◆ face() [2/2]

label & face ( )
inlinenoexcept

Non-const access to the face index.

Definition at line 161 of file tetIndices.H.

References Foam::noexcept.

◆ tetPt() [1/2]

label tetPt ( ) const
inlinenoexcept

Return the characterising tet point index.

Definition at line 166 of file tetIndices.H.

References Foam::noexcept.

Referenced by polyMesh::findTetFacePt(), 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 171 of file tetIndices.H.

References Foam::noexcept.

◆ faceTriIs()

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

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 47 of file tetIndicesI.H.

References Foam::endl(), f(), mesh, triFace(), Foam::Warning, and WarningInFunction.

Referenced by Dual< Type >::add(), Moment< Type >::add(), barycentricToPoint(), Dual< Type >::Dual(), faceTri(), cellPointWeight::findTetrahedron(), cellPointWeight::findTriangle(), Dual< Type >::interpolate(), Moment< Type >::interpolate(), interpolationCellPoint< Type >::interpolate(), Dual< Type >::interpolateGrad(), Moment< Type >::Moment(), oldFaceTri(), oldTet(), and tet().

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

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 88 of file tetIndicesI.H.

References Foam::endl(), f(), mesh, triFace(), 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

The tet geometry for this tet, where point0 is the cell centre.

Definition at line 129 of file tetIndicesI.H.

References faceTriIs(), mesh, and pts.

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

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

◆ oldTet()

Foam::tetPointRef oldTet ( const polyMesh & mesh) const
inline

The tet geometry for this tet (using old positions), where point0 is the cell centre.

Definition at line 144 of file tetIndicesI.H.

References faceTriIs(), mesh, and pts.

Here is the call graph for this function:

◆ faceTri()

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

The triangle geometry for the face for this tet. The normal of the tri points out of the cell.

Definition at line 182 of file tetIndicesI.H.

References faceTriIs(), and mesh.

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

The triangle geometry for the face for this tet (using old positions).

Definition at line 188 of file tetIndicesI.H.

References faceTriIs(), and mesh.

Here is the call graph for this function:

◆ barycentricToPoint()

Foam::point barycentricToPoint ( const polyMesh & mesh,
const barycentric & bary ) const
inline

The x/y/z position for given barycentric coordinates (where point0 is the cell centre).

Definition at line 159 of file tetIndicesI.H.

References Barycentric< Cmpt >::a(), b, Barycentric< Cmpt >::b(), Barycentric< Cmpt >::c(), Barycentric< Cmpt >::d(), faceTriIs(), mesh, pts, Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

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 197 of file tetIndicesI.H.

References b, Foam::diff(), and tetIndices().

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:

◆ operator>>

Istream & operator>> ( Istream & ,
tetIndices &  )
friend

References tetIndices().

◆ operator<<

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

References tetIndices().


The documentation for this class was generated from the following files:
  • src/OpenFOAM/meshes/polyMesh/polyMeshTetDecomposition/tetIndices.H
  • src/OpenFOAM/meshes/polyMesh/polyMeshTetDecomposition/tetIndices.C
  • src/OpenFOAM/meshes/polyMesh/polyMeshTetDecomposition/tetIndicesI.H