AXOM
Axom provides a robust, flexible software infrastructure for the development of multi-physics applications and computational tools.
axom::primal::Triangle< T, NDIMS > Class Template Reference

Represents a triangular geometric shape defined by three points. More...

#include </home/docs/checkouts/readthedocs.org/user_builds/axom/checkouts/latest/src/axom/primal/geometry/Triangle.hpp>

Inheritance diagram for axom::primal::Triangle< T, NDIMS >:

Public Types

using PointType = Point< T, NDIMS >
 
using VectorType = Vector< T, NDIMS >
 
using SphereType = Sphere< T, NDIMS >
 

Public Member Functions

AXOM_HOST_DEVICE Triangle ()
 Default constructor. Creates a degenerate triangle. More...
 
AXOM_HOST_DEVICE Triangle (const PointType &A, const PointType &B, const PointType &C)
 Custom Constructor. Creates a triangle from the 3 points A,B,C. More...
 
AXOM_HOST_DEVICE PointTypeoperator[] (int idx)
 Index operator to get the i^th vertex. More...
 
AXOM_HOST_DEVICE const PointTypeoperator[] (int idx) const
 Index operator to get the i^th vertex. More...
 
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if< TDIM==3, VectorType >::type normal () const
 Returns the normal of the triangle (not normalized) More...
 
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if< TDIM==3, double >::type area () const
 Returns the area of the triangle (3D specialization) More...
 
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if< TDIM==2, double >::type area () const
 Returns the area of the triangle (2D specialization) More...
 
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if< TDIM==2, double >::type signedArea () const
 Returns the signed area of a 2D triangle. More...
 
AXOM_HOST_DEVICE double volume () const
 Returns the volume of the triangle (synonym for area()) More...
 
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if< TDIM==2, double >::type signedVolume () const
 Returns the signed volume of a 2D triangle. More...
 
template<int TDIM = NDIMS>
std::enable_if< TDIM==2, SphereType >::type circumsphere () const
 Returns the circumsphere (circumscribing circle) of the triangle. More...
 
Point< double, 3 > physToBarycentric (const PointType &p, bool skipNormalization=false) const
 Returns the barycentric coordinates of a point within a triangle. More...
 
PointType baryToPhysical (const Point< double, 3 > &bary) const
 Returns the physical coordinates of a barycentric point. More...
 
AXOM_HOST_DEVICE bool degenerate (double eps=1.0e-12) const
 Returns whether the triangle is degenerate. More...
 
bool checkInTriangle (const PointType &p, double eps=1e-8) const
 Returns whether Point P is in the triangle for some 3d Triangle. More...
 
AXOM_HOST_DEVICE double angle (int idx) const
 Computes the request angle corresponding to the given vertex ID. More...
 
std::ostream & print (std::ostream &os) const
 Simple formatted print of a triangle instance. More...
 

Static Public Attributes

static constexpr int DIM = NDIMS
 
static constexpr int NUM_TRI_VERTS = 3
 

Detailed Description

template<typename T, int NDIMS>
class axom::primal::Triangle< T, NDIMS >

Represents a triangular geometric shape defined by three points.

Template Parameters
Tthe coordinate type, e.g., double, float, etc.
NDIMSthe number of dimensions

Member Typedef Documentation

◆ PointType

template<typename T , int NDIMS>
using axom::primal::Triangle< T, NDIMS >::PointType = Point<T, NDIMS>

◆ VectorType

template<typename T , int NDIMS>
using axom::primal::Triangle< T, NDIMS >::VectorType = Vector<T, NDIMS>

◆ SphereType

template<typename T , int NDIMS>
using axom::primal::Triangle< T, NDIMS >::SphereType = Sphere<T, NDIMS>

Constructor & Destructor Documentation

◆ Triangle() [1/2]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE axom::primal::Triangle< T, NDIMS >::Triangle ( )
inline

Default constructor. Creates a degenerate triangle.

◆ Triangle() [2/2]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE axom::primal::Triangle< T, NDIMS >::Triangle ( const PointType A,
const PointType B,
const PointType C 
)
inline

Custom Constructor. Creates a triangle from the 3 points A,B,C.

Parameters
[in]Apoint instance corresponding to vertex A of the triangle.
[in]Bpoint instance corresponding to vertex B of the triangle.
[in]Cpoint instance corresponding to vertex C of the triangle.

Member Function Documentation

◆ operator[]() [1/2]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE PointType& axom::primal::Triangle< T, NDIMS >::operator[] ( int  idx)
inline

Index operator to get the i^th vertex.

Parameters
idxThe index of the desired vertex
Precondition
idx is 0, 1 or 2

References axom::primal::Triangle< T, NDIMS >::NUM_TRI_VERTS, and SLIC_ASSERT.

◆ operator[]() [2/2]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE const PointType& axom::primal::Triangle< T, NDIMS >::operator[] ( int  idx) const
inline

Index operator to get the i^th vertex.

Parameters
idxThe index of the desired vertex
Precondition
idx is 0, 1 or 2

References axom::primal::Triangle< T, NDIMS >::NUM_TRI_VERTS, and SLIC_ASSERT.

◆ normal()

template<typename T , int NDIMS>
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if<TDIM == 3, VectorType>::type axom::primal::Triangle< T, NDIMS >::normal ( ) const
inline

Returns the normal of the triangle (not normalized)

Precondition
This function is only valid when NDIMS = 3
Returns
n triangle normal when NDIMS=3, zero vector otherwise

References axom::primal::Vector< T, NDIMS >::cross_product().

◆ area() [1/2]

template<typename T , int NDIMS>
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if<TDIM == 3, double>::type axom::primal::Triangle< T, NDIMS >::area ( ) const
inline

