tetrahedron< Point, PointRef > Class Template Reference

A tetrahedron primitive. More...

Classes

class  dummyOp
 Dummy. More...
 
class  storeOp
 Store resulting tets. More...
 
class  sumVolOp
 Sum resulting volumes. More...
 

Public Types

enum  { nVertices = 4, nEdges = 6 }
 
typedef Point point_type
 The point type. More...
 
typedef FixedList< tetPoints, 200 > tetIntersectionList
 Storage type for tets originating from intersecting tets. More...
 

Public Member Functions

 tetrahedron (const Point &p0, const Point &p1, const Point &p2, const Point &p3)
 Construct from four points. More...
 
 tetrahedron (const FixedList< Point, 4 > &pts)
 Construct from four points. More...
 
 tetrahedron (const UList< Point > &points, const FixedList< label, 4 > &indices)
 Construct from four points in the list of points. More...
 
 tetrahedron (Istream &)
 Construct from Istream. More...
 
const Pointa () const noexcept
 Return vertex a. More...
 
const Pointb () const noexcept
 Return vertex b. More...
 
const Pointc () const noexcept
 Return vertex c. More...
 
const Pointd () const noexcept
 Return vertex d. More...
 
triPointRef tri (const label facei) const
 Return i-th face. More...
 
vector Sa () const
 Face area normal for side a. More...
 
vector Sb () const
 Face area normal for side b. More...
 
vector Sc () const
 Face area normal for side c. More...
 
vector Sd () const
 Face area normal for side d. More...
 
Point centre () const
 Return centre (centroid) More...
 
scalar mag () const
 Return volume. More...
 
Pair< Pointbox () const
 The enclosing (bounding) box for the tetrahedron. More...
 
treeBoundBox bounds () const
 The bounding box for the tetrahedron. More...
 
Point circumCentre () const
 Return circum-centre. More...
 
scalar circumRadius () const
 Return circum-radius. More...
 
scalar quality () const
 Return quality: Ratio of tetrahedron and circum-sphere volume, scaled so that a regular tetrahedron has a quality of 1. More...
 
Point randomPoint (Random &rndGen) const
 Return a random point in the tetrahedron from a uniform distribution. More...
 
Point barycentricToPoint (const barycentric &bary) const
 Calculate the point from the given barycentric coordinates. More...
 
barycentric pointToBarycentric (const point &pt) const
 Calculate the barycentric coordinates from the given point. More...
 
scalar pointToBarycentric (const point &pt, barycentric &bary) const
 Calculate the barycentric coordinates from the given point. More...
 
pointHit nearestPoint (const point &p) const
 Return nearest point to p on tetrahedron. Is p itself. More...
 
bool inside (const point &pt) const
 Return true if point is inside tetrahedron. More...
 
template<class AboveTetOp , class BelowTetOp >
void sliceWithPlane (const plane &pl, AboveTetOp &aboveOp, BelowTetOp &belowOp) const
 Decompose tet into tets above and below plane. More...
 
void tetOverlap (const tetrahedron< Point, PointRef > &tetB, tetIntersectionList &insideTets, label &nInside, tetIntersectionList &outsideTets, label &nOutside) const
 Decompose tet into tets inside and outside other tet. More...
 
pointHit containmentSphere (const scalar tol) const
 Return (min)containment sphere, i.e. the smallest sphere with. More...
 
void gradNiSquared (scalarField &buffer) const
 Fill buffer with shape function products. More...
 
void gradNiDotGradNj (scalarField &buffer) const
 
void gradNiGradNi (tensorField &buffer) const
 
void gradNiGradNj (tensorField &buffer) const
 

Static Public Member Functions

static Pair< Pointbox (const Point &p0, const Point &p1, const Point &p2, const Point &p3)
 The enclosing (bounding) box for four points. More...
 

Friends

Istreamoperator>> (Istream &, tetrahedron &)
 
Ostreamoperator (Ostream &, const tetrahedron &)
 

Detailed Description

template<class Point, class PointRef>
class Foam::tetrahedron< Point, PointRef >

A tetrahedron primitive.

Ordering of edges needs to be the same for a tetrahedron class, a tetrahedron cell shape model and a tetCell.

Source files

Definition at line 58 of file tetrahedron.H.

Member Typedef Documentation

◆ point_type

typedef Point point_type

The point type.

Definition at line 217 of file tetrahedron.H.

◆ tetIntersectionList

Storage type for tets originating from intersecting tets.

(can possibly be smaller than 200)

Definition at line 224 of file tetrahedron.H.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
nVertices 
nEdges 

