Data Structures | Public Types | Public Member Functions | Static Private Member Functions | Private Attributes | Friends
triangle Class Reference

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

Inheritance diagram for triangle:
Inheritance graph
[legend]

Data Structures

class  dummyOp
 Dummy. More...
 
class  storeOp
 Store resulting tris. More...
 
class  sumAreaOp
 Sum resulting areas. More...
 

Public Types

enum  proxType { NONE, POINT, EDGE }
 Return types for classify. More...
 
typedef FixedList< triPoints, 27 > triIntersectionList
 Storage type for triangles originating from intersecting triangle. More...
 

Public Member Functions

 triangle (const Point &a, const Point &b, const Point &c)
 Construct from three points. More...
 
 triangle (const FixedList< Point, 3 > &)
 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 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...
 
scalar mag () const
 Return scalar magnitude. More...
 
vector normal () const
 Return vector normal. 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=vector::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 randomPoint (cachedRandom &rndGen) const
 Return a random point on the triangle from a uniform. More...
 
scalar barycentric (const point &pt, List< scalar > &bary) const
 Calculate the barycentric coordinates of the given. 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...
 
template<class AboveOp , class BelowOp >
void sliceWithPlane (const plane &pl, 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...
 

Static Private Member Functions

static point planeIntersection (const FixedList< scalar, 3 > &d, const triPoints &t, const label negI, const label posI)
 Helper: calculate intersection point. More...
 
template<class AboveOp , class BelowOp >
static void triSliceWithPlane (const plane &pl, const triPoints &tri, AboveOp &aboveOp, BelowOp &belowOp)
 Helper: slice triangle with plane. More...
 

Private Attributes

PointRef a_
 
PointRef b_
 
PointRef c_
 

Friends

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

Detailed Description

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

Source files

Definition at line 59 of file triangle.H.

Member Typedef Documentation

◆ triIntersectionList

Storage type for triangles originating from intersecting triangle.

with another triangle

Definition at line 91 of file triangle.H.

Member Enumeration Documentation

◆ proxType

enum proxType

Return types for classify.

Enumerator
NONE 
POINT 
EDGE 

Definition at line 94 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 35 of file triangleI.H.

◆ triangle() [2/4]

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

Construct from three points.

◆ triangle() [3/4]

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

Construct from three points in the list of points.

The indices could be from triFace etc.

◆ triangle() [4/4]

triangle ( Istream )
inline

Construct from Istream.

Member Function Documentation

◆ planeIntersection()

Foam::point planeIntersection ( const FixedList< scalar, 3 > &  d,
const triPoints t,
const label  negI,
const label  posI 
)
inlinestaticprivate

Helper: calculate intersection point.

Definition at line 869 of file triangleI.H.

◆ triSliceWithPlane()

void triSliceWithPlane ( const plane pl,
const triPoints tri,
AboveOp &  aboveOp,
BelowOp &  belowOp 
)
inlinestaticprivate

Helper: slice triangle with plane.

Definition at line 35 of file triangle.C.

◆ a()

const Point & a ( ) const
inline

◆ b()

const Point & b ( ) const
inline

◆ c()

const Point & c ( ) const
inline

◆ centre()

Point centre ( ) const
inline

◆ mag()

Foam::scalar mag ( ) const
inline

Return scalar magnitude.

Definition at line 109 of file triangleI.H.

Referenced by pointLinear< Type >::correction(), FreeStream< CloudType >::inflow(), momentOfInertia::massPropertiesShell(), and triangle::sumAreaOp::operator()().

Here is the caller graph for this function:

◆ normal()

Foam::vector normal ( ) const
inline

◆ circumCentre()

Point circumCentre ( ) const
inline

Return circum-centre.

Definition at line 123 of file triangleI.H.

Referenced by main().

Here is the caller graph for this function:

◆ circumRadius()

Foam::scalar circumRadius ( ) const
inline

Return circum-radius.

Definition at line 150 of file triangleI.H.

Referenced by main().

Here is the caller graph for this function:

◆ quality()

Foam::scalar quality ( ) const
inline

Return quality: Ratio of triangle and circum-circle.

area, scaled so that an equilateral triangle has a quality of 1

Definition at line 174 of file triangleI.H.

◆ sweptVol()

Foam::scalar sweptVol ( const triangle t) const
inline

Return swept-volume.

Definition at line 190 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 209 of file triangleI.H.

◆ randomPoint() [1/2]

Point randomPoint ( Random rndGen) const
inline

Return a random point on the triangle from a uniform.

distribution

Definition at line 245 of file triangleI.H.

Referenced by FreeStream< CloudType >::inflow(), and patchInjectionBase::setPositionAndCell().

Here is the caller graph for this function:

◆ randomPoint() [2/2]

Point randomPoint ( cachedRandom rndGen) const
inline

Return a random point on the triangle from a uniform.

distribution

Definition at line 262 of file triangleI.H.

◆ barycentric()

Foam::scalar barycentric ( const point pt,
List< scalar > &  bary 
) const
inline

Calculate the barycentric coordinates of the given.

point, in the same order as a, b, c. Returns the determinant of the solution.

Definition at line 281 of file triangleI.H.

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

Here is the caller 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.

◆ 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.

Referenced by mappedPatchBase::facePoint(), and particle< Type >::hitWallFaces().

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.

Referenced by triSurfaceTools::calcInterpolationWeights(), and face::nearestPointClassify().

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.

Referenced by cellPointWeight::findTriangle(), and tetrahedron::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.

Referenced by ParticleCollector< CloudType >::collectParcelPolygon().

Here is the caller graph for this function:

◆ nearestPoint() [2/2]

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)

◆ sliceWithPlane()

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

Decompose triangle into triangles above and below plane.

Definition at line 113 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 126 of file triangle.C.

Friends And Related Function Documentation

◆ operator>>

Istream& operator>> ( Istream ,
triangle  
)
friend

◆ operator

Ostream& operator ( Ostream ,
const triangle  
)
friend

Field Documentation

◆ a_

PointRef a_
private

Definition at line 142 of file triangle.H.

Referenced by triangle< Foam::Vector, Foam::Vector >::sweptVol().

◆ b_

PointRef b_
private

Definition at line 142 of file triangle.H.

Referenced by triangle< Foam::Vector, Foam::Vector >::sweptVol().

◆ c_

PointRef c_
private

Definition at line 142 of file triangle.H.

Referenced by triangle< Foam::Vector, Foam::Vector >::sweptVol().


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