triangle< Point, PointRef > Class Template Reference

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 Pointa () const
 Return first vertex. More...
 
const Pointb () const
 Return second vertex. More...
 
const Pointc () 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

Istreamoperator>> (Istream &, triangle &)
 
Ostreamoperator (Ostream &, const triangle &)
 

Detailed Description

template<class Point, class PointRef>
class Foam::triangle< Point, PointRef >

A triangle primitive used to calculate face normals and swept volumes.

Source files

Definition at line 79 of file triangle.H.

Member Typedef Documentation

◆ point_type

typedef Point point_type

The point type.

Definition at line 86 of file triangle.H.

◆ triIntersectionList

Storage type for triangles originating from intersecting triangle with another triangle

Definition at line 90 of file triangle.H.

Member Enumeration Documentation

◆ proxType

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.

Constructor & Destructor Documentation

◆ triangle() [1/4]

triangle ( const Point a,
const Point b,
const Point c 
)
inline

Construct from three points.

Definition at line 37 of file triangleI.H.

◆ triangle() [2/4]

triangle ( const FixedList< Point, 3 > &  tri)
inline

Construct from three points.

Definition at line 51 of file triangleI.H.

◆ triangle() [3/4]

triangle ( const UList< Point > &  points,
const FixedList< label, 3 > &  indices 
)
inline

Construct from three points in the list of points.

The indices could be from triFace etc.

Definition at line 63 of file triangleI.H.

◆ triangle() [4/4]

triangle ( Istream is)
inlineexplicit

Construct from Istream.

Definition at line 77 of file triangleI.H.

Member Function Documentation

◆ a()

const Point & a
inline

Return first vertex.

Definition at line 86 of file triangleI.H.

Referenced by triSurfaceTools::calcInterpolationWeights(), offsetSurface::operator()(), and triangle< Point, PointRef >::triangleOverlap().

Here is the caller graph for this function:

◆ b()

const Point & b
inline

Return second vertex.

Definition at line 92 of file triangleI.H.

Referenced by triSurfaceTools::calcInterpolationWeights(), offsetSurface::operator()(), and triangle< Point, PointRef >::triangleOverlap().

Here is the caller graph for this function:

◆ c()

const Point & c
inline

Return third vertex.

Definition at line 98 of file triangleI.H.

Referenced by triSurfaceTools::calcInterpolationWeights(), offsetSurface::operator()(), and triangle< Point, PointRef >::triangleOverlap().

Here is the caller graph for this function:

◆ centre()

Point centre
inline

Return centre (centroid)

Definition at line 105 of file triangleI.H.

Referenced by momentOfInertia::massPropertiesShell(), polyMesh::pointInCell(), and wallBoundedStreamLine::pushIn().

Here is the caller graph for this function:

◆ areaNormal()

Foam::vector areaNormal
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().

Here is the caller graph for this function:

◆ unitNormal()

Foam::vector unitNormal
inline

Return unit normal.

Definition at line 119 of file triangleI.H.

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

Here is the call graph for this function:

◆ FOAM_DEPRECATED_FOR()

FOAM_DEPRECATED_FOR ( 2018-  12,
"areaNormal() or unitNormal()"   
) const
inline

Legacy name for areaNormal().

Deprecated:
(2018-06) Deprecated for new use

Definition at line 215 of file triangle.H.

◆ mag()

Foam::scalar mag
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()().

Here is the caller graph for this function:

◆ circumCentre()

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

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

◆ circumRadius()

Foam::scalar circumRadius
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().

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

◆ quality()

Foam::scalar quality
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().

Here is the call graph for this function:

◆ sweptVol()

Foam::scalar sweptVol ( const triangle< Point, PointRef > &  t) const
inline

Return swept-volume.

Definition at line 198 of file triangleI.H.

◆ inertia()

Foam::tensor inertia ( PointRef  refPt = Zero,
scalar  density = 1.0 
) const
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().

Here is the call graph for this function:

◆ randomPoint()

Point randomPoint ( Random rndGen) const
inline

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

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

◆ barycentricToPoint()

Point barycentricToPoint ( const barycentric2D bary) const
inline

Calculate the point from the given barycentric coordinates.

Definition at line 261 of file triangleI.H.

◆ pointToBarycentric() [1/2]

Foam::barycentric2D pointToBarycentric ( const point pt) const
inline

Calculate the barycentric coordinates from the given point.

Definition at line 271 of file triangleI.H.

Referenced by cellPointWeight::findTriangle(), and offsetSurface::operator()().

Here is the caller graph for this function:

◆ pointToBarycentric() [2/2]

Foam::scalar pointToBarycentric ( const point pt,
barycentric2D bary 
) const
inline

Calculate the barycentric coordinates from the given point.

Returns the determinant.

Definition at line 283 of file triangleI.H.

References Foam::mag().

Here is the call graph for this function:

◆ ray()

Foam::pointHit ray ( const point p,
const vector q,
const intersection::algorithm  alg = intersection::FULL_RAY,
const intersection::direction  dir = intersection::VECTOR 
) const
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.

Here is the call graph for this function:

◆ intersection()

Foam::pointHit intersection ( const point p,
const vector q,
const intersection::algorithm  alg,
const scalar  tol = 0.0 
) const
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().

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

◆ nearestPointClassify()

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

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

◆ nearestPoint() [1/2]

Foam::pointHit nearestPoint ( const point p) const
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().

Here is the caller graph for this function:

◆ classify()

bool classify ( const point p,
label &  nearType,
label &  nearLabel 
) const
inline

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.

◆ nearestPoint() [2/2]

Foam::pointHit nearestPoint ( const linePointRef edge,
pointHit edgePoint 
) const
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().

Here is the call graph for this function:

◆ sign()

int sign ( const point p,
const scalar  tol = SMALL 
) const
inline

The sign for which side of the face plane the point is on.

Uses the supplied tolerance for rounding around zero.

Returns
  • 0: on plane
  • +1: front-side
  • -1: back-side

Definition at line 822 of file triangleI.H.

References p.

◆ sliceWithPlane()

void sliceWithPlane ( const plane pln,
AboveOp &  aboveOp,
BelowOp &  belowOp 
) const
inline

Decompose triangle into triangles above and below plane.

Definition at line 114 of file triangle.C.

◆ triangleOverlap()

void triangleOverlap ( const vector n,
const triangle< Point, PointRef > &  tri,
InsideOp &  insideOp,
OutsideOp &  outsideOp 
) const
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_.

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator>>

Istream & operator>> ( Istream ,
triangle< Point, PointRef > &   
)
friend

◆ operator

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

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