A tetrahedron primitive. More...
#include <tetrahedron.H>

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 Point & | a () const noexcept |
| Return vertex a. | |
| const Point & | b () const noexcept |
| Return vertex b. | |
| const Point & | c () const noexcept |
| Return vertex c. | |
| const Point & | d () 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< Point > | box () 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< Point > | box (const Point &p0, const Point &p1, const Point &p2, const Point &p3) |
| The enclosing (bounding) box for four points. | |
Friends | |
| Istream & | operator>> (Istream &, tetrahedron &) |
| Ostream & | operator (Ostream &, const tetrahedron &) |
A tetrahedron primitive.
Ordering of edges needs to be the same for a tetrahedron class, a tetrahedron cell shape model and a tetCell.
Definition at line 208 of file tetrahedron.H.
| typedef Point point_type |
The point type.
Definition at line 217 of file tetrahedron.H.
| 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.
| anonymous enum |
| Enumerator | |
|---|---|
| nVertices | |
| nEdges | |
Definition at line 304 of file tetrahedron.H.
|
inline |
Construct from four points in the list of points.
Definition at line 98 of file tetrahedronI.H.
References points.
Construct from Istream.
Definition at line 112 of file tetrahedronI.H.
Return vertex a.
Definition at line 351 of file tetrahedron.H.
Referenced by circumRadius(), and Foam::operator<<().

Return vertex b.
Definition at line 356 of file tetrahedron.H.
Referenced by circumRadius(), and Foam::operator<<().

Return vertex c.
Definition at line 361 of file tetrahedron.H.
Referenced by circumRadius(), and Foam::operator<<().

Return vertex d.
Definition at line 366 of file tetrahedron.H.
Referenced by Foam::operator<<().

|
inline |
Return i-th face.
Definition at line 184 of file tetrahedronI.H.
References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.

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


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


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


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


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


Return centre (centroid).
Definition at line 245 of file tetrahedronI.H.
Referenced by nearWallFields::calcAddressing().

|
inline |
Return volume.
Definition at line 252 of file tetrahedronI.H.
Referenced by Dual< Type >::Dual(), gradNiDotGradNj(), gradNiGradNi(), gradNiGradNj(), gradNiSquared(), and Moment< Type >::Moment().

|
inline |
The enclosing (bounding) box for the tetrahedron.
Definition at line 155 of file tetrahedronI.H.
References box().
Referenced by box().


|
inline |
The bounding box for the tetrahedron.
Definition at line 168 of file tetrahedronI.H.
References box().

Return circum-centre.
Definition at line 259 of file tetrahedronI.H.
Referenced by containmentSphere().

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


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


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.

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


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


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

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


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


|
inline |
Decompose tet into tets above and below plane.
Definition at line 1092 of file tetrahedronI.H.
|
inline |
Decompose tet into tets inside and outside other tet.
| Foam::pointHit containmentSphere | ( | const scalar | tol | ) | const |
Return (min)containment sphere, i.e. the smallest sphere with.
all points inside. Returns pointHit with:
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().

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

| void gradNiDotGradNj | ( | scalarField & | buffer | ) | const |
Definition at line 260 of file tetrahedron.C.
References Foam::mag(), mag(), Sa(), Sb(), Sc(), and Sd().

| void gradNiGradNi | ( | tensorField & | buffer | ) | const |
Definition at line 290 of file tetrahedron.C.
References Foam::mag(), mag(), Sa(), Sb(), Sc(), Sd(), and Foam::sqr().

| void gradNiGradNj | ( | tensorField & | buffer | ) | const |
Definition at line 308 of file tetrahedron.C.
References Foam::mag(), mag(), Sa(), Sb(), Sc(), and Sd().

|
friend |
|
friend |