Definition at line 304 of file tetrahedron.H.

Constructor & Destructor Documentation

◆ tetrahedron() [1/4]

tetrahedron ( const Point p0,
const Point p1,
const Point p2,
const Point p3 
)
inline

Construct from four points.

Definition at line 70 of file tetrahedronI.H.

◆ tetrahedron() [2/4]

tetrahedron ( const FixedList< Point, 4 > &  pts)
inline

Construct from four points.

Definition at line 86 of file tetrahedronI.H.

◆ tetrahedron() [3/4]

tetrahedron ( const UList< Point > &  points,
const FixedList< label, 4 > &  indices 
)
inline

Construct from four points in the list of points.

Definition at line 99 of file tetrahedronI.H.

◆ tetrahedron() [4/4]

tetrahedron ( Istream is)
inlineexplicit

Construct from Istream.

Definition at line 112 of file tetrahedronI.H.

Member Function Documentation

◆ a()

const Point& a ( ) const
inlinenoexcept

Return vertex a.

Definition at line 351 of file tetrahedron.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ b()

const Point& b ( ) const
inlinenoexcept

Return vertex b.

Definition at line 356 of file tetrahedron.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ c()

const Point& c ( ) const
inlinenoexcept

Return vertex c.

Definition at line 361 of file tetrahedron.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ d()

const Point& d ( ) const
inlinenoexcept

Return vertex d.

Definition at line 366 of file tetrahedron.H.

Referenced by Foam::operator<<().

Here is the caller graph for this function:

◆ tri()

Foam::triPointRef tri ( const label  facei) const
inline

Return i-th face.

Definition at line 185 of file tetrahedronI.H.

References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ box() [1/2]

Foam::Pair< Point > box ( const Point p0,
const Point p1,
const Point p2,
const Point p3 
)
inlinestatic

The enclosing (bounding) box for four points.

Definition at line 139 of file tetrahedronI.H.

References Foam::max(), Foam::min(), and p0.

Here is the call graph for this function:

◆ Sa()

Foam::vector Sa ( ) const
inline

Face area normal for side a.

Definition at line 217 of file tetrahedronI.H.

References triangle< Point, PointRef >::areaNormal().

Here is the call graph for this function:

◆ Sb()

Foam::vector Sb ( ) const
inline

Face area normal for side b.

Definition at line 224 of file tetrahedronI.H.

References triangle< Point, PointRef >::areaNormal().

Here is the call graph for this function:

◆ Sc()

Foam::vector Sc ( ) const
inline

Face area normal for side c.

Definition at line 231 of file tetrahedronI.H.

References triangle< Point, PointRef >::areaNormal().

Here is the call graph for this function:

◆ Sd()

Foam::vector Sd ( ) const
inline

Face area normal for side d.

Definition at line 238 of file tetrahedronI.H.

References triangle< Point, PointRef >::areaNormal().

Here is the call graph for this function:

◆ centre()

Point centre ( ) const
inline

Return centre (centroid)

Definition at line 245 of file tetrahedronI.H.

Referenced by nearWallFields::calcAddressing().

Here is the caller graph for this function:

◆ mag()

Foam::scalar mag ( ) const
inline

◆ box() [2/2]

Foam::Pair< Point > box ( ) const
inline

The enclosing (bounding) box for the tetrahedron.

Definition at line 155 of file tetrahedronI.H.

Referenced by tetPoints::box().

Here is the caller graph for this function:

◆ bounds()

Foam::treeBoundBox bounds ( ) const
inline

The bounding box for the tetrahedron.

Definition at line 168 of file tetrahedronI.H.

◆ circumCentre()

Point circumCentre ( ) const
inline

Return circum-centre.

Definition at line 259 of file tetrahedronI.H.

References Foam::constant::physicoChemical::b, Foam::constant::universal::c, lambda(), Foam::mag(), Foam::magSqr(), and Foam::constant::physicoChemical::mu.

Here is the call graph for this function:

◆ circumRadius()

Foam::scalar circumRadius ( ) const
inline

Return circum-radius.

Definition at line 286 of file tetrahedronI.H.

References Foam::constant::physicoChemical::b, Foam::constant::universal::c, lambda(), Foam::mag(), Foam::magSqr(), and Foam::constant::physicoChemical::mu.

Here is the call graph for this function:

◆ quality()

Foam::scalar quality ( ) const
inline

Return quality: Ratio of tetrahedron and circum-sphere volume, scaled so that a regular tetrahedron has a quality of 1.

Definition at line 312 of file tetrahedronI.H.

References Foam::mag(), Foam::min(), Foam::pow3(), and Foam::sqrt().

