tetrahedron< Point, PointRef > Class Template Reference

A tetrahedron primitive. More...

Public Types

enum  { nVertices = 4, nEdges = 6 }
 

Public Member Functions

 tetrahedron (const Point &a, const Point &b, const Point &c, const Point &d)
 Construct from points. More...
 
 tetrahedron (const UList< Point > &, const FixedList< label, 4 > &indices)
 Construct from four points in the list of points. More...
 
 tetrahedron (Istream &)
 Construct from Istream. More...
 
const Point & a () const
 Return vertices. More...
 
const Point & b () const
 
const Point & c () const
 
const Point & d () const
 
triPointRef tri (const label facei) const
 Return i-th face. More...
 
vector Sa () const
 Return face normal. More...
 
vector Sb () const
 
vector Sc () const
 
vector Sd () const
 
Point centre () const
 Return centre (centroid) More...
 
scalar mag () const
 Return volume. 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. More...
 
Point randomPoint (Random &rndGen) const
 Return a random point in the tetrahedron from a. 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...
 
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
 
boundBox bounds () const
 Calculate the bounding box. 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 61 of file tetrahedron.H.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
Enumerator
nVertices 
nEdges 

Definition at line 93 of file tetrahedron.H.

Constructor & Destructor Documentation

◆ tetrahedron() [1/3]

tetrahedron ( const Point &  a,
const Point &  b,
const Point &  c,
const Point &  d 
)
inline

Construct from points.

Definition at line 34 of file tetrahedronI.H.

◆ tetrahedron() [2/3]

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

Construct from four points in the list of points.

Definition at line 50 of file tetrahedronI.H.

◆ tetrahedron() [3/3]

tetrahedron ( Istream is)
inline

Construct from Istream.

Definition at line 63 of file tetrahedronI.H.

Member Function Documentation

◆ a()

const Point & a ( ) const
inline

Return vertices.

Definition at line 72 of file tetrahedronI.H.

Referenced by tetOverlapVolume::tetOverlapVolume().

Here is the caller graph for this function:

◆ b()

const Point & b ( ) const
inline

Definition at line 79 of file tetrahedronI.H.

Referenced by tetOverlapVolume::tetOverlapVolume().

Here is the caller graph for this function:

◆ c()

const Point & c ( ) const
inline

Definition at line 86 of file tetrahedronI.H.

Referenced by tetOverlapVolume::tetOverlapVolume().

Here is the caller graph for this function:

◆ d()

const Point & d ( ) const
inline

Definition at line 93 of file tetrahedronI.H.

References tetrahedron< Point, PointRef >::tri().

Referenced by tetOverlapVolume::tetOverlapVolume().

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

◆ tri()

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

Return i-th face.

Definition at line 101 of file tetrahedronI.H.

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

Referenced by tetrahedron< Point, PointRef >::d().

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

◆ Sa()

Foam::vector Sa ( ) const
inline

Return face normal.

Definition at line 135 of file tetrahedronI.H.

◆ Sb()

Foam::vector Sb ( ) const
inline

Definition at line 142 of file tetrahedronI.H.

◆ Sc()

Foam::vector Sc ( ) const
inline

Definition at line 149 of file tetrahedronI.H.

◆ Sd()

Foam::vector Sd ( ) const
inline

Definition at line 156 of file tetrahedronI.H.

◆ centre()

Point centre ( ) const
inline

Return centre (centroid)

Definition at line 163 of file tetrahedronI.H.

◆ mag()

Foam::scalar mag ( ) const
inline

Return volume.

Definition at line 170 of file tetrahedronI.H.

Referenced by Dual< Type >::Dual(), Moment< Type >::Moment(), and AveragingMethod< Foam::Vector >::write().

Here is the caller graph for this function:

◆ circumCentre()

Point circumCentre ( ) const
inline

Return circum-centre.

Definition at line 177 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 204 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 230 of file tetrahedronI.H.

References Foam::mag(), Foam::min(), Foam::pow3(), tetrahedron< Point, PointRef >::randomPoint(), 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 244 of file tetrahedronI.H.

References Foam::barycentric01(), and tetrahedron< Point, PointRef >::barycentricToPoint().

Referenced by tetrahedron< Point, PointRef >::quality().

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

◆ barycentricToPoint()

Point barycentricToPoint ( const barycentric bary) const
inline

Calculate the point from the given barycentric coordinates.

Definition at line 254 of file tetrahedronI.H.

References tetrahedron< Point, PointRef >::pointToBarycentric().

Referenced by interpolation< Foam::Vector >::interpolate(), and tetrahedron< Point, PointRef >::randomPoint().

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 264 of file tetrahedronI.H.

Referenced by tetrahedron< Point, PointRef >::barycentricToPoint(), and 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 276 of file tetrahedronI.H.

References Foam::cmptSum(), Foam::det(), Foam::inv(), Foam::mag(), tetrahedron< Point, PointRef >::nearestPoint(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

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 319 of file tetrahedronI.H.

References PointHit< Point >::distance(), p, and PointHit< Point >::rawPoint().

Referenced by cellPointWeight::findTetrahedron(), and tetrahedron< Point, PointRef >::pointToBarycentric().

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 413 of file tetrahedronI.H.

References IOstream::check(), Foam::mag(), n, Foam::nl, Istream::readBegin(), Istream::readEnd(), and Foam::Zero.

Referenced by polyMeshTetDecomposition::findTet().

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

◆ 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 34 of file tetrahedron.C.

References tetrahedron< Point, PointRef >::gradNiSquared(), Foam::magSqr(), PointHit< Point >::rawPoint(), PointHit< Point >::setDistance(), PointHit< Point >::setHit(), PointHit< Point >::setMiss(), PointHit< Point >::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 246 of file tetrahedron.C.

References tetrahedron< Point, PointRef >::gradNiDotGradNj(), Foam::mag(), and Foam::magSqr().

Referenced by tetrahedron< Point, PointRef >::containmentSphere().

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

◆ gradNiDotGradNj()

void gradNiDotGradNj ( scalarField buffer) const

Definition at line 267 of file tetrahedron.C.

References tetrahedron< Point, PointRef >::gradNiGradNi(), and Foam::mag().

Referenced by tetrahedron< Point, PointRef >::gradNiSquared().

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

◆ gradNiGradNi()

void gradNiGradNi ( tensorField buffer) const

Definition at line 297 of file tetrahedron.C.

References tetrahedron< Point, PointRef >::gradNiGradNj(), Foam::mag(), and Foam::sqr().

Referenced by tetrahedron< Point, PointRef >::gradNiDotGradNj().

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

◆ gradNiGradNj()

void gradNiGradNj ( tensorField buffer) const

Definition at line 315 of file tetrahedron.C.

References Foam::mag().

Referenced by tetrahedron< Point, PointRef >::gradNiGradNi().

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

◆ bounds()

Foam::boundBox bounds ( ) const

Calculate the bounding box.

Definition at line 340 of file tetrahedron.C.

References Foam::constant::physicoChemical::b, Foam::constant::universal::c, Foam::max(), and Foam::min().

Referenced by tetOverlapVolume::cellCellOverlapMinDecomp(), and tetOverlapVolume::cellCellOverlapVolumeMinDecomp().

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 ,
tetrahedron< Point, PointRef > &   
)
friend

◆ operator

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

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