Classes | Public Types | Public Member Functions | Friends | List of all members
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 FixedList< tetPoints, 200 > tetIntersectionList
 Storage type for tets originating from intersecting tets. More...
 

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 randomPoint (cachedRandom &rndGen) const
 Return a random point in the tetrahedron from a. More...
 
scalar barycentric (const point &pt, List< scalar > &bary) const
 Calculate the barycentric coordinates of the given. 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
 

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

Member Typedef Documentation

Storage type for tets originating from intersecting tets.

(can possibly be smaller than 200)

Definition at line 93 of file tetrahedron.H.

Member Enumeration Documentation

anonymous enum
Enumerator
nVertices 
nEdges 

Definition at line 164 of file tetrahedron.H.

Constructor & Destructor Documentation

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

Construct from points.

Definition at line 35 of file tetrahedronI.H.

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

Construct from four points in the list of points.

Definition at line 51 of file tetrahedronI.H.

tetrahedron ( Istream is)
inline

Construct from Istream.

Definition at line 64 of file tetrahedronI.H.

Member Function Documentation

const Point & a ( ) const
inline

Return vertices.

Definition at line 73 of file tetrahedronI.H.

Referenced by particle< Type >::movingTetLambda().

Here is the caller graph for this function:

const Point & b ( ) const
inline

Definition at line 80 of file tetrahedronI.H.

Referenced by particle< Type >::movingTetLambda().

Here is the caller graph for this function:

const Point & c ( ) const
inline

Definition at line 87 of file tetrahedronI.H.

Referenced by particle< Type >::movingTetLambda().

Here is the caller graph for this function:

const Point & d ( ) const
inline

Definition at line 94 of file tetrahedronI.H.

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

Referenced by particle< Type >::movingTetLambda(), and tetrahedron< Point, PointRef >::storeOp::operator()().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Return i-th face.

Definition at line 102 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:

Foam::vector Sa ( ) const
inline

Return face normal.

Definition at line 136 of file tetrahedronI.H.

Referenced by particle< Type >::trackToFace().

Here is the caller graph for this function:

Foam::vector Sb ( ) const
inline

Definition at line 143 of file tetrahedronI.H.

Referenced by particle< Type >::trackToFace().

Here is the caller graph for this function:

Foam::vector Sc ( ) const
inline

Definition at line 150 of file tetrahedronI.H.

Referenced by particle< Type >::trackToFace().

Here is the caller graph for this function:

Foam::vector Sd ( ) const
inline

Definition at line 157 of file tetrahedronI.H.

Referenced by particle< Type >::trackToFace().

Here is the caller graph for this function:

Point centre ( ) const
inline

Return centre (centroid)

Definition at line 164 of file tetrahedronI.H.

Referenced by particle< Type >::findTris(), particle< Type >::hitWallFaces(), and particle< Type >::trackToFace().

Here is the caller graph for this function:

Foam::scalar mag ( ) const
inline

Return volume.

Definition at line 171 of file tetrahedronI.H.

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

Here is the caller graph for this function:

Point circumCentre ( ) const
inline

Return circum-centre.

Definition at line 178 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:

Foam::scalar circumRadius ( ) const
inline

Return circum-radius.

Definition at line 205 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:

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

References Foam::mag(), Foam::min(), Foam::pow3(), tetrahedron< Point, PointRef >::randomPoint(), and Foam::sqrt().

Referenced by polyMeshTetDecomposition::findBasePoint(), and polyMeshTetDecomposition::findSharedBasePoint().

Here is the call graph for this function:

Here is the caller graph for this function:

Point randomPoint ( Random rndGen) const
inline

Return a random point in the tetrahedron from a.

uniform distribution

Definition at line 245 of file tetrahedronI.H.

References s(), and Random::scalar01().

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

Here is the call graph for this function:

Here is the caller graph for this function:

Point randomPoint ( cachedRandom rndGen) const
inline

Return a random point in the tetrahedron from a.

uniform distribution

Definition at line 281 of file tetrahedronI.H.

References tetrahedron< Point, PointRef >::barycentric(), s(), and cachedRandom::sample01().

Here is the call graph for this function:

Foam::scalar barycentric ( const point pt,
List< scalar > &  bary 
) const
inline

Calculate the barycentric coordinates of the given.

point, in the same order as a, b, c, d. Returns the determinant of the solution.

Definition at line 317 of file tetrahedronI.H.

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

Referenced by cellPointWeight::findTetrahedron(), interpolationCellPoint< Type >::interpolate(), tetrahedron< Point, PointRef >::randomPoint(), and Dual< Type >::~Dual().

Here is the call graph for this function:

Here is the caller graph for this function:

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

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

if inside.

Definition at line 362 of file tetrahedronI.H.

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

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

Here is the call graph for this function:

Here is the caller graph for this function:

bool inside ( const point pt) const
inline

Return true if point is inside tetrahedron.

Definition at line 456 of file tetrahedronI.H.

References Foam::mag(), n, and Foam::Zero.

Referenced by polyMeshTetDecomposition::findTet().

Here is the call graph for this function:

Here is the caller graph for this function:

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

Decompose tet into tets above and below plane.

Definition at line 985 of file tetrahedronI.H.

References token::BEGIN_LIST, token::END_LIST, Foam::nl, Istream::readBegin(), and token::SPACE.

Referenced by tetrahedron< Point, PointRef >::storeOp::operator()(), and tetOverlapVolume::tetOverlapVolume().

Here is the call graph for this function:

Here is the caller graph for this function:

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.

Definition at line 34 of file tetrahedron.C.

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

Here is the call graph for this function:

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 124 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().

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

Here is the call graph for this function:

Here is the caller graph for this function:

void gradNiSquared ( scalarField buffer) const

Fill buffer with shape function products.

Definition at line 336 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:

void gradNiDotGradNj ( scalarField buffer) const

Definition at line 357 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:

void gradNiGradNi ( tensorField buffer) const

Definition at line 387 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:

void gradNiGradNj ( tensorField buffer) const

Definition at line 405 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:

Friends And Related Function Documentation

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

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