A tetrahedron primitive. More...
Data Structures | |
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 |
Static Private Member Functions | |
static point | planeIntersection (const FixedList< scalar, 4 > &, const tetPoints &, const label, const label) |
template<class TetOp > | |
static void | decomposePrism (const FixedList< point, 6 > &points, TetOp &op) |
template<class AboveTetOp , class BelowTetOp > | |
static void | tetSliceWithPlane (const plane &pl, const tetPoints &tet, AboveTetOp &aboveOp, BelowTetOp &belowOp) |
Private Attributes | |
PointRef | a_ |
PointRef | b_ |
PointRef | c_ |
PointRef | d_ |
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 62 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 93 of file tetrahedron.H.
anonymous enum |
Enumerator | |
---|---|
nVertices | |
nEdges |
Definition at line 164 of file tetrahedron.H.
|
inline |
Construct from points.
Definition at line 35 of file tetrahedronI.H.
|
inline |
Construct from four points in the list of points.
|
inline |
Construct from Istream.
|
inlinestaticprivate |
Definition at line 580 of file tetrahedronI.H.
References tetrahedron::d().
Referenced by tetrahedron::tetSliceWithPlane().
Definition at line 596 of file tetrahedronI.H.
References points.
Referenced by tetrahedron::tetSliceWithPlane().
|
inlinestaticprivate |
Definition at line 611 of file tetrahedronI.H.
References Foam::abort(), tetrahedron::d(), tetrahedron::decomposePrism(), Foam::FatalError, FatalErrorInFunction, forAll, plane::normal(), p, tetrahedron::planeIntersection(), and plane::refPoint().
Referenced by tetrahedron::sliceWithPlane().
|
inline |
Return vertices.
Definition at line 73 of file tetrahedronI.H.
Referenced by volumeOptimizer::evaluateFunc(), volumeOptimizer::evaluateGradientsExact(), volumeOptimizer::evaluateStabilisationFactor(), particle< Type >::movingTetLambda(), Foam::help::pointInTetrahedron(), and simplexSmoother::simplexSmoother().
|
inline |
Definition at line 80 of file tetrahedronI.H.
Referenced by volumeOptimizer::evaluateFunc(), volumeOptimizer::evaluateGradientsExact(), volumeOptimizer::evaluateStabilisationFactor(), particle< Type >::movingTetLambda(), Foam::help::pointInTetrahedron(), and simplexSmoother::simplexSmoother().
|
inline |
Definition at line 87 of file tetrahedronI.H.
Referenced by volumeOptimizer::evaluateFunc(), volumeOptimizer::evaluateGradientsExact(), volumeOptimizer::evaluateStabilisationFactor(), particle< Type >::movingTetLambda(), Foam::help::pointInTetrahedron(), and simplexSmoother::simplexSmoother().
|
inline |
Definition at line 94 of file tetrahedronI.H.
Referenced by volumeOptimizer::evaluateFunc(), volumeOptimizer::evaluateGradientsExact(), volumeOptimizer::evaluateStabilisationFactor(), particle< Type >::movingTetLambda(), tetrahedron::planeIntersection(), Foam::help::pointInTetrahedron(), and tetrahedron::tetSliceWithPlane().
|
inline |
Return i-th face.
Definition at line 102 of file tetrahedronI.H.
References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.
|
inline |
Return face normal.
Definition at line 136 of file tetrahedronI.H.
References triangle::normal().
Referenced by particle< Type >::trackToFace().
|
inline |
Definition at line 143 of file tetrahedronI.H.
References triangle::normal().
Referenced by particle< Type >::trackToFace().
|
inline |
Definition at line 150 of file tetrahedronI.H.
References triangle::normal().
Referenced by particle< Type >::trackToFace().
|
inline |
Definition at line 157 of file tetrahedronI.H.
References triangle::normal().
Referenced by particle< Type >::trackToFace().
|
inline |
Return centre (centroid)
Definition at line 164 of file tetrahedronI.H.
Referenced by particle< Type >::findTris(), particle< Type >::hitWallFaces(), primitiveMesh::makeCellCentresAndVols(), tetOverlapVolume::sumMomentOp::operator()(), and particle< Type >::trackToFace().
|
inline |
Return volume.
Definition at line 171 of file tetrahedronI.H.
Referenced by Foam::help::angleBetweenFaces(), Foam::polyMeshGenChecks::checkCellPartTetrahedra(), meshSurfaceCheckEdgeTypes::classifyEdges(), triSurfaceClassifyEdges::classifyEdgesTypes(), Dual< Type >::Dual(), volumeOptimizer::evaluateFunc(), volumeOptimizer::evaluateGradientsExact(), volumeOptimizer::evaluateStabilisationFactor(), DSMCCloud< DSMCParcel< ParcelType > >::initialise(), Foam::help::isSharedEdgeConvex(), partTet::mag(), main(), primitiveMesh::makeCellCentresAndVols(), Moment< Type >::Moment(), tetOverlapVolume::sumMomentOp::operator()(), tetrahedron::sumVolOp::operator()(), tetOverlapVolume::tetTetOverlap(), and AveragingMethod< Foam::Vector >::write().
|
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.
Referenced by partTet::crcmCentre(), and main().
|
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.
Referenced by partTet::crcmRadius(), and main().
|
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(), and Foam::sqrt().
Referenced by polyMeshTetDecomposition::findBasePoint(), and polyMeshTetDecomposition::findSharedBasePoint().
Return a random point in the tetrahedron from a.
uniform distribution
Definition at line 245 of file tetrahedronI.H.
References rndGen(), s(), and cachedRandom::scalar01().
Referenced by DSMCCloud< DSMCParcel< ParcelType > >::initialise().
|
inline |
Return a random point in the tetrahedron from a.
uniform distribution
Definition at line 281 of file tetrahedronI.H.
References rndGen(), s(), and cachedRandom::sample01().
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(), List::setSize(), Vector< Cmpt >::x(), Vector< Cmpt >::y(), and Vector< Cmpt >::z().
Referenced by cellPointWeight::findTetrahedron(), interpolationCellPoint< Foam::Vector >::interpolate(), Dual< Type >::tetGeometry(), and Implicit< CloudType >::velocityCorrection().
|
inline |
Return nearest point to p on tetrahedron. Is p itself.
if inside.
Definition at line 362 of file tetrahedronI.H.
References PointHit::distance(), triangle::nearestPoint(), p, and PointHit::rawPoint().
Referenced by cellPointWeight::findTetrahedron().
|
inline |
Return true if point is inside tetrahedron.
Definition at line 456 of file tetrahedronI.H.
References Foam::mag(), n, and Vector< scalar >::zero.
Referenced by polyMeshTetDecomposition::findTet().
|
inline |
Decompose tet into tets above and below plane.
Definition at line 985 of file tetrahedronI.H.
References tetrahedron::a_, tetrahedron::b_, tetrahedron::c_, tetrahedron::d_, and tetrahedron::tetSliceWithPlane().
Referenced by tetOverlapVolume::tetTetOverlap().
|
inline |
Decompose tet into tets inside and outside other tet.
Definition at line 34 of file tetrahedron.C.
References tetrahedron::a_, tetrahedron::b_, tetrahedron::c_, and tetrahedron::d_.
Referenced by main().
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 130 of file tetrahedron.C.
References Foam::magSqr(), PointHit::rawPoint(), PointHit::setDistance(), PointHit::setHit(), PointHit::setMiss(), PointHit::setPoint(), and Foam::sqrt().
void gradNiSquared | ( | scalarField & | buffer | ) | const |
Fill buffer with shape function products.
Definition at line 335 of file tetrahedron.C.
References Foam::mag(), and Foam::magSqr().
void gradNiDotGradNj | ( | scalarField & | buffer | ) | const |
Definition at line 356 of file tetrahedron.C.
References Foam::mag().
void gradNiGradNi | ( | tensorField & | buffer | ) | const |
Definition at line 386 of file tetrahedron.C.
References Foam::mag(), and Foam::sqr().
void gradNiGradNj | ( | tensorField & | buffer | ) | const |
Definition at line 404 of file tetrahedron.C.
References Foam::mag().
|
friend |
|
friend |
|
private |
Definition at line 133 of file tetrahedron.H.
Referenced by tetrahedron::sliceWithPlane(), and tetrahedron::tetOverlap().
|
private |
Definition at line 133 of file tetrahedron.H.
Referenced by tetrahedron::sliceWithPlane(), and tetrahedron::tetOverlap().
|
private |
Definition at line 133 of file tetrahedron.H.
Referenced by tetrahedron::sliceWithPlane(), and tetrahedron::tetOverlap().
|
private |
Definition at line 133 of file tetrahedron.H.
Referenced by tetrahedron::sliceWithPlane(), and tetrahedron::tetOverlap().
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.