A triangle primitive used to calculate face normals and swept volumes. More...
Classes | |
class | dummyOp |
Dummy. More... | |
class | storeOp |
Store resulting tris. More... | |
class | sumAreaOp |
Sum resulting areas. More... | |
Public Types | |
enum | proxType { NONE = 0 , POINT , EDGE } |
The proximity classifications. More... | |
typedef Point | point_type |
The point type. More... | |
typedef FixedList< triPoints, 27 > | triIntersectionList |
Public Member Functions | |
triangle (const Point &a, const Point &b, const Point &c) | |
Construct from three points. More... | |
triangle (const FixedList< Point, 3 > &tri) | |
Construct from three points. More... | |
triangle (const UList< Point > &points, const FixedList< label, 3 > &indices) | |
Construct from three points in the list of points. More... | |
triangle (Istream &is) | |
Construct from Istream. More... | |
const Point & | a () const |
Return first vertex. More... | |
const Point & | b () const |
Return second vertex. More... | |
const Point & | c () const |
Return third vertex. More... | |
Point | centre () const |
Return centre (centroid) More... | |
vector | areaNormal () const |
The area normal - with magnitude equal to area of triangle. More... | |
vector | unitNormal () const |
Return unit normal. More... | |
FOAM_DEPRECATED_FOR (2018-12, "areaNormal() or unitNormal()") vector normal() const | |
Legacy name for areaNormal(). More... | |
scalar | mag () const |
Return scalar magnitude. More... | |
Point | circumCentre () const |
Return circum-centre. More... | |
scalar | circumRadius () const |
Return circum-radius. More... | |
scalar | quality () const |
Return quality: Ratio of triangle and circum-circle. More... | |
scalar | sweptVol (const triangle &t) const |
Return swept-volume. More... | |
tensor | inertia (PointRef refPt=Zero, scalar density=1.0) const |
Return the inertia tensor, with optional reference. More... | |
Point | randomPoint (Random &rndGen) const |
Return a random point on the triangle from a uniform. More... | |
Point | barycentricToPoint (const barycentric2D &bary) const |
Calculate the point from the given barycentric coordinates. More... | |
barycentric2D | pointToBarycentric (const point &pt) const |
Calculate the barycentric coordinates from the given point. More... | |
scalar | pointToBarycentric (const point &pt, barycentric2D &bary) const |
Calculate the barycentric coordinates from the given point. More... | |
pointHit | ray (const point &p, const vector &q, const intersection::algorithm=intersection::FULL_RAY, const intersection::direction dir=intersection::VECTOR) const |
Return point intersection with a ray. More... | |
pointHit | intersection (const point &p, const vector &q, const intersection::algorithm alg, const scalar tol=0.0) const |
Fast intersection with a ray. More... | |
pointHit | nearestPointClassify (const point &p, label &nearType, label &nearLabel) const |
Find the nearest point to p on the triangle and classify it: More... | |
pointHit | nearestPoint (const point &p) const |
Return nearest point to p on triangle. More... | |
bool | classify (const point &p, label &nearType, label &nearLabel) const |
Classify nearest point to p in triangle plane. More... | |
pointHit | nearestPoint (const linePointRef &edge, pointHit &edgePoint) const |
Return nearest point to line on triangle. Returns hit if. More... | |
int | sign (const point &p, const scalar tol=SMALL) const |
The sign for which side of the face plane the point is on. More... | |
template<class AboveOp , class BelowOp > | |
void | sliceWithPlane (const plane &pln, AboveOp &aboveOp, BelowOp &belowOp) const |
Decompose triangle into triangles above and below plane. More... | |
template<class InsideOp , class OutsideOp > | |
void | triangleOverlap (const vector &n, const triangle< Point, PointRef > &tri, InsideOp &insideOp, OutsideOp &outsideOp) const |
Decompose triangle into triangles inside and outside. More... | |
Friends | |
Istream & | operator>> (Istream &, triangle &) |
Ostream & | operator (Ostream &, const triangle &) |
A triangle primitive used to calculate face normals and swept volumes.
Definition at line 79 of file triangle.H.
typedef Point point_type |
The point type.
Definition at line 86 of file triangle.H.
typedef FixedList<triPoints, 27> triIntersectionList |
Storage type for triangles originating from intersecting triangle with another triangle
Definition at line 90 of file triangle.H.
enum proxType |
The proximity classifications.
Enumerator | |
---|---|
NONE | Unknown proximity. |
POINT | Close to point. |
EDGE | Close to edge. |
Definition at line 93 of file triangle.H.
Construct from three points.
Definition at line 37 of file triangleI.H.
Construct from three points.
Definition at line 51 of file triangleI.H.
Construct from three points in the list of points.
The indices could be from triFace etc.
Definition at line 63 of file triangleI.H.
Construct from Istream.
Definition at line 77 of file triangleI.H.
|
inline |
Return first vertex.
Definition at line 86 of file triangleI.H.
Referenced by triSurfaceTools::calcInterpolationWeights(), offsetSurface::operator()(), and triangle< Point, PointRef >::triangleOverlap().
|
inline |
Return second vertex.
Definition at line 92 of file triangleI.H.
Referenced by triSurfaceTools::calcInterpolationWeights(), offsetSurface::operator()(), and triangle< Point, PointRef >::triangleOverlap().
|
inline |
Return third vertex.
Definition at line 98 of file triangleI.H.
Referenced by triSurfaceTools::calcInterpolationWeights(), offsetSurface::operator()(), and triangle< Point, PointRef >::triangleOverlap().
|
inline |
Return centre (centroid)
Definition at line 105 of file triangleI.H.
Referenced by momentOfInertia::massPropertiesShell(), polyMesh::pointInCell(), and wallBoundedStreamLine::pushIn().
|
inline |
The area normal - with magnitude equal to area of triangle.
Definition at line 112 of file triangleI.H.
Referenced by Foam::calcEdgeNormalFromFace(), polyMeshGeometry::checkFaceTwist(), primitiveMeshGeometry::checkFaceTwist(), polyMeshGeometry::checkTriangleTwist(), polyMesh::pointInCell(), tetrahedron< Point, PointRef >::Sa(), tetrahedron< Point, PointRef >::Sb(), tetrahedron< Point, PointRef >::Sc(), and tetrahedron< Point, PointRef >::Sd().
|
inline |
Return unit normal.
Definition at line 119 of file triangleI.H.
References Foam::mag(), n, s(), and Foam::Zero.
|
inline |
Legacy name for areaNormal().
Definition at line 215 of file triangle.H.
|
inline |
Return scalar magnitude.
Definition at line 128 of file triangleI.H.
Referenced by pointLinear< Type >::correction(), FreeStream< CloudType >::inflow(), momentOfInertia::massPropertiesShell(), and triangle< Point, PointRef >::sumAreaOp::operator()().
|
inline |
Return circum-centre.
Definition at line 135 of file triangleI.H.
References Foam::mag().
Referenced by tetrahedron< Point, PointRef >::containmentSphere(), and Foam::edgeMeshTools::featureProximity().
|
inline |
Return circum-radius.
Definition at line 162 of file triangleI.H.
References Foam::mag(), Foam::max(), Foam::min(), and Foam::sqrt().
Referenced by Foam::edgeMeshTools::featureProximity().
|
inline |
Return quality: Ratio of triangle and circum-circle.
area, scaled so that an equilateral triangle has a quality of 1
Definition at line 183 of file triangleI.H.
References Foam::mag(), Foam::sqr(), and Foam::sqrt().
Return swept-volume.
Definition at line 198 of file triangleI.H.
|
inline |
Return the inertia tensor, with optional reference.
point and density specification
Definition at line 217 of file triangleI.H.
References Foam::I, Foam::mag(), and Tensor< Cmpt >::T().
Return a random point on the triangle from a uniform.
distribution
Definition at line 254 of file triangleI.H.
References Foam::barycentric2D01(), and rndGen.
Referenced by FreeStream< CloudType >::inflow(), and patchInjectionBase::setPositionAndCell().
|
inline |
Calculate the point from the given barycentric coordinates.
Definition at line 261 of file triangleI.H.
|
inline |
Calculate the barycentric coordinates from the given point.
Definition at line 271 of file triangleI.H.
Referenced by cellPointWeight::findTriangle(), and offsetSurface::operator()().
|
inline |
Calculate the barycentric coordinates from the given point.
Returns the determinant.
Definition at line 283 of file triangleI.H.
References Foam::mag().
|
inline |
Return point intersection with a ray.
For a hit, the distance is signed. Positive number represents the point in front of triangle. In case of miss pointHit is set to nearest point on triangle and its distance to the distance between the original point and the plane intersection point
Definition at line 322 of file triangleI.H.
References intersection::CONTACT_SPHERE, intersection::FULL_RAY, intersection::HALF_RAY, PointHit< PointType >::hit(), Foam::mag(), Foam::min(), n, p, intersection::planarTol(), PointHit< PointType >::rawPoint(), PointHit< PointType >::setDistance(), PointHit< PointType >::setHit(), PointHit< PointType >::setMiss(), PointHit< PointType >::setPoint(), intersection::VECTOR, and intersection::VISIBLE.
|
inline |
Fast intersection with a ray.
For a hit, the pointHit.distance() is the line parameter t : intersection=p+t*q. Only defined for VISIBLE, FULL_RAY or HALF_RAY. tol increases the virtual size of the triangle by a relative factor.
Definition at line 438 of file triangleI.H.
References Foam::det(), intersection::FULL_RAY, intersection::HALF_RAY, intersection::VISIBLE, and Foam::Zero.
Referenced by mappedPatchBase::facePoint(), and triangleFuncs::intersectBb().
Foam::pointHit nearestPointClassify | ( | const point & | p, |
label & | nearType, | ||
label & | nearLabel | ||
) | const |
Find the nearest point to p on the triangle and classify it:
+ near point (nearType=POINT, nearLabel=0, 1, 2) + near edge (nearType=EDGE, nearLabel=0, 1, 2) Note: edges are counted from starting vertex so e.g. edge 2 is from f[2] to f[0]
Definition at line 521 of file triangleI.H.
References cp, Foam::mag(), and p.
Referenced by triSurfaceTools::calcInterpolationWeights(), and face::nearestPointClassify().
|
inline |
Return nearest point to p on triangle.
Definition at line 670 of file triangleI.H.
References p.
Referenced by wallBoundedStreamLine::findNearestTet(), cellPointWeight::findTriangle(), and tetrahedron< Point, PointRef >::nearestPoint().
Classify nearest point to p in triangle plane.
w.r.t. triangle edges and points. Returns inside (true)/outside (false).
Definition at line 684 of file triangleI.H.
References p.
|
inline |
Return nearest point to line on triangle. Returns hit if.
point is inside triangle. Sets edgePoint to point on edge (hit if nearest is inside line)
Definition at line 696 of file triangleI.H.
References PointHit< PointType >::distance(), intersection::FULL_RAY, PointHit< PointType >::hit(), PointHit< PointType >::hitPoint(), Foam::ln(), Foam::mag(), PointHit< PointType >::missPoint(), PointHit< PointType >::setDistance(), PointHit< PointType >::setHit(), PointHit< PointType >::setMiss(), and PointHit< PointType >::setPoint().
|
inline |
The sign for which side of the face plane the point is on.
Uses the supplied tolerance for rounding around zero.
Definition at line 822 of file triangleI.H.
References p.
|
inline |
Decompose triangle into triangles above and below plane.
Definition at line 114 of file triangle.C.
|
inline |
Decompose triangle into triangles inside and outside.
(with respect to user provided normal) other triangle.
Definition at line 127 of file triangle.C.
References triangle< Point, PointRef >::a(), triangle< Point, PointRef >::b(), triangle< Point, PointRef >::c(), Foam::mag(), n, triangle< Point, PointRef >::storeOp::nTris_, s(), and triangle< Point, PointRef >::storeOp::tris_.