triangle< Point, PointRef > Class Template Reference

A triangle primitive used to calculate face areas and swept volumes. More...

Public Types

enum  proxType { NONE , POINT , EDGE }
 Return types for classify. More...
 

Public Member Functions

 triangle (const Point &a, const Point &b, const Point &c)
 Construct from three points. More...
 
 triangle (const UList< Point > &, const FixedList< label, 3 > &indices)
 Construct from three points in the list of points. More...
 
 triangle (Istream &)
 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 area () const
 Return vector area. More...
 
scalar mag () const
 Return scalar magnitude. More...
 
vector normal () const
 Return unit normal. More...
 
Tuple2< Point, scalar > circumCircle () const
 Return the circum centre and 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=vector::zero, scalar density=1.0) const
 Return the inertia tensor, with optional reference. More...
 
Point randomPoint (randomGenerator &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::algorithm::fullRay, const intersection::direction dir=intersection::direction::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...
 

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 areas and swept volumes.

Source files

Definition at line 79 of file triangle.H.

Member Enumeration Documentation

◆ proxType

enum proxType

Return types for classify.

Enumerator
NONE 
POINT 
EDGE 

Definition at line 89 of file triangle.H.

Constructor & Destructor Documentation

◆ triangle() [1/3]

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

Construct from three points.

Definition at line 33 of file triangleI.H.

◆ triangle() [2/3]

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 47 of file triangleI.H.

◆ triangle() [3/3]

triangle ( Istream is)
inline

Construct from Istream.

Definition at line 61 of file triangleI.H.

Member Function Documentation

◆ a()

const Point & a
inline

Return first vertex.

Definition at line 70 of file triangleI.H.

Referenced by triSurfaceTools::calcInterpolationWeights().

Here is the caller graph for this function:

◆ b()

const Point & b
inline

Return second vertex.

Definition at line 76 of file triangleI.H.

Referenced by triSurfaceTools::calcInterpolationWeights().

Here is the caller graph for this function:

◆ c()

const Point & c
inline

Return third vertex.

Definition at line 82 of file triangleI.H.

Referenced by triSurfaceTools::calcInterpolationWeights().

Here is the caller graph for this function:

◆ centre()

Point centre
inline

Return centre (centroid)

Definition at line 89 of file triangleI.H.

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

Here is the caller graph for this function:

◆ area()

◆ mag()

Foam::scalar mag
inline

Return scalar magnitude.

Definition at line 103 of file triangleI.H.

References Foam::mag().

Referenced by pointLinear< Type >::correction(), FreeStream< CloudType >::inflow(), momentOfInertia::massPropertiesShell(), and PatchFlowRateInjection< CloudType >::setPositionAndCell().

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

◆ normal()

Foam::vector normal
inline

Return unit normal.

Definition at line 110 of file triangleI.H.

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

Here is the call graph for this function:

◆ circumCircle()

Foam::Tuple2< Point, Foam::scalar > circumCircle
inline

Return the circum centre and radius.

Definition at line 120 of file triangleI.H.

References Foam::constant::physicoChemical::c1, Foam::constant::physicoChemical::c2, Foam::mag(), Foam::max(), Foam::sqrt(), and VectorSpace< Form, Cmpt, Ncmpts >::uniform().

Here is the call 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 151 of file triangleI.H.

References Foam::mag(), Foam::min(), 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 160 of file triangleI.H.

◆ inertia()

Foam::tensor inertia ( PointRef  refPt = vector::zero,
scalar  density = 1.0 
) const
inline

Return the inertia tensor, with optional reference.

point and density specification

Definition at line 179 of file triangleI.H.

References Foam::I, Foam::mag(), VectorSpace< Form, Cmpt, Mrows *Ncols >::one, Foam::fvm::S(), and Tensor< Cmpt >::T().

Here is the call graph for this function:

◆ randomPoint()

Point randomPoint ( randomGenerator rndGen) const
inline

Return a random point on the triangle from a uniform.

distribution

Definition at line 216 of file triangleI.H.

References Foam::barycentric2D01(), and rndGen().

Referenced by FreeStream< CloudType >::inflow().

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 226 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 236 of file triangleI.H.

Referenced by cellPointWeight::findTriangle().

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 248 of file triangleI.H.

References Foam::mag().

Here is the call graph for this function:

◆ ray()

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 287 of file triangleI.H.

References intersection::contactSphere, intersection::fullRay, intersection::halfRay, PointHit< Point >::hit(), Foam::mag(), Foam::min(), n, p, intersection::planarTol(), PointHit< Point >::rawPoint(), PointHit< Point >::setDistance(), PointHit< Point >::setHit(), PointHit< Point >::setMiss(), PointHit< Point >::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 algorithm::visible, algorithm::fullRay or algorithm::halfRay. tol increases the virtual size of the triangle by a relative factor.

Definition at line 408 of file triangleI.H.

References Foam::det(), intersection::fullRay, intersection::halfRay, intersection::visible, and Foam::Zero.

Referenced by mappedInternalPatchBase::samplePoints().

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 495 of file triangleI.H.

References Foam::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 644 of file triangleI.H.

References p.

Referenced by 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 658 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 670 of file triangleI.H.

References PointHit< Point >::distance(), intersection::fullRay, PointHit< Point >::hit(), PointHit< Point >::hitPoint(), Foam::ln(), Foam::mag(), Foam::minDist(), PointHit< Point >::missPoint(), PointHit< Point >::setDistance(), PointHit< Point >::setHit(), PointHit< Point >::setMiss(), and PointHit< Point >::setPoint().

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: