Loading...
Searching...
No Matches
tetrahedron< Point, PointRef > Class Template Reference

A tetrahedron primitive. More...

#include <tetrahedron.H>

Inheritance diagram for tetrahedron< Point, PointRef >:

Classes

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

Public Types

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

Public Member Functions

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

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 208 of file tetrahedron.H.

Member Typedef Documentation

◆ point_type

template<class Point, class PointRef>
typedef Point point_type

The point type.

Definition at line 217 of file tetrahedron.H.

◆ tetIntersectionList

template<class Point, class PointRef>
typedef FixedList<tetPoints, 200> 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

template<class Point, class PointRef>
anonymous enum
Enumerator
nVertices 
nEdges 

Definition at line 304 of file tetrahedron.H.

Constructor & Destructor Documentation

◆ tetrahedron() [1/4]

template<class Point, class PointRef>
tetrahedron ( const Point & p0,
const Point & p1,
const Point & p2,
const Point & p3 )
inline

Construct from four points.

Definition at line 69 of file tetrahedronI.H.

References p0.

◆ tetrahedron() [2/4]

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

Construct from four points.

Definition at line 85 of file tetrahedronI.H.

References pts.

◆ tetrahedron() [3/4]

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

Construct from four points in the list of points.

Definition at line 98 of file tetrahedronI.H.

References points.

◆ tetrahedron() [4/4]

template<class Point, class PointRef>
tetrahedron ( Istream & is)
inlineexplicit

Construct from Istream.

Definition at line 112 of file tetrahedronI.H.

Member Function Documentation

◆ a()

template<class Point, class PointRef>
const Point & a ( ) const
inlinenoexcept

Return vertex a.

Definition at line 351 of file tetrahedron.H.

Referenced by circumRadius(), and Foam::operator<<().

Here is the caller graph for this function:

◆ b()

template<class Point, class PointRef>
const Point & b ( ) const
inlinenoexcept

Return vertex b.

Definition at line 356 of file tetrahedron.H.

Referenced by circumRadius(), and Foam::operator<<().

Here is the caller graph for this function:

◆ c()

template<class Point, class PointRef>
const Point & c ( ) const
inlinenoexcept

Return vertex c.

Definition at line 361 of file tetrahedron.H.

Referenced by circumRadius(), and Foam::operator<<().

Here is the caller graph for this function:

◆ d()

template<class Point, class PointRef>
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()

template<class Point, class PointRef>
Foam::triPointRef tri ( const label facei) const
inline

Return i-th face.

Definition at line 184 of file tetrahedronI.H.

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

Here is the call graph for this function:

◆ box() [1/2]

template<class Point, class PointRef>
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 138 of file tetrahedronI.H.

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

Referenced by bounds().

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

◆ Sa()

template<class Point, class PointRef>
Foam::vector Sa ( ) const
inline

Face area normal for side a.

Definition at line 217 of file tetrahedronI.H.

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

Referenced by gradNiDotGradNj(), gradNiGradNi(), gradNiGradNj(), and gradNiSquared().

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

◆ Sb()

template<class Point, class PointRef>
Foam::vector Sb ( ) const
inline

Face area normal for side b.

Definition at line 224 of file tetrahedronI.H.

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

Referenced by gradNiDotGradNj(), gradNiGradNi(), gradNiGradNj(), and gradNiSquared().

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

◆ Sc()

template<class Point, class PointRef>
Foam::vector Sc ( ) const
inline

Face area normal for side c.

Definition at line 231 of file tetrahedronI.H.

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

Referenced by gradNiDotGradNj(), gradNiGradNi(), gradNiGradNj(), and gradNiSquared().

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

◆ Sd()

template<class Point, class PointRef>
Foam::vector Sd ( ) const
inline

Face area normal for side d.

Definition at line 238 of file tetrahedronI.H.

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

Referenced by gradNiDotGradNj(), gradNiGradNi(), gradNiGradNj(), and gradNiSquared().

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

◆ centre()

template<class Point, class PointRef>
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()

template<class Point, class PointRef>
Foam::scalar mag ( ) const
inline

Return volume.

Definition at line 252 of file tetrahedronI.H.

Referenced by Dual< Type >::Dual(), gradNiDotGradNj(), gradNiGradNi(), gradNiGradNj(), gradNiSquared(), and Moment< Type >::Moment().

Here is the caller graph for this function:

◆ box() [2/2]

template<class Point, class PointRef>
Foam::Pair< Point > box ( ) const
inline

The enclosing (bounding) box for the tetrahedron.

Definition at line 155 of file tetrahedronI.H.

References box().

Referenced by box().

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