Referenced by polyMeshTetDecomposition::minQuality().

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

◆ randomPoint()

Point randomPoint ( Random rndGen) const
inline

Return a random point in the tetrahedron from a uniform distribution.

Definition at line 326 of file tetrahedronI.H.

References Foam::barycentric01(), and rndGen.

Here is the call graph for this function:

◆ barycentricToPoint()

Point barycentricToPoint ( const barycentric bary) const
inline

Calculate the point from the given barycentric coordinates.

Definition at line 336 of file tetrahedronI.H.

References Barycentric< Cmpt >::a(), Barycentric< Cmpt >::b(), Barycentric< Cmpt >::c(), and Barycentric< Cmpt >::d().

Referenced by interpolation< Foam::vector >::interpolate().

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

◆ pointToBarycentric() [1/2]

Foam::barycentric pointToBarycentric ( const point pt) const
inline

Calculate the barycentric coordinates from the given point.

Definition at line 346 of file tetrahedronI.H.

Referenced by cellPointWeight::findTetrahedron().

Here is the caller graph for this function:

◆ pointToBarycentric() [2/2]

Foam::scalar pointToBarycentric ( const point pt,
barycentric bary 
) const
inline

Calculate the barycentric coordinates from the given point.

Returns the determinant.

Definition at line 358 of file tetrahedronI.H.

References Foam::det(), Foam::inv(), and Foam::mag().

Here is the call graph for this function:

◆ nearestPoint()

Foam::pointHit nearestPoint ( const point p) const
inline

Return nearest point to p on tetrahedron. Is p itself.

if inside.

Definition at line 401 of file tetrahedronI.H.

References triangle< Point, PointRef >::areaNormal(), PointHit< PointType >::distance(), triangle< Point, PointRef >::nearestPoint(), p, and PointHit< PointType >::point().

Referenced by cellPointWeight::findTetrahedron().

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

◆ inside()

bool inside ( const point pt) const
inline

Return true if point is inside tetrahedron.

Definition at line 515 of file tetrahedronI.H.

References p, and triangle< Point, PointRef >::unitNormal().

Referenced by polyMeshTetDecomposition::findTet().

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

◆ sliceWithPlane()

void sliceWithPlane ( const plane pl,
AboveTetOp &  aboveOp,
BelowTetOp &  belowOp 
) const
inline

Decompose tet into tets above and below plane.

Definition at line 1108 of file tetrahedronI.H.

◆ tetOverlap()

void tetOverlap ( const tetrahedron< Point, PointRef > &  tetB,
tetIntersectionList insideTets,
label &  nInside,
tetIntersectionList outsideTets,
label &  nOutside 
) const
inline

Decompose tet into tets inside and outside other tet.

◆ containmentSphere()

Foam::pointHit containmentSphere ( const scalar  tol) const

Return (min)containment sphere, i.e. the smallest sphere with.

all points inside. Returns pointHit with:

  • hit : if sphere is equal to circumsphere (biggest sphere)
  • point : centre of sphere
  • distance : radius of sphere
  • eligiblemiss: false Tol (small compared to 1, e.g. 1e-9) is used to determine whether point is inside: mag(pt - ctr) < (1+tol)*radius.

Definition at line 28 of file tetrahedron.C.

References Foam::magSqr(), PointHit< PointType >::point(), PointHit< PointType >::setDistance(), PointHit< PointType >::setHit(), PointHit< PointType >::setMiss(), PointHit< PointType >::setPoint(), and Foam::sqrt().

Here is the call graph for this function:

◆ gradNiSquared()

void gradNiSquared ( scalarField buffer) const

Fill buffer with shape function products.

Definition at line 240 of file tetrahedron.C.

References Foam::mag(), and Foam::magSqr().

Here is the call graph for this function:

◆ gradNiDotGradNj()

void gradNiDotGradNj ( scalarField buffer) const

Definition at line 261 of file tetrahedron.C.

References Foam::mag().

Here is the call graph for this function:

◆ gradNiGradNi()

void gradNiGradNi ( tensorField buffer) const

Definition at line 291 of file tetrahedron.C.

References Foam::mag(), and Foam::sqr().

Here is the call graph for this function:

◆ gradNiGradNj()

void gradNiGradNj ( tensorField buffer) const

Definition at line 309 of file tetrahedron.C.

References Foam::mag().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator>>

Istream& operator>> ( Istream ,
tetrahedron< Point, PointRef > &   
)
friend

◆ operator

Ostream& operator ( Ostream ,
const tetrahedron< Point, PointRef > &   
)
friend

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