Returns the area of the triangle (3D specialization)

References axom::primal::Triangle< T, NDIMS >::normal().

◆ area() [2/2]

template<typename T , int NDIMS>
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if<TDIM == 2, double>::type axom::primal::Triangle< T, NDIMS >::area ( ) const
inline

Returns the area of the triangle (2D specialization)

References axom::utilities::abs(), and axom::primal::Triangle< T, NDIMS >::signedArea().

◆ signedArea()

template<typename T , int NDIMS>
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if<TDIM == 2, double>::type axom::primal::Triangle< T, NDIMS >::signedArea ( ) const
inline

Returns the signed area of a 2D triangle.

The area is positive when the vertices are oriented counter-clockwise.

Note
This function is only available for triangles in 2D since signed areas don't make sense for 3D triangles

References axom::numerics::determinant().

◆ volume()

template<typename T , int NDIMS>
AXOM_HOST_DEVICE double axom::primal::Triangle< T, NDIMS >::volume ( ) const
inline

Returns the volume of the triangle (synonym for area())

See also
area()

References axom::primal::Triangle< T, NDIMS >::area().

◆ signedVolume()

template<typename T , int NDIMS>
template<int TDIM = NDIMS>
AXOM_HOST_DEVICE std::enable_if<TDIM == 2, double>::type axom::primal::Triangle< T, NDIMS >::signedVolume ( ) const
inline

Returns the signed volume of a 2D triangle.

References axom::primal::Triangle< T, NDIMS >::signedArea().

◆ circumsphere()

template<typename T , int NDIMS>
template<int TDIM = NDIMS>
std::enable_if<TDIM == 2, SphereType>::type axom::primal::Triangle< T, NDIMS >::circumsphere ( ) const
inline

Returns the circumsphere (circumscribing circle) of the triangle.

Derived from formula on https://mathworld.wolfram.com/Circumcircle.html but uses observation that radius can be derived from distance to a vertex

Note
This function is only available for triangles in 2D

References axom::numerics::determinant().

◆ physToBarycentric()

template<typename T , int NDIMS>
Point<double, 3> axom::primal::Triangle< T, NDIMS >::physToBarycentric ( const PointType p,
bool  skipNormalization = false 
) const
inline

Returns the barycentric coordinates of a point within a triangle.

Parameters
[in]pThe point at which we want to compute barycentric coordinates
[in]skipNormalizationDetermines if the result should be normalized by the area of the triangle. One might want to skip this if they were only interested in the relative weights rather than their actual amounts
Precondition
The point lies in this triangle's plane.
Postcondition
The barycentric coordinates sum to 1 when skipNormalization is false Otherwise, the sum of coordinates will be proportional to area of the triangle

Algorithm adapted from Real Time Collision Detection by Christer Ericson.

References axom::primal::abs(), axom::primal::Triangle< T, NDIMS >::area(), axom::primal::Vector< T, NDIMS >::cross_product(), axom::primal::Triangle< T, NDIMS >::DIM, axom::utilities::isNearlyEqual(), and SLIC_CHECK.

◆ baryToPhysical()

template<typename T , int NDIMS>
PointType axom::primal::Triangle< T, NDIMS >::baryToPhysical ( const Point< double, 3 > &  bary) const
inline

Returns the physical coordinates of a barycentric point.

Parameters
[in]baryBarycentric coordinates relative to this triangle
Returns
Physical point represented by bary

References axom::primal::Point< T, NDIMS >::array(), axom::utilities::isNearlyEqual(), axom::primal::Triangle< T, NDIMS >::NUM_TRI_VERTS, and SLIC_CHECK_MSG.

◆ degenerate()

template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool axom::primal::Triangle< T, NDIMS >::degenerate ( double  eps = 1.0e-12) const
inline

Returns whether the triangle is degenerate.

Returns
true iff the triangle is degenerate (0 area)
See also
primal::Point

References axom::primal::Triangle< T, NDIMS >::area(), and axom::utilities::isNearlyEqual().

◆ checkInTriangle()

template<typename T , int NDIMS>
bool axom::primal::Triangle< T, NDIMS >::checkInTriangle ( const PointType p,
double  eps = 1e-8 
) const
inline

Returns whether Point P is in the triangle for some 3d Triangle.

Returns
true iff P is in the triangle
See also
primal::Point

References axom::utilities::isNearlyEqual(), and axom::primal::Triangle< T, NDIMS >::physToBarycentric().

◆ angle()

template<typename T , int NDIMS>
AXOM_HOST_DEVICE double axom::primal::Triangle< T, NDIMS >::angle ( int  idx) const
inline

Computes the request angle corresponding to the given vertex ID.

Parameters
[in]idxthe index of the corresponding vertex
Returns
alpha the incidence angle in the range [0, pi].
Precondition
idx >= 0 && idx < NUM_TRI_VERTS

References axom::utilities::clampVal(), axom::numerics::dot_product(), and SLIC_ASSERT.

◆ print()

template<typename T , int NDIMS>
std::ostream& axom::primal::Triangle< T, NDIMS >::print ( std::ostream &  os) const
inline

Simple formatted print of a triangle instance.

Parameters
osThe output stream to write to
Returns
A reference to the modified ostream

Member Data Documentation

◆ DIM

template<typename T , int NDIMS>
constexpr int axom::primal::Triangle< T, NDIMS >::DIM = NDIMS
staticconstexpr

◆ NUM_TRI_VERTS

template<typename T , int NDIMS>
constexpr int axom::primal::Triangle< T, NDIMS >::NUM_TRI_VERTS = 3
staticconstexpr

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