◆ bounds()

template<class Point, class PointRef>
Foam::treeBoundBox bounds ( ) const
inline

The bounding box for the tetrahedron.

Definition at line 168 of file tetrahedronI.H.

References box().

Here is the call graph for this function:

◆ circumCentre()

template<class Point, class PointRef>
Point circumCentre ( ) const
inline

Return circum-centre.

Definition at line 259 of file tetrahedronI.H.

Referenced by containmentSphere().

Here is the caller graph for this function:

◆ circumRadius()

template<class Point, class PointRef>
Foam::scalar circumRadius ( ) const
inline

Return circum-radius.

Definition at line 286 of file tetrahedronI.H.

References a(), b(), c(), lambda(), Foam::mag(), Foam::magSqr(), and mu.

Referenced by quality().

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

◆ quality()

template<class Point, class PointRef>
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 circumRadius(), Foam::mag(), Foam::min(), Foam::pow3(), and Foam::sqrt().

Referenced by polyMeshTetDecomposition::checkFaceTets(), isoSurfaceTopo::isoSurfaceTopo(), and polyMeshTetDecomposition::minQuality().

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

◆ randomPoint()

template<class Point, class PointRef>
Point randomPoint ( Random & rndGen) const
inline

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

Definition at line 325 of file tetrahedronI.H.

References Foam::barycentric01(), barycentricToPoint(), and rndGen.

Here is the call graph for this function:

◆ barycentricToPoint()

template<class Point, class PointRef>
Point barycentricToPoint ( const barycentric & bary) const
inline

Calculate the point from the given barycentric coordinates.

Definition at line 335 of file tetrahedronI.H.

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

Referenced by interpolation< Type >::interpolate(), and randomPoint().

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

◆ pointToBarycentric() [1/2]

template<class Point, class PointRef>
Foam::barycentric pointToBarycentric ( const point & pt) const
inline

Calculate the barycentric coordinates from the given point.

Definition at line 345 of file tetrahedronI.H.

References pointToBarycentric().

Referenced by cellPointWeight::findTetrahedron(), and pointToBarycentric().

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

◆ pointToBarycentric() [2/2]

template<class Point, class PointRef>
Foam::scalar pointToBarycentric ( const point & pt,
barycentric & bary ) const
inline

Calculate the barycentric coordinates from the given point.

Returns the determinant.

Definition at line 357 of file tetrahedronI.H.

References Foam::det(), Foam::inv(), Foam::mag(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().

Here is the call graph for this function:

◆ nearestPoint()

template<class Point, class PointRef>
Foam::pointHit nearestPoint ( const point & p) const
inline

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

if inside.

Definition at line 400 of file tetrahedronI.H.

References triangle< Point, PointRef >::areaNormal(), PointHit< PointType >::distance(), inside(), 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()

template<class Point, class PointRef>
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(), and nearestPoint().

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

◆ sliceWithPlane()

template<class Point, class PointRef>
template<class AboveTetOp, class BelowTetOp>
void sliceWithPlane ( const plane & pl,
AboveTetOp & aboveOp,
BelowTetOp & belowOp ) const
inline

Decompose tet into tets above and below plane.

Definition at line 1092 of file tetrahedronI.H.

◆ tetOverlap()

template<class Point, class PointRef>
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()

template<class Point, class PointRef>
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 27 of file tetrahedron.C.

References circumCentre(), 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()

template<class Point, class PointRef>
void gradNiSquared ( scalarField & buffer) const

Fill buffer with shape function products.

Definition at line 239 of file tetrahedron.C.

References Foam::mag(), mag(), Foam::magSqr(), Sa(), Sb(), Sc(), and Sd().

Here is the call graph for this function:

◆ gradNiDotGradNj()

template<class Point, class PointRef>
void gradNiDotGradNj ( scalarField & buffer) const

Definition at line 260 of file tetrahedron.C.

References Foam::mag(), mag(), Sa(), Sb(), Sc(), and Sd().

Here is the call graph for this function:

◆ gradNiGradNi()

template<class Point, class PointRef>
void gradNiGradNi ( tensorField & buffer) const

Definition at line 290 of file tetrahedron.C.

References Foam::mag(), mag(), Sa(), Sb(), Sc(), Sd(), and Foam::sqr().

Here is the call graph for this function:

◆ gradNiGradNj()

template<class Point, class PointRef>
void gradNiGradNj ( tensorField & buffer) const

Definition at line 308 of file tetrahedron.C.

References Foam::mag(), mag(), Sa(), Sb(), Sc(), and Sd().

Here is the call graph for this function:

◆ operator>>

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

◆ operator

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

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