AXOM
Axom provides a robust, flexible software infrastructure for the development of multi-physics applications and computational tools.
axom::primal Namespace Reference

Classes

class  BezierCurve
 Represents a Bezier curve defined by an array of control points. More...
 
class  BezierPatch
 Represents a 3D Bezier patch defined by a 2D array of control points. More...
 
class  BoundingBox
 
class  CurvedPolygon
 Represents a polygon with curved edges defined by BezierCurves. More...
 
class  Hexahedron
 Represents an hexahedral geometric shape defined by eight points. More...
 
class  KnotVector
 Represents the knot vector for a B-Spline/NURBS Curve/Surface. More...
 
class  NumericArray
 A simple statically sized array of data with component-wise operators. More...
 
struct  NonChar
 Type trait to avoid outputting chars when a value is expected This avoids unintentionally outputting system beeps. More...
 
struct  NonChar< char >
 
struct  NonChar< unsigned char >
 
class  NURBSCurve
 Represents a NURBS curve defined by an array of control points, weights and knots. More...
 
class  NURBSPatch
 Represents a 3D NURBS patch defined by a 2D array of control points. More...
 
class  Octahedron
 Represents an octahedral geometric shape defined by six points. More...
 
class  OrientedBoundingBox
 OrientedBoundingBox represents an oriented bounding box defined by its centroid, axes, and extents. More...
 
class  Plane
 Defines an oriented plane within a 2-D or 3-D Euclidean space and provides associated operators, such as, projection, signed distance, orientation, etc. More...
 
class  Point
 The point class represents a point, \( p \in \mathcal{R}^d \) . It provides access methods to set and query the point coordinates. More...
 
class  Polygon
 Represents a polygon defined by an array of points. More...
 
class  Polyhedron
 Represents a convex polyhedron defined by an array of points and optionally their neighbors (a point is a neighbor of another point if there is an edge between them) More...
 
class  NeighborCollection
 Represents a collection of neighbor relations between vertices. More...
 
class  Quadrilateral
 Represents a quadrilateral geometric shape defined by four points. More...
 
class  Ray
 
class  Segment
 
class  Sphere
 Defines an oriented Sphere in 2-D (i.e., a circle) or 3-D given by its center, \( \mathcal{X} \) and radius \( \mathcal{R} \). The Sphere object provides associated operations on a sphere, such as, signed distance and orientation. More...
 
class  Tetrahedron
 Represents a tetrahedral geometric shape defined by four points. More...
 
class  Triangle
 Represents a triangular geometric shape defined by three points. More...
 
class  Vector
 Represents a vector, \( v \in \mathcal{R}^d \). It provides access methods for setting and querying the vector components as well as vector math operators, e.g., adding, subtracting, dot_product and cross_product. More...
 

Typedefs

Pre-defined point types
using Point2D = Point< double, 2 >
 
using Point3D = Point< double, 3 >
 
Pre-defined Vector types
using Vector2D = Vector< double, 2 >
 
using Vector3D = Vector< double, 3 >
 

Enumerations

enum  OrientationResult { ON_BOUNDARY , ON_POSITIVE_SIDE , ON_NEGATIVE_SIDE }
 Enumerates possible return values for orientation tests. More...
 
enum class  PolygonArray { Dynamic , Static }
 

Functions

template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const BezierCurve< T, NDIMS > &bCurve)
 Overloaded output operator for Bezier Curves. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const BezierPatch< T, NDIMS > &bPatch)
 Overloaded output operator for Bezier Patches. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const CurvedPolygon< T, NDIMS > &poly)
 Overloaded output operator for polygons. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Hexahedron< T, NDIMS > &hex)
 Free functions implementing Hexahedron's operators. More...
 
template<typename T >
std::ostream & operator<< (std::ostream &os, const KnotVector< T > &kVec)
 Overloaded output operator for knot vectors. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const NURBSCurve< T, NDIMS > &bCurve)
 Overloaded output operator for NURBS Curves. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const NURBSPatch< T, NDIMS > &nPatch)
 Overloaded output operator for NURBS Patches. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Octahedron< T, NDIMS > &oct)
 Free functions implementing Octahedron's operators. More...
 
template<typename T , int NDIMS, axom::primal::PolygonArray ARRAY_TYPE, int MAX_VERTS>
std::ostream & operator<< (std::ostream &os, const Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS > &poly)
 Overloaded output operator for polygons. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Polyhedron< T, NDIMS > &poly)
 Overloaded output operator for polyhedrons. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Quadrilateral< T, NDIMS > &quad)
 Overloaded output operator for quadrilaterals. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Ray< T, NDIMS > &ray)
 Overloaded output operator for rays. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Tetrahedron< T, NDIMS > &tet)
 Free functions implementing Tetrahedron's operators. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Triangle< T, NDIMS > &tri)
 Overloaded output operator for triangles. More...
 
template<typename T >
Polygon< T, 3 > clip (const Triangle< T, 3 > &tri, const BoundingBox< T, 3 > &bbox)
 Clips a 3D triangle against an axis-aligned bounding box in 3D. More...
 
template<typename T , axom::primal::PolygonArray ARRAY_TYPE, int MAX_VERTS>
AXOM_SUPPRESS_HD_WARN AXOM_HOST_DEVICE Polygon< T, 2, ARRAY_TYPE, MAX_VERTS > clip (const Polygon< T, 2, ARRAY_TYPE, MAX_VERTS > &subjectPolygon, const Polygon< T, 2, ARRAY_TYPE, MAX_VERTS > &clipPolygon, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 2D subject polygon against a clip polygon in 2D, returning their geometric intersection as a polygon. More...
 
template<typename T >
AXOM_HOST_DEVICE Polyhedron< T, 3 > clip (const Hexahedron< T, 3 > &hex, const Tetrahedron< T, 3 > &tet, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 3D hexahedron against a tetrahedron in 3D, returning the geometric intersection of the hexahedron and the tetrahedron as a polyhedron. More...
 
template<typename T >
AXOM_HOST_DEVICE Polyhedron< T, 3 > clip (const Tetrahedron< T, 3 > &tet, const Hexahedron< T, 3 > &hex, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 3D hexahedron against a tetrahedron in 3D, returning the geometric intersection of the hexahedron and the tetrahedron as a polyhedron. More...
 
template<typename T >
AXOM_HOST_DEVICE Polyhedron< T, 3 > clip (const Octahedron< T, 3 > &oct, const Tetrahedron< T, 3 > &tet, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 3D octahedron against a tetrahedron in 3D, returning the geometric intersection of the octahedron and the tetrahedron as a polyhedron. More...
 
template<typename T >
AXOM_HOST_DEVICE Polyhedron< T, 3 > clip (const Tetrahedron< T, 3 > &tet, const Octahedron< T, 3 > &oct, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 3D octahedron against a tetrahedron in 3D, returning the geometric intersection of the octahedron and the tetrahedron as a polyhedron. More...
 
template<typename T >
AXOM_HOST_DEVICE Polyhedron< T, 3 > clip (const Tetrahedron< T, 3 > &tet1, const Tetrahedron< T, 3 > &tet2, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 3D tetrahedron against another tetrahedron in 3D, returning the geometric intersection as a polyhedron. More...
 
template<typename T >
AXOM_HOST_DEVICE Polyhedron< T, 3 > clip (const Tetrahedron< T, 3 > &tet, const Plane< T, 3 > &plane, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 3D tetrahedron against the half-space defined by a plane and returns the resulting polyhedron. More...
 
template<typename T >
AXOM_HOST_DEVICE Polyhedron< T, 3 > clip (const Plane< T, 3 > &plane, const Tetrahedron< T, 3 > &tet, double eps=1.e-10, bool tryFixOrientation=false)
 Clips a 3D tetrahedron against the half-space defined by a plane and returns the resulting polyhedron. More...
 
template<typename T , int NDIMS>
OrientedBoundingBox< T, NDIMS > merge_boxes (const OrientedBoundingBox< T, NDIMS > &l, const OrientedBoundingBox< T, NDIMS > &r)
 Creates an oriented bounding box which contains the passed in OBBs. More...
 
template<typename T , int NDIMS>
BoundingBox< T, NDIMS > merge_boxes (const BoundingBox< T, NDIMS > &l, const BoundingBox< T, NDIMS > &r)
 Constructor. Creates a bounding box which contains the passed in bounding boxes. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox< T, NDIMS > compute_bounding_box (const Triangle< T, NDIMS > &tri)
 Creates a bounding box around a Triangle. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox< T, NDIMS > compute_bounding_box (const Quadrilateral< T, NDIMS > &quad)
 Creates a bounding box around a Quadrilateral. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox< T, NDIMS > compute_bounding_box (const Octahedron< T, NDIMS > &oct)
 Creates a bounding box around an Octahedron. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox< T, NDIMS > compute_bounding_box (const Hexahedron< T, NDIMS > &hex)
 Creates a bounding box around a Hexahedron. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox< T, NDIMS > compute_bounding_box (const Polyhedron< T, NDIMS > &poly)
 Creates a bounding box around a Polyhedron. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox< T, NDIMS > compute_bounding_box (const Tetrahedron< T, NDIMS > &tet)
 Creates a bounding box around a Tetrahedron. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox< T, NDIMS > compute_bounding_box (const Polygon< T, NDIMS > &poly)
 Creates a bounding box around a Polygon. More...
 
template<typename T >
sector_area (const primal::BezierCurve< T, 2 > &curve)
 Calculates the sector area of a planar, nonrational Bezier Curve. More...
 
template<typename T >
primal::Point< T, 2 > sector_centroid (const primal::BezierCurve< T, 2 > &curve)
 Calculates the sector centroid of a planar, nonrational Bezier Curve. More...
 
template<typename T >
area (const primal::CurvedPolygon< T, 2 > &poly, double tol=1e-8)
 Returns the area enclosed by the CurvedPolygon. More...
 
template<typename T >
primal::Point< T, 2 > centroid (const primal::CurvedPolygon< T, 2 > &poly, double tol=1e-8)
 Returns the centroid of the CurvedPolygon. More...
 
template<typename Lambda , typename T , int NDIMS>
double evaluate_scalar_line_integral (const primal::CurvedPolygon< T, NDIMS > cpoly, Lambda &&scalar_integrand, int npts)
 Evaluate a line integral along the boundary of a CurvedPolygon object on a scalar field. More...
 
template<typename Lambda , typename T , int NDIMS>
double evaluate_scalar_line_integral (const primal::BezierCurve< T, NDIMS > &c, Lambda &&scalar_integrand, int npts)
 Evaluate a line integral on a single Bezier curve on a scalar field. More...
 
template<typename Lambda , typename T , int NDIMS>
double evaluate_vector_line_integral (const primal::CurvedPolygon< T, NDIMS > cpoly, Lambda &&vector_integrand, int npts)
 Evaluate a line integral along the boundary of a CurvedPolygon object on a vector field. More...
 
template<typename Lambda , typename T , int NDIMS>
double evaluate_vector_line_integral (const primal::BezierCurve< T, NDIMS > &c, Lambda &&vector_integrand, int npts)
 Evaluate a line integral on a single Bezier curve on a vector field. More...
 
template<class Lambda , typename T >
double evaluate_area_integral (const primal::CurvedPolygon< T, 2 > cpoly, Lambda &&integrand, int npts_Q, int npts_P=0)
 Evaluate an integral on the interior of a CurvedPolygon object. More...
 
template<typename T >
bool in_curved_polygon (const Point< T, 2 > &query, const CurvedPolygon< T, 2 > &cpoly, bool useNonzeroRule=true, double edge_tol=1e-8, double EPS=1e-8)
 Robustly determine if query point is interior to a curved polygon. More...
 
template<typename T >
bool in_polygon (const Point< T, 2 > &query, const Polygon< T, 2 > &poly, bool includeBoundary=false, bool useNonzeroRule=true, double edge_tol=1e-8)
 Determines containment for a point in a polygon. More...
 
template<typename T >
bool in_polyhedron (const Point< T, 3 > &query, const Polyhedron< T, 3 > &poly, bool includeBoundary=false, bool useNonzeroRule=true, double edge_tol=1e-8, double EPS=1e-8)
 Determines containment for a point in a polygon. More...
 
template<typename T >
bool in_sphere (const Point< T, 2 > &q, const Point< T, 2 > &p0, const Point< T, 2 > &p1, const Point< T, 2 > &p2, double EPS=1e-8)
 Tests whether a query point lies inside a 2D triangle's circumcircle. More...
 
template<typename T >
bool in_sphere (const Point< T, 2 > &q, const Triangle< T, 2 > &tri, double EPS=1e-8)
 Tests whether a query point lies inside a 2D triangle's circumcircle. More...
 
template<typename T >
bool in_sphere (const Point< T, 3 > &q, const Point< T, 3 > &p0, const Point< T, 3 > &p1, const Point< T, 3 > &p2, const Point< T, 3 > &p3, double EPS=1e-8)
 Tests whether a query point lies inside a 3D tetrahedron's circumsphere. More...
 
template<typename T >
bool in_sphere (const Point< T, 3 > &q, const Tetrahedron< T, 3 > &tet, double EPS=1e-8)
 Tests whether a query point lies inside a 3D tetrahedron's circumsphere. More...
 
template<typename T >
AXOM_HOST_DEVICEintersection_volume (const Hexahedron< T, 3 > &hex, const Tetrahedron< T, 3 > &tet, double eps=1.e-10, bool tryFixOrientation=false)
 Finds the absolute (unsigned) intersection volume between a hexahedron and a tetrahedron. More...
 
template<typename T >
AXOM_HOST_DEVICEintersection_volume (const Tetrahedron< T, 3 > &tet, const Hexahedron< T, 3 > &hex, double eps=1.e-10, bool tryFixOrientation=false)
 Finds the absolute (unsigned) intersection volume between a tetrahedron and a hexahedron. More...
 
template<typename T >
AXOM_HOST_DEVICEintersection_volume (const Octahedron< T, 3 > &oct, const Tetrahedron< T, 3 > &tet, double eps=1.e-10, bool tryFixOrientation=false)
 Finds the absolute (unsigned) intersection volume between a octahedron and a tetrahedron. More...
 
template<typename T >
AXOM_HOST_DEVICEintersection_volume (const Tetrahedron< T, 3 > &tet, const Octahedron< T, 3 > &oct, double eps=1.e-10, bool tryFixOrientation=false)
 Finds the absolute (unsigned) intersection volume between a tetrahedron and a octahedron. More...
 
template<typename T >
AXOM_HOST_DEVICEintersection_volume (const Tetrahedron< T, 3 > &tet1, const Tetrahedron< T, 3 > &tet2, double eps=1.e-10, bool tryFixOrientation=false)
 Finds the absolute (unsigned) intersection volume between a tetrahedron and another tetrahedron. More...
 
template<typename T >
bool is_convex (const Polygon< T, 2 > &poly, double EPS=1e-8)
 Determines if a polygon defined by ordered vertices is convex. More...
 
template<typename T >
int orientation (const Point< T, 3 > &p, const Triangle< T, 3 > &tri, double EPS=1e-9)
 Computes the orientation of a point p with respect to an oriented triangle tri. More...
 
template<typename T >
int orientation (const Point< T, 2 > &p, const Segment< T, 2 > &seg, double EPS=1e-9)
 Computes the orientation of a point p with respect to an oriented segment. More...
 
template<typename Tp , int NDIMS = 3>
void split (const Octahedron< Tp, NDIMS > &oct, axom::Array< Tetrahedron< Tp, NDIMS >> &out)
 Splits an Octahedron into eight Tetrahedrons. More...
 
double squared_distance (const double *A, const double *B, int N)
 Computes the squared distance from point A to point B, represented by arrays of length N. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE double squared_distance (const Point< T, NDIMS > &A, const Point< T, NDIMS > &B)
 Computes the squared distance from point A to point B. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE double squared_distance (const Point< T, NDIMS > &P, const BoundingBox< T, NDIMS > &B)
 Computes the minimum squared distance from a query point, P, to a given axis-aligned bounding box B. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE double squared_distance (const BoundingBox< T, NDIMS > &A, const BoundingBox< T, NDIMS > &B)
 Computes the minimum squared distance between 2 axis-aligned boxes. More...
 
template<typename T , int NDIMS>
double squared_distance (const Point< T, NDIMS > &P, const Segment< T, NDIMS > &S)
 Computes the minimum squared distance from a query point, P, to the given segment, S. More...
 
template<typename T , int NDIMS>
double squared_distance (const Point< T, NDIMS > &P, const Triangle< T, NDIMS > &tri)
 Computes the minimum squared distance from a query point, P, to the closest point on the given triangle. More...
 
Forward Declared Overloaded Operators
template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool operator== (const BoundingBox< T, NDIMS > &lhs, const BoundingBox< T, NDIMS > &rhs)
 Equality comparison operator for bounding boxes. Two bounding boxes are equal when they have the same bounds. More...
 
template<typename T , int NDIMS>
bool operator!= (const BoundingBox< T, NDIMS > &lhs, const BoundingBox< T, NDIMS > &rhs)
 Inequality comparison operator for bounding boxes. Two bounding boxes are unequal when they have different bounds. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const BoundingBox< T, NDIMS > &bb)
 Overloaded output operator for bounding boxes. More...
 
template<typename T , int SIZE>
AXOM_HOST_DEVICE bool operator== (const NumericArray< T, SIZE > &lhs, const NumericArray< T, SIZE > &rhs)
 Checks if two numeric arrays are component-wise equal. More...
 
template<typename T , int SIZE>
bool operator!= (const NumericArray< T, SIZE > &lhs, const NumericArray< T, SIZE > &rhs)
 Checks if two numeric arrays are not component-wise equal. More...
 
template<typename T , int SIZE>
AXOM_HOST_DEVICE NumericArray< T, SIZE > operator+ (const NumericArray< T, SIZE > &lhs, const NumericArray< T, SIZE > &rhs)
 Performs component-wise addition of two numeric arrays. More...
 
template<typename T , int SIZE>
AXOM_HOST_DEVICE NumericArray< T, SIZE > operator- (const NumericArray< T, SIZE > &lhs, const NumericArray< T, SIZE > &rhs)
 Performs component-wise subtraction of two numeric arrays. More...
 
template<typename T , int SIZE>
NumericArray< T, SIZE > operator- (const NumericArray< T, SIZE > &arr)
 Unary negation of a numeric array instance. More...
 
template<typename T , int SIZE>
NumericArray< T, SIZE > operator* (const NumericArray< T, SIZE > &arr, double scalar)
 Scalar multiplication a numeric array; Scalar on rhs. More...
 
template<typename T , int SIZE>
NumericArray< T, SIZE > operator* (double scalar, const NumericArray< T, SIZE > &arr)
 Scalar multiplication a numeric array; Scalar on lhs. More...
 
template<typename T , int SIZE>
AXOM_HOST_DEVICE NumericArray< T, SIZE > operator* (const NumericArray< T, SIZE > &lhs, const NumericArray< T, SIZE > &rhs)
 Component-wise multiplication of NumericArrays. More...
 
template<typename T , int SIZE>
NumericArray< T, SIZE > operator/ (const NumericArray< T, SIZE > &lhs, const NumericArray< T, SIZE > &rhs)
 Component-wise division of NumericArrays. More...
 
template<typename T , int SIZE>
NumericArray< T, SIZE > operator/ (const NumericArray< T, SIZE > &arr, double scalar)
 Scalar division of NumericArray; Scalar on rhs. More...
 
template<typename T , int SIZE>
NumericArray< T, SIZE > abs (const NumericArray< T, SIZE > &arr)
 Coordinate-wise absolute value on the NumericArray. More...
 
template<typename T , int SIZE>
std::ostream & operator<< (std::ostream &os, const NumericArray< T, SIZE > &arr)
 Overloaded output operator for numeric arrays. More...
 
template<typename T , int NDIMS>
bool operator== (const OrientedBoundingBox< T, NDIMS > &lhs, const OrientedBoundingBox< T, NDIMS > &rhs)
 Equality comparison operator for oriented bounding boxes. More...
 
template<typename T , int NDIMS>
bool operator!= (const OrientedBoundingBox< T, NDIMS > &lhs, const OrientedBoundingBox< T, NDIMS > &rhs)
 Inequality comparison operator for oriented bounding boxes. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const OrientedBoundingBox< T, NDIMS > &box)
 Overloaded output operator for oriented bounding boxes. More...
 
template<typename T , int NDIMS>
OrientedBoundingBox< T, NDIMS > compute_oriented_bounding_box (const Point< T, NDIMS > *pts, int n)
 Creates a bounding box which contains the collection of passed in points. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Plane< T, NDIMS > &p)
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool operator== (const Point< T, NDIMS > &lhs, const Point< T, NDIMS > &rhs)
 Equality comparison operator for points. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool operator!= (const Point< T, NDIMS > &lhs, const Point< T, NDIMS > &rhs)
 Inequality comparison operator for points. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Point< T, NDIMS > &pt)
 Overloaded output operator for points. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool operator== (const Segment< T, NDIMS > &lhs, const Segment< T, NDIMS > &rhs)
 Equality comparison operator for Segment. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool operator!= (const Segment< T, NDIMS > &lhs, const Segment< T, NDIMS > &rhs)
 Inequality comparison operator for Segment. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Segment< T, NDIMS > &seg)
 Overloaded output operator for Segment. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Sphere< T, NDIMS > &s)
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > operator+ (const Vector< T, NDIMS > &A, const Vector< T, NDIMS > &B)
 Adds vectors A, B and stores the result into a new vector C. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Point< T, NDIMS > operator+ (const Point< T, NDIMS > &P, const Vector< T, NDIMS > &V)
 Adds vector V to point P and stores the result into a new point. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Point< T, NDIMS > operator+ (const Vector< T, NDIMS > &V, const Point< T, NDIMS > &P)
 Adds vector V to point P and stores the result into a new point. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Point< T, NDIMS > operator- (const Point< T, NDIMS > &P, const Vector< T, NDIMS > &V)
 Subtracts vector V from point P and stores the result into a new point. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > operator- (const Vector< T, NDIMS > &A, const Vector< T, NDIMS > &B)
 Subtracts vectors A, B and stores the result into a new vector C. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > operator- (const Point< T, NDIMS > &h, const Point< T, NDIMS > &t)
 Subtracts Point t from Point h, yielding a vector. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > operator- (const Vector< T, NDIMS > &vec1)
 Unary negation of a vector instance. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > operator* (const Vector< T, NDIMS > &vec, const T scalar)
 Scalar multiplication of vector; Scalar on rhs. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > operator* (const T scalar, const Vector< T, NDIMS > &vec)
 Scalar multiplication of vector; Scalar on lhs. More...
 
template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > operator/ (const Vector< T, NDIMS > &vec, const T scalar)
 Scalar division of vector; Scalar on rhs. More...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Vector< T, NDIMS > &vec)
 Overloaded output operator for vectors. More...
 
Forward-declared helper functions
template<typename T >
AXOM_HOST_DEVICE Plane< T, 2 > make_plane (const Point< T, 2 > &x1, const Point< T, 2 > &x2)
 Constructs a Plane in 2D that goes through the specified points. More...
 
template<typename T >
AXOM_HOST_DEVICE Plane< T, 3 > make_plane (const Point< T, 3 > &x1, const Point< T, 3 > &x2, const Point< T, 3 > &x3)
 Constructs a Plane in 3D that goes through the specified points. More...
 
Triangle Intersection Routines
template<typename T >
AXOM_HOST_DEVICE bool intersect (const Triangle< T, 3 > &t1, const Triangle< T, 3 > &t2, bool includeBoundary=false, double EPS=1E-08)
 Tests if 3D Triangles t1 and t2 intersect. More...
 
template<typename T >
bool intersect (const Triangle< T, 2 > &t1, const Triangle< T, 2 > &t2, bool includeBoundary=false, double EPS=1E-08)
 Tests if 2D Triangles t1 and t2 intersect. More...
 
template<typename T >
bool intersect (const Triangle< T, 3 > &tri, const BoundingBox< T, 3 > &bb)
 Determines if a triangle and a bounding box intersect. More...
 
template<typename T >
bool intersect (const Triangle< T, 3 > &tri, const Ray< T, 3 > &ray)
 Determines if a 3D triangle intersects a 3D ray. More...
 
template<typename T >
bool intersect (const Triangle< T, 3 > &tri, const Ray< T, 3 > &ray, T &t)
 Determines if a 3D triangle intersects a 3D ray. More...
 
template<typename T >
bool intersect (const Triangle< T, 3 > &tri, const Ray< T, 3 > &ray, T &t, Point< double, 3 > &p)
 Determines if a 3D triangle intersects a 3D ray. More...
 
template<typename T >
bool intersect (const Triangle< T, 3 > &tri, const Segment< T, 3 > &seg)
 Determines if a 3D triangle intersects a 3D segment. More...
 
template<typename T >
bool intersect (const Triangle< T, 3 > &tri, const Segment< T, 3 > &seg, T &t)
 Determines if a 3D triangle intersects a 3D segment. More...
 
template<typename T >
bool intersect (const Triangle< T, 3 > &tri, const Segment< T, 3 > &seg, T &t, Point< double, 3 > &p)
 Determines if a 3D triangle intersects a 3D segment. More...
 
Ray Intersection Routines
template<typename T >
bool intersect (const Ray< T, 2 > &R, const Segment< T, 2 > &S, T &ray_param, T &seg_param, const T EPS=1e-8)
 Computes the intersection of the given ray, R, with the segment, S. More...
 
template<typename T >
bool intersect (const Ray< T, 2 > &R, const Segment< T, 2 > &S, T &ray_param)
 Computes the intersection of the given ray, R, with the segment, S. More...
 
template<typename T >
bool intersect (const Ray< T, 2 > &R, const Segment< T, 2 > &S, Point< T, 2 > &ip, const T EPS=1e-8)
 Computes the intersection of the given ray, R, with the segment, S. More...
 
template<typename T , int DIM>
AXOM_HOST_DEVICE bool intersect (const Ray< T, DIM > &R, const BoundingBox< T, DIM > &bb, Point< T, DIM > &ip)
 Computes the intersection of the given ray, R, with the Box, bb. More...
 
Segment-BoundingBox Intersection Routines
template<typename T , int DIM>
bool intersect (const Segment< T, DIM > &S, const BoundingBox< T, DIM > &bb, T &tmin, T &tmax, const double &EPS=1e-8)
 Computes the intersection of the given segment, S, with the Box, bb. If an intersection is found, output parameter ip contains an intersection point. More...
 
template<typename T , int DIM>
bool intersect (const Segment< T, DIM > &S, const BoundingBox< T, DIM > &bb, Point< T, DIM > &ip, const double &EPS=1e-8)
 This variant returns a point within the intersection as an OUT parameters. More...
 
template<typename T , int DIM>
bool intersect (const Segment< T, DIM > &S, const BoundingBox< T, DIM > &bb)
 
Axis-Aligned Bounding Box Intersection Routines
template<typename T , int DIM>
AXOM_HOST_DEVICE bool intersect (const BoundingBox< T, DIM > &bb1, const BoundingBox< T, DIM > &bb2)
 Determines if two axis aligned bounding boxes intersect. More...
 
Sphere Intersection Routines
template<typename T , int DIM>
bool intersect (const Sphere< T, DIM > &s1, const Sphere< T, DIM > &s2, double TOL=1.e-9)
 Determines if two spheres intersect. More...
 
Oriented Bounding Box Intersection Routines
template<typename T >
bool intersect (const OrientedBoundingBox< T, 1 > &b1, const OrientedBoundingBox< T, 1 > &b2)
 
template<typename T >
bool intersect (const OrientedBoundingBox< T, 2 > &b1, const OrientedBoundingBox< T, 2 > &b2)
 Determines if a 2D OBB intersects a 2D OBB. More...
 
template<typename T >
bool intersect (const OrientedBoundingBox< T, 3 > &b1, const OrientedBoundingBox< T, 3 > &b2, double EPS=1E-4)
 Determines if a 3D OBB intersects a 3D OBB. More...
 
Bezier Curve Intersection Routines
template<typename T >
bool intersect (const BezierCurve< T, 2 > &c1, const BezierCurve< T, 2 > &c2, std::vector< T > &sp, std::vector< T > &tp, double tol=1E-8)
 Tests if two Bezier Curves c1 and c2 intersect. More...
 
Plane Intersection Routines
template<typename T >
AXOM_HOST_DEVICE bool intersect (const Plane< T, 3 > &p, const BoundingBox< T, 3 > &bb, bool checkOverlaps=false, double EPS=1E-08)
 Determines if a 3D plane intersects a 3D bounding box. By default (checkOverlaps is false), checks if |s| <= r, where "s" is the distance of the bounding box center to the plane, and "r" is the projected radius of the bounding box along the line parallel to the plane normal and going through the box center. If checkOverlaps is true, checks if |s| < r, where the bounding box overlaps both half spaces of the plane. More...
 
template<typename T , int DIM>
AXOM_HOST_DEVICE bool intersect (const Plane< T, DIM > &plane, const Segment< T, DIM > &seg, T &t)
 Determines if a plane intersects a segment. More...
 
template<typename T >
AXOM_HOST_DEVICE bool intersect (const Plane< T, 3 > &p, const Tetrahedron< T, 3 > &tet, Polygon< T, 3 > &intersection)
 Determines if a 3D plane intersects a tetrahedron. More...
 
Winding number operations between 2D points and primitives
template<typename T >
double winding_number (const Point< T, 2 > &q, const Segment< T, 2 > &s, double edge_tol=1e-8)
 
template<typename T >
int winding_number (const Point< T, 2 > &q, const Triangle< T, 2 > &tri, bool includeBoundary=false, double edge_tol=1e-8)
 
template<typename T >
int winding_number (const Point< T, 2 > &R, const Polygon< T, 2 > &P, bool &isOnEdge, bool includeBoundary, double edge_tol)
 Computes the winding number for a 2D point wrt a 2D polygon. More...
 
template<typename T >
int winding_number (const Point< T, 2 > &R, const Polygon< T, 2 > &P, bool includeBoundary=false, double edge_tol=1e-8)
 Computes the winding number for a 2D point wrt a 2D polygon. More...
 
template<typename T >
double winding_number (const Point< T, 2 > &q, const BezierCurve< T, 2 > &c, double edge_tol=1e-8, double EPS=1e-8)
 Computes the GWN for a 2D point wrt a 2D Bezier curve. More...
 
template<typename T >
double winding_number (const Point< T, 2 > &q, const NURBSCurve< T, 2 > &n, double edge_tol=1e-8, double EPS=1e-8)
 Computes the GWN for a 2D point wrt a 2D NURBS curve. More...
 
template<typename T >
double winding_number (const Point< T, 2 > &q, const CurvedPolygon< T, 2 > &cpoly, double edge_tol=1e-8, double EPS=1e-8)
 Computes the GWN for a 2D point wrt to a 2D curved polygon. More...
 
template<typename T >
double winding_number (const Point< T, 2 > &q, const axom::Array< BezierCurve< T, 2 >> &carray, double edge_tol=1e-8, double EPS=1e-8)
 Computes the GWN for a 2D point wrt to a collection of 2D Bezier curves. More...
 
template<typename T >
double winding_number (const Point< T, 2 > &q, const axom::Array< NURBSCurve< T, 2 >> &narray, double edge_tol=1e-8, double EPS=1e-8)
 Computes the GWN for a 2D point wrt to a collection of 2D NURBS curves. More...
 
Winding number operations between 3D points and primitives
template<typename T >
double winding_number (const Point< T, 3 > &q, const Triangle< T, 3 > &tri, bool &isOnFace, const double edge_tol=1e-8, const double EPS=1e-8)
 Computes the GWN for a 3D point wrt a 3D triangle. More...
 
template<typename T >
double winding_number (const Point< T, 3 > &q, const Triangle< T, 3 > &tri, const double edge_tol=1e-8, const double EPS=1e-8)
 Computes the GWN for a 3D point wrt a 3D triangle. More...
 
template<typename T >
double winding_number (const Point< T, 3 > &q, const Polygon< T, 3 > &poly, bool &isOnFace, const double edge_tol=1e-8, const double EPS=1e-8)
 Computes the GWN for a 3D point wrt a 3D planar polygon. More...
 
template<typename T >
double winding_number (const Point< T, 3 > &q, const Polygon< T, 3 > &poly, const double edge_tol=1e-8, const double EPS=1e-8)
 Computes the GWN for a 3D point wrt a 3D planar polygon. More...
 
template<typename T >
int winding_number (const Point< T, 3 > &query, const Polyhedron< T, 3 > &poly, bool includeBoundary=false, double edge_tol=1e-8, double EPS=1e-8)
 Computes the winding number for a 3D point wrt a 3D convex polyhedron. More...
 

Typedef Documentation

◆ Point2D

using axom::primal::Point2D = typedef Point<double, 2>

◆ Point3D

using axom::primal::Point3D = typedef Point<double, 3>

◆ Vector2D

using axom::primal::Vector2D = typedef Vector<double, 2>

◆ Vector3D

using axom::primal::Vector3D = typedef Vector<double, 3>

Enumeration Type Documentation

◆ OrientationResult

Enumerates possible return values for orientation tests.

Enumerator
ON_BOUNDARY 

primitive is on the boundary of a primitive

ON_POSITIVE_SIDE 

primitive is on the positive side of a primitive

ON_NEGATIVE_SIDE 

primitive is on the negative side of a primitive

◆ PolygonArray

The polygon can have a dynamic or static array to store vertices

Enumerator
Dynamic 
Static 

Function Documentation

◆ operator<<() [1/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const BezierCurve< T, NDIMS > &  bCurve 
)

Overloaded output operator for Bezier Curves.

Free functions related to BezierCurve.

References axom::primal::BezierCurve< T, NDIMS >::print().

◆ operator<<() [2/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const BezierPatch< T, NDIMS > &  bPatch 
)

Overloaded output operator for Bezier Patches.

Free functions related to BezierPatch.

References axom::primal::BezierPatch< T, NDIMS >::print().

◆ operator==() [1/5]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool axom::primal::operator== ( const BoundingBox< T, NDIMS > &  lhs,
const BoundingBox< T, NDIMS > &  rhs 
)

Equality comparison operator for bounding boxes. Two bounding boxes are equal when they have the same bounds.

Free functions implementing comparison and arithmetic operators.

References axom::primal::BoundingBox< T, NDIMS >::getMax(), and axom::primal::BoundingBox< T, NDIMS >::getMin().

◆ operator!=() [1/5]

template<typename T , int NDIMS>
bool axom::primal::operator!= ( const BoundingBox< T, NDIMS > &  lhs,
const BoundingBox< T, NDIMS > &  rhs 
)

Inequality comparison operator for bounding boxes. Two bounding boxes are unequal when they have different bounds.

◆ operator<<() [3/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const BoundingBox< T, NDIMS > &  bb 
)

Overloaded output operator for bounding boxes.

References axom::primal::BoundingBox< T, NDIMS >::print().

◆ operator<<() [4/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const CurvedPolygon< T, NDIMS > &  poly 
)

Overloaded output operator for polygons.

Free functions implementing Polygon's operators.

References axom::primal::CurvedPolygon< T, NDIMS >::print().

◆ operator<<() [5/22]

template<typename T , int NDIMS>
std::ostream& axom::primal::operator<< ( std::ostream &  os,
const Hexahedron< T, NDIMS > &  hex 
)

Free functions implementing Hexahedron's operators.

References axom::primal::Hexahedron< T, NDIMS >::print().

◆ operator<<() [6/22]

template<typename T >
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const KnotVector< T > &  kVec 
)

Overloaded output operator for knot vectors.

Free functions related to KnotVector.

References axom::primal::KnotVector< T >::print().

◆ operator==() [2/5]

template<typename T , int SIZE>
AXOM_HOST_DEVICE bool axom::primal::operator== ( const NumericArray< T, SIZE > &  lhs,
const NumericArray< T, SIZE > &  rhs 
)

Checks if two numeric arrays are component-wise equal.

Free functions implementing comparison and arithmetic operators.

Parameters
[in]lhsnumeric array instance on the left-hand side.
[in]rhsnumeric array instance on the right-hand side.
Returns
status true if lhs==rhs, otherwise, false.

◆ operator!=() [2/5]

template<typename T , int SIZE>
bool axom::primal::operator!= ( const NumericArray< T, SIZE > &  lhs,
const NumericArray< T, SIZE > &  rhs 
)

Checks if two numeric arrays are not component-wise equal.

Parameters
[in]lhsnumeric array instance on the left-hand side.
[in]rhsnumeric array instance on the right-hand side.
Returns
status true if lhs!=rhs, otherwise, false.

◆ operator+() [1/4]

template<typename T , int SIZE>
AXOM_HOST_DEVICE NumericArray< T, SIZE > axom::primal::operator+ ( const NumericArray< T, SIZE > &  lhs,
const NumericArray< T, SIZE > &  rhs 
)
inline

Performs component-wise addition of two numeric arrays.

Parameters
[in]lhsnumeric array instance on the left-hand side.
[in]rhsnumeric array instance on the right-hand side.
Returns
C resulting numeric array from the component-wise addition.

◆ operator-() [1/6]

template<typename T , int SIZE>
AXOM_HOST_DEVICE NumericArray< T, SIZE > axom::primal::operator- ( const NumericArray< T, SIZE > &  lhs,
const NumericArray< T, SIZE > &  rhs 
)
inline

Performs component-wise subtraction of two numeric arrays.

Parameters
[in]lhsnumeric array instance on the left-hand side.
[in]rhsnumeric array instance on the right-hand side.
Returns
C resulting numeric array from component-wise subtraction.

◆ operator-() [2/6]

template<typename T , int SIZE>
NumericArray< T, SIZE > axom::primal::operator- ( const NumericArray< T, SIZE > &  arr)
inline

Unary negation of a numeric array instance.

Parameters
[in]arrnumeric array instance on the left-hand side.
Returns
C resulting numeric array from unary negation.

◆ operator*() [1/5]

template<typename T , int SIZE>
NumericArray< T, SIZE > axom::primal::operator* ( const NumericArray< T, SIZE > &  arr,
double  scalar 
)
inline

Scalar multiplication a numeric array; Scalar on rhs.

Parameters
[in]arrnumeric array instance.
[in]scalaruser-supplied scalar.
Returns
C resutling numeric array, \( \ni: C_i = scalar*arr_i, \forall i\)

◆ operator*() [2/5]

template<typename T , int SIZE>
NumericArray< T, SIZE > axom::primal::operator* ( double  scalar,
const NumericArray< T, SIZE > &  arr 
)
inline

Scalar multiplication a numeric array; Scalar on lhs.

Parameters
[in]scalaruser-supplied scalar.
[in]arrnumeric array instance.
Returns
C resulting numeric array, \( \ni: C_i = scalar*arr_i, \forall i\)

◆ operator*() [3/5]

template<typename T , int SIZE>
AXOM_HOST_DEVICE NumericArray< T, SIZE > axom::primal::operator* ( const NumericArray< T, SIZE > &  lhs,
const NumericArray< T, SIZE > &  rhs 
)
inline

Component-wise multiplication of NumericArrays.

Parameters
[in]lhsnumeric array instance on the left-hand side.
[in]rhsnumeric array instance on the right-hand side.
Returns
C resulting numeric array, \( \ni: C_i = lhs_i * rhs_i, \forall i\)

◆ operator/() [1/3]

template<typename T , int SIZE>
NumericArray< T, SIZE > axom::primal::operator/ ( const NumericArray< T, SIZE > &  lhs,
const NumericArray< T, SIZE > &  rhs 
)
inline

Component-wise division of NumericArrays.

Parameters
[in]lhsnumeric array instance on the left-hand side.
[in]rhsnumeric array instance on the right-hand side.
Returns
C resulting numeric array, \( \ni: C_i = lhs_i / rhs_i, \forall i\)
Precondition
\( rhs_i != 0.0, \forall i \)

◆ operator/() [2/3]

template<typename T , int SIZE>
NumericArray< T, SIZE > axom::primal::operator/ ( const NumericArray< T, SIZE > &  arr,
double  scalar 
)
inline

Scalar division of NumericArray; Scalar on rhs.

Parameters
[in]arrnumeric array instance
[in]scalaruser-supplied scalar
Returns
C resulting numeric array, \( \ni: C_i = arr_i/scalar, \forall i\)
Precondition
scalar != 0.0

◆ abs()

template<typename T , int SIZE>
NumericArray< T, SIZE > axom::primal::abs ( const NumericArray< T, SIZE > &  arr)
inline

Coordinate-wise absolute value on the NumericArray.

Parameters
[in]arrnumeric array instance
Precondition
std::abs is defined for template type T
Returns
A NumericArray whose coordinates are the absolute value of arr

References axom::utilities::abs().

◆ operator<<() [7/22]

template<typename T , int SIZE>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const NumericArray< T, SIZE > &  arr 
)

Overloaded output operator for numeric arrays.

Parameters
[in]osC++ output stream
[in]arrnumeric array instance.

References axom::primal::NumericArray< T, SIZE >::print().

◆ operator<<() [8/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const NURBSCurve< T, NDIMS > &  bCurve 
)

Overloaded output operator for NURBS Curves.

Free functions related to NURBSCurve.

References axom::primal::NURBSCurve< T, NDIMS >::print().

◆ operator<<() [9/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const NURBSPatch< T, NDIMS > &  nPatch 
)

Overloaded output operator for NURBS Patches.

Free functions related to NURBSPatch.

References axom::primal::NURBSPatch< T, NDIMS >::print().

◆ operator<<() [10/22]

template<typename T , int NDIMS>
std::ostream& axom::primal::operator<< ( std::ostream &  os,
const Octahedron< T, NDIMS > &  oct 
)

Free functions implementing Octahedron's operators.

References axom::primal::Octahedron< T, NDIMS >::print().

◆ operator==() [3/5]

template<typename T , int NDIMS>
bool axom::primal::operator== ( const OrientedBoundingBox< T, NDIMS > &  lhs,
const OrientedBoundingBox< T, NDIMS > &  rhs 
)

Equality comparison operator for oriented bounding boxes.

Free functions implementing comparison and arithmetic operators.

Two oriented bounding boxes are equal when they have the same bounds and axes.

References axom::primal::OrientedBoundingBox< T, NDIMS >::contains().

◆ operator!=() [3/5]

template<typename T , int NDIMS>
bool axom::primal::operator!= ( const OrientedBoundingBox< T, NDIMS > &  lhs,
const OrientedBoundingBox< T, NDIMS > &  rhs 
)

Inequality comparison operator for oriented bounding boxes.

Two bounding boxes are unequal when they have different bounds

◆ operator<<() [11/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const OrientedBoundingBox< T, NDIMS > &  box 
)

Overloaded output operator for oriented bounding boxes.

References axom::primal::OrientedBoundingBox< T, NDIMS >::print().

◆ compute_oriented_bounding_box()

template<typename T , int NDIMS>
OrientedBoundingBox< T, NDIMS > axom::primal::compute_oriented_bounding_box ( const Point< T, NDIMS > *  pts,
int  n 
)

Creates a bounding box which contains the collection of passed in points.

Call srand() to initialize the random number generator before using this function.

Parameters
[in]ptsarray of points
[in]nnumber of points
Note
if n <= 0, invokes default constructor

◆ operator<<() [12/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Plane< T, NDIMS > &  p 
)

◆ make_plane() [1/2]

template<typename T >
AXOM_HOST_DEVICE Plane< T, 2 > axom::primal::make_plane ( const Point< T, 2 > &  x1,
const Point< T, 2 > &  x2 
)

Constructs a Plane in 2D that goes through the specified points.

Parameters
[in]x1coordinates of the first point.
[in]x2coordinates of the second point.
Returns
plane the plane/line defined by the two points
Precondition
x1 should not be equal to x2

References axom::primal::Vector< T, NDIMS >::is_zero(), and SLIC_CHECK_MSG.

◆ make_plane() [2/2]

template<typename T >
AXOM_HOST_DEVICE Plane< T, 3 > axom::primal::make_plane ( const Point< T, 3 > &  x1,
const Point< T, 3 > &  x2,
const Point< T, 3 > &  x3 
)

Constructs a Plane in 3D that goes through the specified points.

Parameters
[in]x1coordinates of the first point.
[in]x2coordinates of the second point.
[in]x3coordinates of the third point.
Returns
plane the plane defined by the three points
Precondition
The user-supplied points, x1, x2, x3 should not be collinear.

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

◆ operator==() [4/5]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool axom::primal::operator== ( const Point< T, NDIMS > &  lhs,
const Point< T, NDIMS > &  rhs 
)

Equality comparison operator for points.

◆ operator!=() [4/5]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool axom::primal::operator!= ( const Point< T, NDIMS > &  lhs,
const Point< T, NDIMS > &  rhs 
)

Inequality comparison operator for points.

Inequality operator for points.

◆ operator<<() [13/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Point< T, NDIMS > &  pt 
)

Overloaded output operator for points.

Free functions implementing Point's operators.

References axom::primal::Point< T, NDIMS >::print().

◆ operator<<() [14/22]

template<typename T , int NDIMS, axom::primal::PolygonArray ARRAY_TYPE, int MAX_VERTS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS > &  poly 
)

Overloaded output operator for polygons.

Free functions implementing Polygon's operators.

References axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::print().

◆ operator<<() [15/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Polyhedron< T, NDIMS > &  poly 
)

Overloaded output operator for polyhedrons.

Free functions implementing Polyhedron's operators.

References axom::primal::Polyhedron< T, NDIMS >::print().

◆ operator<<() [16/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Quadrilateral< T, NDIMS > &  quad 
)

Overloaded output operator for quadrilaterals.

Free functions implementing Quadrilateral's operators.

References axom::primal::Quadrilateral< T, NDIMS >::print().

◆ operator<<() [17/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Ray< T, NDIMS > &  ray 
)

Overloaded output operator for rays.

Free functions implementing Ray's operators.

References axom::primal::Ray< T, NDIMS >::print().

◆ operator==() [5/5]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool axom::primal::operator== ( const Segment< T, NDIMS > &  lhs,
const Segment< T, NDIMS > &  rhs 
)

Equality comparison operator for Segment.

Equality comparison operator for segments.

◆ operator!=() [5/5]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE bool axom::primal::operator!= ( const Segment< T, NDIMS > &  lhs,
const Segment< T, NDIMS > &  rhs 
)

Inequality comparison operator for Segment.

Inequality operator for segments.

◆ operator<<() [18/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Segment< T, NDIMS > &  seg 
)

Overloaded output operator for Segment.

Free functions implementing Segments's operators.

References axom::primal::Segment< T, NDIMS >::print().

◆ operator<<() [19/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Sphere< T, NDIMS > &  s 
)

◆ operator<<() [20/22]

template<typename T , int NDIMS>
std::ostream& axom::primal::operator<< ( std::ostream &  os,
const Tetrahedron< T, NDIMS > &  tet 
)

Free functions implementing Tetrahedron's operators.

References axom::primal::Tetrahedron< T, NDIMS >::print().

◆ operator<<() [21/22]

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

Overloaded output operator for triangles.

Free functions implementing Triangle's operators.

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

◆ operator+() [2/4]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > axom::primal::operator+ ( const Vector< T, NDIMS > &  A,
const Vector< T, NDIMS > &  B 
)
inline

Adds vectors A, B and stores the result into a new vector C.

Parameters
[in]Avector on the left-hand side.
[in]Bvector on the right-hand side.
Returns
C resulting vector, \( C_i = A_i + B_i \forall i \)

◆ operator+() [3/4]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Point< T, NDIMS > axom::primal::operator+ ( const Point< T, NDIMS > &  P,
const Vector< T, NDIMS > &  V 
)

Adds vector V to point P and stores the result into a new point.

Parameters
[in]Ppoint on the left-hand side.
[in]Vvector on the right-hand side.
Returns
resulting point, \( p'_i = p_i + v_i \forall i \)

References axom::primal::Point< T, NDIMS >::array(), and axom::primal::Vector< T, NDIMS >::array().

◆ operator+() [4/4]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Point< T, NDIMS > axom::primal::operator+ ( const Vector< T, NDIMS > &  V,
const Point< T, NDIMS > &  P 
)

Adds vector V to point P and stores the result into a new point.

Parameters
[in]Vvector on the left-hand side.
[in]Ppoint on the right-hand side.
Returns
resulting point, \( p'_i = v_i + p_i \forall i \)

◆ operator-() [3/6]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Point< T, NDIMS > axom::primal::operator- ( const Point< T, NDIMS > &  P,
const Vector< T, NDIMS > &  V 
)

Subtracts vector V from point P and stores the result into a new point.

Parameters
[in]Ppoint on the left-hand side.
[in]Vvector on the right-hand side.
Returns
resulting point, \( p'_i = p_i + v_i \forall i \)

References axom::primal::Point< T, NDIMS >::array(), and axom::primal::Vector< T, NDIMS >::array().

◆ operator-() [4/6]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > axom::primal::operator- ( const Vector< T, NDIMS > &  A,
const Vector< T, NDIMS > &  B 
)
inline

Subtracts vectors A, B and stores the result into a new vector C.

Parameters
[in]Avector on the left-hand side.
[in]Bvector on the right-hand side.
Returns
C resulting vector, \( C_i = A_i - B_i \forall i \)

◆ operator-() [5/6]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > axom::primal::operator- ( const Point< T, NDIMS > &  h,
const Point< T, NDIMS > &  t 
)

Subtracts Point t from Point h, yielding a vector.

Parameters
[in]hthe head of the resulting vector
[in]tthe tail of the resulting vector
Returns
resulting vector, \( V_i = h_i - t_i \forall i \)

◆ operator-() [6/6]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > axom::primal::operator- ( const Vector< T, NDIMS > &  vec1)
inline

Unary negation of a vector instance.

Parameters
[in]vec1vector instance to negate.
Returns
C resulting vector from unary negation.

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

◆ operator*() [4/5]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > axom::primal::operator* ( const Vector< T, NDIMS > &  vec,
const T  scalar 
)
inline

Scalar multiplication of vector; Scalar on rhs.

Free functions involving vectors.

Parameters
[in]vecvector instance.
[in]scalaruser-supplied scalar.
Returns
C resulting vector, \( \ni: C_i = scalar*vec_i, \forall i\)

◆ operator*() [5/5]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > axom::primal::operator* ( const T  scalar,
const Vector< T, NDIMS > &  vec 
)
inline

Scalar multiplication of vector; Scalar on lhs.

Parameters
[in]scalaruser-supplied scalar.
[in]vecvector instance.
Returns
C resulting vector, \( \ni: C_i = scalar*vec_i, \forall i\)

◆ operator/() [3/3]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Vector< T, NDIMS > axom::primal::operator/ ( const Vector< T, NDIMS > &  vec,
const T  scalar 
)
inline

Scalar division of vector; Scalar on rhs.

Parameters
[in]vecvector instance
[in]scalaruser-supplied scalar.
Returns
C resulting vector, \( \ni: C_i = vec_i / scalar, \forall i\)
Precondition
scalar != 0.0

◆ operator<<() [22/22]

template<typename T , int NDIMS>
std::ostream & axom::primal::operator<< ( std::ostream &  os,
const Vector< T, NDIMS > &  vec 
)

Overloaded output operator for vectors.

Parameters
[in]osC++ output stream
[in]vecvector instance

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

◆ clip() [1/9]

template<typename T >
Polygon<T, 3> axom::primal::clip ( const Triangle< T, 3 > &  tri,
const BoundingBox< T, 3 > &  bbox 
)

Clips a 3D triangle against an axis-aligned bounding box in 3D.

Parameters
[in]triThe triangle to clip
[in]bboxThe bounding box to clip against
Returns
A planar polygon of the triangle clipped against the bounding box. If the triangle is completely outside the bounding box, the returned polygon is empty (i.e. it has no vertices).
Note
Using a specialization of the Sutherland-Hodgeman clipping algorithm for axis aligned planes

References axom::primal::BoundingBox< T, NDIMS >::contains(), axom::primal::BoundingBox< T, NDIMS >::getMax(), axom::primal::BoundingBox< T, NDIMS >::getMin(), axom::primal::BoundingBox< T, NDIMS >::intersectsWith(), and axom::utilities::swap().

◆ clip() [2/9]

template<typename T , axom::primal::PolygonArray ARRAY_TYPE, int MAX_VERTS>
AXOM_SUPPRESS_HD_WARN AXOM_HOST_DEVICE Polygon<T, 2, ARRAY_TYPE, MAX_VERTS> axom::primal::clip ( const Polygon< T, 2, ARRAY_TYPE, MAX_VERTS > &  subjectPolygon,
const Polygon< T, 2, ARRAY_TYPE, MAX_VERTS > &  clipPolygon,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 2D subject polygon against a clip polygon in 2D, returning their geometric intersection as a polygon.

This function clips the subject polygon by the planes obtained from the clip polygon's edges (normals point inward). Clipping the subject polygon by each plane gives the polygon above that plane. Clipping the polygon by a plane involves finding new vertices at the intersection of the polygon edges and the plane, and removing vertices from the polygon that are below the plane.

Parameters
[in]subjectPolygonThe subject polygon
[in]clipPolygonThe clip polygon
[in]epsThe tolerance for plane point orientation. Defaults to 1.e-10.
[in]tryFixOrientationIf true, takes each shape with a negative signed area and swaps the order of some vertices in that shape to try to obtain a nonnegative signed area. Defaults to false.
Returns
A polygon of the subject polygon clipped against the clip polygon.
Note
Function is based off the Sutherland–Hodgman algorithm.
Warning
Polygons with static array types must have enough vertices preallocated for the output polygon. It is mandatory that MAX_VERTS >= subjectPolygon.numVertices() + clipPolygon.numVertices() for the output polygon with the largest possible vertex count. Otherwise, if there is not enough preallocated vertices, output polygon will have missing vertices.
See also
axom::primal::Polygon::addVertex(), axom::StaticArray::push_back() for behavior when there is not enough preallocated vertices.
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned Polygon will have a non-positive and/or unexpected area.
If tryFixOrientation flag is false and some of the shapes have a negative signed area, the returned Polygon will have a non-positive and/or unexpected area.

◆ clip() [3/9]

template<typename T >
AXOM_HOST_DEVICE Polyhedron<T, 3> axom::primal::clip ( const Hexahedron< T, 3 > &  hex,
const Tetrahedron< T, 3 > &  tet,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 3D hexahedron against a tetrahedron in 3D, returning the geometric intersection of the hexahedron and the tetrahedron as a polyhedron.

This function clips the hexahedron by the 4 planes obtained from the tetrahedron's faces (normals point inward). Clipping the hexahedron/polyhedron by each plane gives the polyhedron above that plane. Clipping the polyhedron by a plane involves finding new vertices at the intersection of the polyhedron edges and the plane, removing vertices from the polyhedron that are below the plane, and redefining the neighbors for each vertex (a vertex is a neighbor of another vertex if there is an edge between them).

Parameters
[in]hexThe hexahedron to clip
[in]tetThe tetrahedron to clip against
[in]epsThe epsilon value
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
A polyhedron of the hexahedron clipped against the tetrahedron.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned Polyhedron will have a non-positive and/or unexpected volume.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned Polyhedron will have a non-positive and/or unexpected volume.

◆ clip() [4/9]

template<typename T >
AXOM_HOST_DEVICE Polyhedron<T, 3> axom::primal::clip ( const Tetrahedron< T, 3 > &  tet,
const Hexahedron< T, 3 > &  hex,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 3D hexahedron against a tetrahedron in 3D, returning the geometric intersection of the hexahedron and the tetrahedron as a polyhedron.

This function clips the hexahedron by the 4 planes obtained from the tetrahedron's faces (normals point inward). Clipping the hexahedron/polyhedron by each plane gives the polyhedron above that plane. Clipping the polyhedron by a plane involves finding new vertices at the intersection of the polyhedron edges and the plane, removing vertices from the polyhedron that are below the plane, and redefining the neighbors for each vertex (a vertex is a neighbor of another vertex if there is an edge between them).

Parameters
[in]tetThe tetrahedron to clip against
[in]hexThe hexahedron to clip
[in]epsThe epsilon value
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
A polyhedron of the hexahedron clipped against the tetrahedron.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned Polyhedron will have a non-positive and/or unexpected volume.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned Polyhedron will have a non-positive and/or unexpected volume.

References clip().

◆ clip() [5/9]

template<typename T >
AXOM_HOST_DEVICE Polyhedron<T, 3> axom::primal::clip ( const Octahedron< T, 3 > &  oct,
const Tetrahedron< T, 3 > &  tet,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 3D octahedron against a tetrahedron in 3D, returning the geometric intersection of the octahedron and the tetrahedron as a polyhedron.

This function clips the octahedron by the 4 planes obtained from the tetrahedron's faces (normals point inward). Clipping the octahedron/polyhedron by each plane gives the polyhedron above that plane. Clipping the polyhedron by a plane involves finding new vertices at the intersection of the polyhedron edges and the plane, removing vertices from the polyhedron that are below the plane, and redefining the neighbors for each vertex (a vertex is a neighbor of another vertex if there is an edge between them).

Parameters
[in]octThe octahedron to clip
[in]tetThe tetrahedron to clip against
[in]epsThe epsilon value
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
A polyhedron of the octahedron clipped against the tetrahedron.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned Polyhedron will have a non-positive and/or unexpected volume.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned Polyhedron will have a non-positive and/or unexpected volume.

◆ clip() [6/9]

template<typename T >
AXOM_HOST_DEVICE Polyhedron<T, 3> axom::primal::clip ( const Tetrahedron< T, 3 > &  tet,
const Octahedron< T, 3 > &  oct,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 3D octahedron against a tetrahedron in 3D, returning the geometric intersection of the octahedron and the tetrahedron as a polyhedron.

This function clips the octahedron by the 4 planes obtained from the tetrahedron's faces (normals point inward). Clipping the octahedron/polyhedron by each plane gives the polyhedron above that plane. Clipping the polyhedron by a plane involves finding new vertices at the intersection of the polyhedron edges and the plane, removing vertices from the polyhedron that are below the plane, and redefining the neighbors for each vertex (a vertex is a neighbor of another vertex if there is an edge between them).

Parameters
[in]tetThe tetrahedron to clip against
[in]octThe octahedron to clip
[in]tetThe tetrahedron to clip against
[in]epsThe epsilon value
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
A polyhedron of the octahedron clipped against the tetrahedron.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned Polyhedron will have a non-positive and/or unexpected volume.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned Polyhedron will have a non-positive and/or unexpected volume.

References clip().

◆ clip() [7/9]

template<typename T >
AXOM_HOST_DEVICE Polyhedron<T, 3> axom::primal::clip ( const Tetrahedron< T, 3 > &  tet1,
const Tetrahedron< T, 3 > &  tet2,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 3D tetrahedron against another tetrahedron in 3D, returning the geometric intersection as a polyhedron.

This function clips the first tetrahedron by the 4 planes obtained from the second tetrahedron's faces (normals point inward). Clipping the tetrahedron/polyhedron by each plane gives the polyhedron above that plane. Clipping the polyhedron by a plane involves finding new vertices at the intersection of the polyhedron edges and the plane, removing vertices from the polyhedron that are below the plane, and redefining the neighbors for each vertex (a vertex is a neighbor of another vertex if there is an edge between them).

Parameters
[in]tet1The tetrahedron to clip
[in]tet2The tetrahedron to clip against
[in]epsThe epsilon value
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
A polyhedron of the tetrahedron clipped against the other tetrahedron.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned Polyhedron will have a non-positive and/or unexpected volume.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned Polyhedron will have a non-positive and/or unexpected volume.

◆ clip() [8/9]

template<typename T >
AXOM_HOST_DEVICE Polyhedron<T, 3> axom::primal::clip ( const Tetrahedron< T, 3 > &  tet,
const Plane< T, 3 > &  plane,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 3D tetrahedron against the half-space defined by a plane and returns the resulting polyhedron.

This function clips a tetrahedron against the half-space defined by a plane. This involves finding new vertices at the intersection of the polyhedron edges and the plane, removing vertices from the polyhedron that are below the plane, and redefining the neighbors for each vertex (a vertex is a neighbor of another vertex if there is an edge between them).

Parameters
[in]tetThe tetrahedron to clip
[in]planeThe plane defining the half-space used to clip the tetrahedron
[in]epsThe tolerance for plane point orientation
[in]tryFixOrientationIf true and the tetrahedron has a negative signed volume, swaps the order of some vertices in the tetrathedron to try to obtain a nonnegative signed volume. Defaults to false.
Returns
The polyhedron obtained from clipping a tetrahedron against the half-space defined by a plane.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.
Warning
tryFixOrientation flag does not guarantee the tetrahedron's vertex order will be valid. It is the responsiblity of the caller to pass a tetrahedron with a valid vertex order. Otherwise, the returned polyhedron will have a non-positive and/or unexpected volume.
If the tryFixOrientation flag is false and the tetrahedron has a negative signed volume, the returned polyhedron will have a non-positive and/or unexpected volume.

◆ clip() [9/9]

template<typename T >
AXOM_HOST_DEVICE Polyhedron<T, 3> axom::primal::clip ( const Plane< T, 3 > &  plane,
const Tetrahedron< T, 3 > &  tet,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Clips a 3D tetrahedron against the half-space defined by a plane and returns the resulting polyhedron.

This function clips a tetrahedron against the half-space defined by a plane. This involves finding new vertices at the intersection of the polyhedron edges and the plane, removing vertices from the polyhedron that are below the plane, and redefining the neighbors for each vertex (a vertex is a neighbor of another vertex if there is an edge between them).

Parameters
[in]planeThe plane defining the half-space used to clip the tetrahedron
[in]tetThe tetrahedron to clip
[in]epsThe tolerance for plane point orientation
[in]tryFixOrientationIf true and the tetrahedron has a negative signed volume, swaps the order of some vertices in the tetrathedron to try to obtain a nonnegative signed volume. Defaults to false.
Returns
The polyhedron obtained from clipping a tetrahedron against the half-space defined by a plane.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.
Warning
tryFixOrientation flag does not guarantee the tetrahedron's vertex order will be valid. It is the responsiblity of the caller to pass a tetrahedron with a valid vertex order. Otherwise, the returned polyhedron will have a non-positive and/or unexpected volume.
If the tryFixOrientation flag is false and the tetrahedron has a negative signed volume, the returned polyhedron will have a non-positive and/or unexpected volume.

References clip().

◆ merge_boxes() [1/2]

template<typename T , int NDIMS>
OrientedBoundingBox<T, NDIMS> axom::primal::merge_boxes ( const OrientedBoundingBox< T, NDIMS > &  l,
const OrientedBoundingBox< T, NDIMS > &  r 
)

Creates an oriented bounding box which contains the passed in OBBs.

Call srand() to initialize the random number generator before using this function.

Parameters
[in]lleft obb
[in]rright obb

References axom::primal::OrientedBoundingBox< T, NDIMS >::contains(), and axom::primal::OrientedBoundingBox< T, NDIMS >::vertices().

◆ merge_boxes() [2/2]

template<typename T , int NDIMS>
BoundingBox<T, NDIMS> axom::primal::merge_boxes ( const BoundingBox< T, NDIMS > &  l,
const BoundingBox< T, NDIMS > &  r 
)

Constructor. Creates a bounding box which contains the passed in bounding boxes.

Parameters
[in]lleft bb
[in]rright bb

References axom::primal::BoundingBox< T, NDIMS >::addBox().

◆ compute_bounding_box() [1/7]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box ( const Triangle< T, NDIMS > &  tri)

Creates a bounding box around a Triangle.

Parameters
[in]triThe Triangle

◆ compute_bounding_box() [2/7]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box ( const Quadrilateral< T, NDIMS > &  quad)

Creates a bounding box around a Quadrilateral.

Parameters
[in]quadThe Quadrilateral

◆ compute_bounding_box() [3/7]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box ( const Octahedron< T, NDIMS > &  oct)

Creates a bounding box around an Octahedron.

Parameters
[in]octThe Octahedron

◆ compute_bounding_box() [4/7]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box ( const Hexahedron< T, NDIMS > &  hex)

Creates a bounding box around a Hexahedron.

Parameters
[in]hexThe Hexahedron

◆ compute_bounding_box() [5/7]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box ( const Polyhedron< T, NDIMS > &  poly)

◆ compute_bounding_box() [6/7]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box ( const Tetrahedron< T, NDIMS > &  tet)

Creates a bounding box around a Tetrahedron.

Parameters
[in]tetThe tetrahedron

◆ compute_bounding_box() [7/7]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box ( const Polygon< T, NDIMS > &  poly)

Creates a bounding box around a Polygon.

Parameters
[in]polyThe polygon

References axom::primal::BoundingBox< T, NDIMS >::addPoint(), and axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::numVertices().

◆ sector_area()

template<typename T >
T axom::primal::sector_area ( const primal::BezierCurve< T, 2 > &  curve)

Calculates the sector area of a planar, nonrational Bezier Curve.

The sector area is the area between the curve and the origin. The equation and derivation is described in: Ueda, K. "Signed area of sectors between spline curves and the origin" IEEE International Conference on Information Visualization, 1999.

References axom::primal::BezierCurve< T, NDIMS >::getOrder(), axom::primal::BezierCurve< T, NDIMS >::isRational(), and SLIC_ASSERT.

◆ sector_centroid()

template<typename T >
primal::Point<T, 2> axom::primal::sector_centroid ( const primal::BezierCurve< T, 2 > &  curve)

Calculates the sector centroid of a planar, nonrational Bezier Curve.

This is the centroid of the region between the curve and the origin. The equation and derivation are generalizations of: Ueda, K. "Signed area of sectors between spline curves and the origin" IEEE International Conference on Information Visualization, 1999.

References axom::primal::BezierCurve< T, NDIMS >::getOrder(), axom::primal::BezierCurve< T, NDIMS >::isRational(), and SLIC_ASSERT.

◆ area()

template<typename T >
T axom::primal::area ( const primal::CurvedPolygon< T, 2 > &  poly,
double  tol = 1e-8 
)

◆ centroid()

template<typename T >
primal::Point<T, 2> axom::primal::centroid ( const primal::CurvedPolygon< T, 2 > &  poly,
double  tol = 1e-8 
)

◆ evaluate_scalar_line_integral() [1/2]

template<typename Lambda , typename T , int NDIMS>
double axom::primal::evaluate_scalar_line_integral ( const primal::CurvedPolygon< T, NDIMS >  cpoly,
Lambda &&  scalar_integrand,
int  npts 
)

Evaluate a line integral along the boundary of a CurvedPolygon object on a scalar field.

The line integral is evaluated on each curve in the CurvedPolygon, and added together to represent the total integral. The Polygon need not be connected. Uses 1D Gaussian quadrature generated by MFEM.

Parameters
[in]cpolythe CurvedPolygon object
[in]scalar_integrandthe lambda function representing the integrand. Must accept a Point<T, NDIM> as input and return a double
[in]nptsthe number of quadrature points to evaluate the line integral on each edge of the CurvedPolygon
Returns
the value of the integral

References axom::primal::CurvedPolygon< T, NDIMS >::numEdges(), and axom::mint::SEGMENT.

◆ evaluate_scalar_line_integral() [2/2]

template<typename Lambda , typename T , int NDIMS>
double axom::primal::evaluate_scalar_line_integral ( const primal::BezierCurve< T, NDIMS > &  c,
Lambda &&  scalar_integrand,
int  npts 
)

Evaluate a line integral on a single Bezier curve on a scalar field.

Parameters
[in]cthe Bezier curve object
[in]scalar_integrandthe lambda function representing the integrand. Must accept a Point<T, NDIMS> as input, and return a double.
[in]nptsthe number of quadrature nodes
Returns
the value of the integral

References axom::mint::SEGMENT.

◆ evaluate_vector_line_integral() [1/2]

template<typename Lambda , typename T , int NDIMS>
double axom::primal::evaluate_vector_line_integral ( const primal::CurvedPolygon< T, NDIMS >  cpoly,
Lambda &&  vector_integrand,
int  npts 
)

Evaluate a line integral along the boundary of a CurvedPolygon object on a vector field.

The line integral is evaluated on each curve in the CurvedPolygon, and added together to represent the total integral. The Polygon need not be connected. Uses 1D Gaussian quadrature generated by MFEM.

Parameters
[in]cpolythe CurvedPolygon object
[in]vector_integrandthe lambda function representing the integrand. Must accept a Point<T, NDIM> as input and return a Vector<double, NDIM>
[in]nptsthe number of quadrature points to evaluate the line integral on each edge of the CurvedPolygon
Returns
the value of the integral

References axom::primal::CurvedPolygon< T, NDIMS >::numEdges(), and axom::mint::SEGMENT.

◆ evaluate_vector_line_integral() [2/2]

template<typename Lambda , typename T , int NDIMS>
double axom::primal::evaluate_vector_line_integral ( const primal::BezierCurve< T, NDIMS > &  c,
Lambda &&  vector_integrand,
int  npts 
)

Evaluate a line integral on a single Bezier curve on a vector field.

Parameters
[in]cthe Bezier curve object
[in]vector_integrandthe lambda function representing the integrand. Must accept a Point<T, NDIMS> as input, and return a Vector<double, NDIMS>.
[in]nptsthe number of quadrature nodes
Returns
the value of the integral

References axom::mint::SEGMENT.

◆ evaluate_area_integral()

template<class Lambda , typename T >
double axom::primal::evaluate_area_integral ( const primal::CurvedPolygon< T, 2 >  cpoly,
Lambda &&  integrand,
int  npts_Q,
int  npts_P = 0 
)

Evaluate an integral on the interior of a CurvedPolygon object.

See above definition for details.

Parameters
[in]csthe array of Bezier curve objects that bound the region
[in]integrandthe lambda function representing the integrand. Must accept a 2D point as input and return a double
[in]npts_Qthe number of quadrature points to evaluate the line integral
[in]npts_Pthe number of quadrature points to evaluate the antiderivative
Returns
the value of the integral

References axom::utilities::min(), axom::primal::CurvedPolygon< T, NDIMS >::numEdges(), and axom::mint::SEGMENT.

◆ in_curved_polygon()

template<typename T >
bool axom::primal::in_curved_polygon ( const Point< T, 2 > &  query,
const CurvedPolygon< T, 2 > &  cpoly,
bool  useNonzeroRule = true,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Robustly determine if query point is interior to a curved polygon.

Parameters
[in]queryThe query point to test
[in]cpolyThe CurvedPolygon object to test for containment
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Determines containment using the (rounded) winding number with respect to the given curved polygon. This algorithm is robust, as the winding number is rounded in the final in/out determination. Different protocols determine containment from the winding number differently. Nonzero Rule: If the winding number is nonzero, the point is interior. Even/Odd rule: If the winding number is odd, it is interior. Exterior otherwise.

Returns
A boolean value indicating containment.

References winding_number().

◆ in_polygon()

template<typename T >
bool axom::primal::in_polygon ( const Point< T, 2 > &  query,
const Polygon< T, 2 > &  poly,
bool  includeBoundary = false,
bool  useNonzeroRule = true,
double  edge_tol = 1e-8 
)

Determines containment for a point in a polygon.

Parameters
[in]queryThe query point to test
[in]polyThe Polygon object to test for containment
[in]includeBoundaryIf true, points on the boundary are considered interior.
[in]useNonzeroRuleIf false, use even/odd protocol for inclusion
[in]edge_tolThe distance at which a point is considered on the boundary

Determines containment using the winding number with respect to the given polygon. Different protocols determine containment from the winding number differently. Nonzero Rule: If the winding number is nonzero, the point is interior. Even/Odd rule: If the winding number is odd, it is interior. Exterior otherwise.

Returns
boolean value indicating containment.

References winding_number().

◆ in_polyhedron()

template<typename T >
bool axom::primal::in_polyhedron ( const Point< T, 3 > &  query,
const Polyhedron< T, 3 > &  poly,
bool  includeBoundary = false,
bool  useNonzeroRule = true,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Determines containment for a point in a polygon.

Parameters
[in]queryThe query point to test
[in]polyThe Polyhedron object to test for containment
[in]includeBoundaryIf true, points on the boundary are considered interior.
[in]useNonzeroRuleIf false, use even/odd protocol for inclusion
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSThe tolerance level for collinearity

Determines containment using the winding number with respect to the given polygon. Different protocols determine containment from the winding number differently. Nonzero Rule: If the winding number is nonzero, the point is interior. Even/Odd rule: If the winding number is odd, it is interior. Exterior otherwise.

Returns
boolean value indicating containment.

References winding_number().

◆ in_sphere() [1/4]

template<typename T >
bool axom::primal::in_sphere ( const Point< T, 2 > &  q,
const Point< T, 2 > &  p0,
const Point< T, 2 > &  p1,
const Point< T, 2 > &  p2,
double  EPS = 1e-8 
)
inline

Tests whether a query point lies inside a 2D triangle's circumcircle.

A triangle's circumcircle is the unique circle (i.e. a 2-sphere) that passes through each of its three vertices.

Parameters
[in]qthe query point
[in]p0the first vertex of the triangle
[in]p1the second vertex of the triangle
[in]p2the third vertex of the triangle
[in]EPStolerance for determining if q is on the boundary. Default: 1e-8.
Returns
true if the point is inside the circumcircle, false if it is on the circle's boundary or outside the circle

References axom::numerics::determinant(), and axom::utilities::isNearlyEqual().

◆ in_sphere() [2/4]

template<typename T >
bool axom::primal::in_sphere ( const Point< T, 2 > &  q,
const Triangle< T, 2 > &  tri,
double  EPS = 1e-8 
)
inline

Tests whether a query point lies inside a 2D triangle's circumcircle.

Parameters
[in]qthe query point
[in]trithe triangle
[in]EPStolerance for determining if q is on the boundary. Default: 1e-8.
See also
in_sphere

References in_sphere().

◆ in_sphere() [3/4]

template<typename T >
bool axom::primal::in_sphere ( const Point< T, 3 > &  q,
const Point< T, 3 > &  p0,
const Point< T, 3 > &  p1,
const Point< T, 3 > &  p2,
const Point< T, 3 > &  p3,
double  EPS = 1e-8 
)
inline

Tests whether a query point lies inside a 3D tetrahedron's circumsphere.

A tetrahedron's circumsphere is the unique sphere that passes through each of its four vertices.

Parameters
[in]qthe query point
[in]p0the first vertex of the tetrahedron
[in]p1the second vertex of the tetrahedron
[in]p2the third vertex of the tetrahedron
[in]p3the fourth vertex of the tetrahedron
[in]EPStolerance for determining if q is on the boundary. Default: 1e-8.
Returns
true if the point is inside the circumsphere, false if it is on the sphere's boundary or outside the sphere

References axom::numerics::determinant(), and axom::utilities::isNearlyEqual().

◆ in_sphere() [4/4]

template<typename T >
bool axom::primal::in_sphere ( const Point< T, 3 > &  q,
const Tetrahedron< T, 3 > &  tet,
double  EPS = 1e-8 
)
inline

Tests whether a query point lies inside a 3D tetrahedron's circumsphere.

Parameters
[in]qthe query point
[in]tetthe tetrahedron
[in]EPStolerance for determining if q is on the boundary. Default: 1e-8.
See also
in_sphere

References in_sphere().

◆ intersect() [1/25]

template<typename T >
AXOM_HOST_DEVICE bool axom::primal::intersect ( const Triangle< T, 3 > &  t1,
const Triangle< T, 3 > &  t2,
bool  includeBoundary = false,
double  EPS = 1E-08 
)

Tests if 3D Triangles t1 and t2 intersect.

Parameters
[in]t1The first triangle
[in]t2The second triangle
[in]includeBoundaryIndicates if boundaries should be considered when detecting intersections (default: false)
[in]EPSTolerance for determining intersections (default: 1E-8)
Returns
status true iff t1 intersects with t2, otherwise, false.

If parameter includeBoundary is false (default), this function will return true if the interior of t1 intersects the interior of t2. To include triangle boundaries in intersections, specify includeBoundary as true.

◆ intersect() [2/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 2 > &  t1,
const Triangle< T, 2 > &  t2,
bool  includeBoundary = false,
double  EPS = 1E-08 
)

Tests if 2D Triangles t1 and t2 intersect.

Parameters
[in]t1The first triangle
[in]t2The second triangle
[in]includeBoundaryIndicates if boundaries should be considered when detecting intersections (default: false)
[in]EPSTolerance for determining intersections (default: 1E-8)
Returns
status true iff t1 intersects with t2, otherwise, false.

If parameter includeBoundary is false (default), this function will return true if the interior of t1 intersects the interior of t2. To include triangle boundaries in intersections, specify includeBoundary as true.

◆ intersect() [3/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 3 > &  tri,
const BoundingBox< T, 3 > &  bb 
)

Determines if a triangle and a bounding box intersect.

Parameters
[in]triuser-supplied triangle (with three vertices).
[in]bbuser-supplied axis aligned bounding box.
Returns
true iff tri intersects with bb, otherwise, false.

◆ intersect() [4/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 3 > &  tri,
const Ray< T, 3 > &  ray 
)

Determines if a 3D triangle intersects a 3D ray.

Parameters
[in]triA 3D triangle
[in]rayA 3D ray
Returns
true iff tri intersects with ray, otherwise, false.

◆ intersect() [5/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 3 > &  tri,
const Ray< T, 3 > &  ray,
T &  t 
)

Determines if a 3D triangle intersects a 3D ray.

Parameters
[in]triA 3D triangle
[in]rayA 3D ray
[out]tIntersection point of tri and R, w.r.t. parametrization of R
Note
If there is an intersection, the intersection point is: R.at(t)
Returns
true iff tri intersects with ray, otherwise, false.

◆ intersect() [6/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 3 > &  tri,
const Ray< T, 3 > &  ray,
T &  t,
Point< double, 3 > &  p 
)

Determines if a 3D triangle intersects a 3D ray.

Parameters
[in]triA 3D triangle
[in]rayA 3D ray
[out]tIntersection point of tri and R, w.r.t. parametrization of R
[out]pIntersection point of tri and R, in barycentric coordinates relative to tri.
Note
If there is an intersection, the intersection point is: R.at(t)
Returns
true iff tri intersects with ray, otherwise, false.
Note
t and p only valid when function returns true

References axom::primal::Point< T, NDIMS >::array().

◆ intersect() [7/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 3 > &  tri,
const Segment< T, 3 > &  seg 
)

Determines if a 3D triangle intersects a 3D segment.

Parameters
[in]triA 3D triangle
[in]segA 3D line segment
Returns
true iff tri intersects with seg, otherwise, false.

◆ intersect() [8/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 3 > &  tri,
const Segment< T, 3 > &  seg,
T &  t 
)

Determines if a 3D triangle intersects a 3D segment.

Parameters
[in]triA 3D triangle
[in]segA 3D line segment
[out]tIntersection point of tri and seg, w.r.t. seg's parametrization
Returns
true iff tri intersects with seg, otherwise, false.

◆ intersect() [9/25]

template<typename T >
bool axom::primal::intersect ( const Triangle< T, 3 > &  tri,
const Segment< T, 3 > &  seg,
T &  t,
Point< double, 3 > &  p 
)

Determines if a 3D triangle intersects a 3D segment.

Parameters
[in]triA 3D triangle
[in]segA 3D line segment
[out]tIntersection point of tri and seg, w.r.t. seg's parametrization
[out]pIntersection point of tri and R, in barycentric coordinates relative to tri.
Note
If there is an intersection, the intersection point pt is: pt = seg.source() + t * ( seg.dest() - seg.target() )
Returns
true iff tri intersects with seg, otherwise, false.
Note
t and p only valid when function returns true

References axom::primal::Point< T, NDIMS >::array().

◆ intersect() [10/25]

template<typename T >
bool axom::primal::intersect ( const Ray< T, 2 > &  R,
const Segment< T, 2 > &  S,
T &  ray_param,
T &  seg_param,
const T  EPS = 1e-8 
)

Computes the intersection of the given ray, R, with the segment, S.

Parameters
[in]Rthe specified ray
[in]Sthe segment to check
[out]ray_paramparametric coordinate of intersection along R, valid only if status=true.
[out]seg_paramparametric coordinate of intersection along S, valid only if status=true.
[in]EPStolerance for intersection tests
Returns
status true iff R intersects with S, otherwise, false.
See also
primal::Ray
primal::Segment

◆ intersect() [11/25]

template<typename T >
bool axom::primal::intersect ( const Ray< T, 2 > &  R,
const Segment< T, 2 > &  S,
T &  ray_param 
)

Computes the intersection of the given ray, R, with the segment, S.

Parameters
[in]Rthe specified ray
[in]Sthe segment to check
[out]ray_paramparametric coordinate of intersection along R, valid only if status=true
Note
If you need to specify a tolerance for the intersection tests, please use the overload of this function with two [OUT] parameters (ray_param and seg_param)
Returns
status true iff R intersects with S, otherwise, false.
See also
primal::Ray
primal::Segment

References intersect().

◆ intersect() [12/25]

template<typename T >
bool axom::primal::intersect ( const Ray< T, 2 > &  R,
const Segment< T, 2 > &  S,
Point< T, 2 > &  ip,
const T  EPS = 1e-8 
)

Computes the intersection of the given ray, R, with the segment, S.

Parameters
[in]Rthe specified ray
[in]Sthe segment to check
[out]ipthe intersection point on S, valid only if status=true.
[in]EPStolerance for intersection tests
Returns
status true iff R intersects with S, otherwise, false.
See also
primal::Ray
primal::Segment
primal::Point

References axom::primal::Ray< T, NDIMS >::at().

◆ intersect() [13/25]

template<typename T , int DIM>
AXOM_HOST_DEVICE bool axom::primal::intersect ( const Ray< T, DIM > &  R,
const BoundingBox< T, DIM > &  bb,
Point< T, DIM > &  ip 
)

Computes the intersection of the given ray, R, with the Box, bb.

Parameters
[in]Rthe specified ray
[in]bbthe user-supplied axis-aligned bounding box
[out]ipthe intersection point where R intersects bb.
Returns
status true iff bb intersects with R, otherwise, false.
See also
primal::Ray
primal::Segment
primal::BoundingBox
Note
Computes Ray Box intersection using the slab method from pg 180 of Real Time Collision Detection by Christer Ericson.

◆ intersect() [14/25]

template<typename T , int DIM>
bool axom::primal::intersect ( const Segment< T, DIM > &  S,
const BoundingBox< T, DIM > &  bb,
T &  tmin,
T &  tmax,
const double &  EPS = 1e-8 
)

Computes the intersection of the given segment, S, with the Box, bb. If an intersection is found, output parameter ip contains an intersection point.

Returns
status true iff bb intersects with S, otherwise, false.
Note
The intersection between segment S and box bb intersect, will, in general, be a along a (1D) subset of segment S. One variant of this function returns the two parametric coordinates of the intersections along S found while determining if there is a valid intersection. Another variant returns an intersection point along S Specifically, it is the point of smallest parametric coordinate that is contained in bb (i.e. with parameter tmin). These are only valid when the function returns true

Computes Segment-Box intersection using the slab method from pg 180 of Real Time Collision Detection by Christer Ericson. This variant returns the two parametric coordinates of the intersection segment as OUT parameters

References axom::primal::Segment< T, NDIMS >::length().

◆ intersect() [15/25]

template<typename T , int DIM>
bool axom::primal::intersect ( const Segment< T, DIM > &  S,
const BoundingBox< T, DIM > &  bb,
Point< T, DIM > &  ip,
const double &  EPS = 1e-8 
)

This variant returns a point within the intersection as an OUT parameters.

References axom::primal::Segment< T, NDIMS >::at(), and intersect().

◆ intersect() [16/25]

template<typename T , int DIM>
bool axom::primal::intersect ( const Segment< T, DIM > &  S,
const BoundingBox< T, DIM > &  bb 
)

References intersect().

◆ intersect() [17/25]

template<typename T , int DIM>
AXOM_HOST_DEVICE bool axom::primal::intersect ( const BoundingBox< T, DIM > &  bb1,
const BoundingBox< T, DIM > &  bb2 
)

Determines if two axis aligned bounding boxes intersect.

Parameters
[in]bb1user-supplied axis aligned bounding box.
[in]bb2user-supplied axis aligned bounding box.
Returns
true iff bb1 intersects with bb2, otherwise, false.

References axom::primal::BoundingBox< T, NDIMS >::intersectsWith().

◆ intersect() [18/25]

template<typename T , int DIM>
bool axom::primal::intersect ( const Sphere< T, DIM > &  s1,
const Sphere< T, DIM > &  s2,
double  TOL = 1.e-9 
)

Determines if two spheres intersect.

Parameters
[in]s1user-supplied sphere object to check for intersection.
[in]s2user-supplied sphere object to check for intersection.
[in]TOLtolerance used for intersection check (optional)
Note
If TOL is not supplied, the default is 1.e-9.
Returns
status true iff s1 intersects with s2, otherwise, false.

References axom::primal::Sphere< T, NDIMS >::intersectsWith().

◆ intersect() [19/25]

template<typename T >
bool axom::primal::intersect ( const OrientedBoundingBox< T, 1 > &  b1,
const OrientedBoundingBox< T, 1 > &  b2 
)

◆ intersect() [20/25]

template<typename T >
bool axom::primal::intersect ( const OrientedBoundingBox< T, 2 > &  b1,
const OrientedBoundingBox< T, 2 > &  b2 
)

Determines if a 2D OBB intersects a 2D OBB.

Parameters
[in]b1A 2D OrientedBoundingBox
[in]b2A 2D OrientedBoundingBox
Returns
true iff b1 intersects with b2, otherwise, false.

◆ intersect() [21/25]

template<typename T >
bool axom::primal::intersect ( const OrientedBoundingBox< T, 3 > &  b1,
const OrientedBoundingBox< T, 3 > &  b2,
double  EPS = 1E-4 
)

Determines if a 3D OBB intersects a 3D OBB.

Parameters
[in]b1A 3D OrientedBoundingBox
[in]b2A 3D OrientedBoundingBox
[in]EPSerror tolerance for intersection
Returns
true iff b1 intersects with b2, otherwise, false.

◆ intersect() [22/25]

template<typename T >
bool axom::primal::intersect ( const BezierCurve< T, 2 > &  c1,
const BezierCurve< T, 2 > &  c2,
std::vector< T > &  sp,
std::vector< T > &  tp,
double  tol = 1E-8 
)

Tests if two Bezier Curves c1 and c2 intersect.

Returns
status true iff c1 intersects c2, otherwise false.
Parameters
[in]c1the first BezierCurve, parametrized in [0,1)
[in]c2the second BezierCurve, parametrized in [0,1)
[out]spvector of parameter space intersection points for c1
[out]tpvector of parameter space intersection points for c2
[in]toltolerance parameter for determining if a curve can be approximated by a line segment.
Returns
True if the curves intersect, false otherwise. Intersection parameters are stored in sp and tp

Finds all intersection points between the two curves.

Note
This function assumes two dimensional curves in a plane.
This function assumes that the curves are in general position. Specifically, we assume that all intersections are at points and that the curves don't overlap.
This function assumes the all intersections have multiplicity one, i.e. there are no points at which the curves and their derivatives both intersect. Thus, the function does not find tangencies.
This function assumes that the curves are half-open, i.e. they contain their first endpoint, but not their last endpoint. Thus, the curves do not intersect at \( s==1 \) or at \( t==1 \).

References axom::primal::BezierCurve< T, NDIMS >::getOrder().

◆ intersect() [23/25]

template<typename T >
AXOM_HOST_DEVICE bool axom::primal::intersect ( const Plane< T, 3 > &  p,
const BoundingBox< T, 3 > &  bb,
bool  checkOverlaps = false,
double  EPS = 1E-08 
)

Determines if a 3D plane intersects a 3D bounding box. By default (checkOverlaps is false), checks if |s| <= r, where "s" is the distance of the bounding box center to the plane, and "r" is the projected radius of the bounding box along the line parallel to the plane normal and going through the box center. If checkOverlaps is true, checks if |s| < r, where the bounding box overlaps both half spaces of the plane.

Parameters
[in]pA 3D plane
[in]bbA 3D bounding box
[in]checkOverlapsIf true, checks if bounding box overlaps both halfspaces of the plane. Otherwise, overlap of both halfspaces is not guaranteed. Default is false.
[in]EPStolerance parameter for determining if "s" is just within min/max of "r".
Returns
true iff plane intersects with bounding box, otherwise, false.
Note
Uses method from pg 164 of Real Time Collision Detection by Christer Ericson.

◆ intersect() [24/25]

template<typename T , int DIM>
AXOM_HOST_DEVICE bool axom::primal::intersect ( const Plane< T, DIM > &  plane,
const Segment< T, DIM > &  seg,
T &  t 
)

Determines if a plane intersects a segment.

Parameters
[in]planeA plane
[in]segA line segment
[out]tIntersection point of plane and seg, w.r.t. seg's parametrization
Note
If there is an intersection, the intersection point pt is: pt = seg.at(t)
Returns
true iff plane intersects with seg, otherwise, false.
Note
t is only valid when function returns true
Uses method from pg 176 of Real Time Collision Detection by Christer Ericson.

◆ intersect() [25/25]

template<typename T >
AXOM_HOST_DEVICE bool axom::primal::intersect ( const Plane< T, 3 > &  p,
const Tetrahedron< T, 3 > &  tet,
Polygon< T, 3 > &  intersection 
)

Determines if a 3D plane intersects a tetrahedron.

Parameters
[in]pA 3D plane
[in]tetA 3D tetrahedron
[out]intersectionA polygon containing the intersection.
Returns
true if plane intersects with tetrahedron, otherwise, false.
Note
If no intersection is found, the output polygon will be empty. If the plane intersects at a tetrahedron vertex, the polygon will contain duplicated points.

◆ intersection_volume() [1/5]

template<typename T >
AXOM_HOST_DEVICE T axom::primal::intersection_volume ( const Hexahedron< T, 3 > &  hex,
const Tetrahedron< T, 3 > &  tet,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Finds the absolute (unsigned) intersection volume between a hexahedron and a tetrahedron.

Parameters
[in]hexThe hexahedron
[in]tetThe tetrahedron
[in]epsThe tolerance for determining the intersection
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
Intersection volume between the hexahedron and tetrahedron
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned volume may be zero and/or unexpected.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned volume of intersection may be zero and/or unexpected.

References clip().

◆ intersection_volume() [2/5]

template<typename T >
AXOM_HOST_DEVICE T axom::primal::intersection_volume ( const Tetrahedron< T, 3 > &  tet,
const Hexahedron< T, 3 > &  hex,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Finds the absolute (unsigned) intersection volume between a tetrahedron and a hexahedron.

Parameters
[in]hexThe tetrahedron
[in]tetThe hexahedron
[in]epsThe tolerance for determining the intersection
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
Intersection volume between the tetrahedron and hexahedron
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned volume may be zero and/or unexpected.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned volume of intersection may be zero and/or unexpected.

References intersection_volume().

◆ intersection_volume() [3/5]

template<typename T >
AXOM_HOST_DEVICE T axom::primal::intersection_volume ( const Octahedron< T, 3 > &  oct,
const Tetrahedron< T, 3 > &  tet,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Finds the absolute (unsigned) intersection volume between a octahedron and a tetrahedron.

Parameters
[in]octThe octahedron
[in]tetThe tetrahedron
[in]epsThe tolerance for determining the intersection
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
Intersection volume between the octahedron and tetrahedron
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned volume may be zero and/or unexpected.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned volume of intersection may be zero and/or unexpected.

References clip().

◆ intersection_volume() [4/5]

template<typename T >
AXOM_HOST_DEVICE T axom::primal::intersection_volume ( const Tetrahedron< T, 3 > &  tet,
const Octahedron< T, 3 > &  oct,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Finds the absolute (unsigned) intersection volume between a tetrahedron and a octahedron.

Parameters
[in]octThe tetrahedron
[in]tetThe octahedron
[in]epsThe tolerance for determining the intersection
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
Intersection volume between the tetrahedron and octahedron
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned volume may be zero and/or unexpected.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned volume of intersection may be zero and/or unexpected.

References intersection_volume().

◆ intersection_volume() [5/5]

template<typename T >
AXOM_HOST_DEVICE T axom::primal::intersection_volume ( const Tetrahedron< T, 3 > &  tet1,
const Tetrahedron< T, 3 > &  tet2,
double  eps = 1.e-10,
bool  tryFixOrientation = false 
)

Finds the absolute (unsigned) intersection volume between a tetrahedron and another tetrahedron.

Parameters
[in]tet1The tetrahedron
[in]tet2The other tetrahedron
[in]epsThe tolerance for determining the intersection
[in]tryFixOrientationIf true, takes each shape with a negative signed volume and swaps the order of some vertices in that shape to try to obtain a nonnegative signed volume. Defaults to false.
Returns
Intersection volume between the tetrahedra
Warning
tryFixOrientation flag does not guarantee the shapes' vertex orders will be valid. It is the responsiblity of the caller to pass shapes with a valid vertex order. Otherwise, if the shapes have invalid vertex orders, the returned volume may be zero and/or unexpected.
If tryFixOrientation flag is false and some of the shapes have a negative signed volume, the returned volume of intersection may be zero and/or unexpected.

References clip().

◆ is_convex()

template<typename T >
bool axom::primal::is_convex ( const Polygon< T, 2 > &  poly,
double  EPS = 1e-8 
)

Determines if a polygon defined by ordered vertices is convex.

Parameters
[in]polyThe polygon

A 2D polygon is convex when every line that does not contain an edge intersects the shape at most twice. Checks whether for each pair of vertices P[i-1]P[i+1], the point P[i] and (P[0] or P[n]) lie on the same side of the line connecting them.

Algorithm adapted from: Y. L. Ma, W.T. Hewitt. "Point inversion and projection for NURBS curve and surface: Control polygon approach" Computer Aided Geometric Design 20(2):79-99, May 2003.

Note
Only defined in 2D
Returns
A boolean value indicating convexity

References axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::numVertices(), ON_BOUNDARY, and orientation().

◆ orientation() [1/2]

template<typename T >
int axom::primal::orientation ( const Point< T, 3 > &  p,
const Triangle< T, 3 > &  tri,
double  EPS = 1e-9 
)
inline

Computes the orientation of a point p with respect to an oriented triangle tri.

Parameters
[in]pthe query point
[in]trian oriented triangle
[in]EPSa tolerance for determining if p and tri are coplanar
Returns
The orientation of p with respect to tri
Note
The triangle lies in a plane that divides space into the positive half-space and the negative half-space. The triangle's normal vector points in the positive half-space. The return value of this routine can be one of the following:
  • ON_BOUNDARY, when p is coplanar with tri (within tolerance EPS)
  • ON_POSITIVE_SIDE, when p lies in the positive half-space
  • ON_NEGATIVE_SIDE, when p lies in the negative half-space
See also
OrientationResult

References axom::numerics::determinant(), axom::utilities::isNearlyEqual(), ON_BOUNDARY, ON_NEGATIVE_SIDE, and ON_POSITIVE_SIDE.

◆ orientation() [2/2]

template<typename T >
int axom::primal::orientation ( const Point< T, 2 > &  p,
const Segment< T, 2 > &  seg,
double  EPS = 1e-9 
)
inline

Computes the orientation of a point p with respect to an oriented segment.

Parameters
[in]pthe query point
[in]segan oriented segment
[in]EPSa tolerance for determining if p and seg are coplanar
Returns
The orientation of p with respect to seg
Note
The segment lies in a plane that divides space into the positive half-space and the negative half-space. The segment's normal vector points in the positive half-space. The return value can be one of the following:
  • ON_BOUNDARY, when p is coplanar with seg (within tolerance EPS)
  • ON_POSITIVE_SIDE, when p lies in the positive half-space
  • ON_NEGATIVE_SIDE, when p lies in the negative half-space
See also
OrientationResult

References axom::numerics::determinant(), axom::utilities::isNearlyEqual(), ON_BOUNDARY, ON_NEGATIVE_SIDE, and ON_POSITIVE_SIDE.

◆ split()

template<typename Tp , int NDIMS = 3>
void axom::primal::split ( const Octahedron< Tp, NDIMS > &  oct,
axom::Array< Tetrahedron< Tp, NDIMS >> &  out 
)

Splits an Octahedron into eight Tetrahedrons.

Template Parameters
Tpthe coordinate type, such double or float
NDIMSthe number of spatial dimensions (must be 3).
Parameters
[in]octThe Octahedron to split
[out]outThe axom::Array of Tetrahedron objects; the fragments of oct are appended to out.
Precondition
NDIMS == 3

The tets are produced by putting a vertex at the centroid of the oct and drawing an edge from each vertex to the centroid.

References SLIC_ASSERT.

◆ squared_distance() [1/6]

double axom::primal::squared_distance ( const double *  A,
const double *  B,
int  N 
)
inline

Computes the squared distance from point A to point B, represented by arrays of length N.

Parameters
[in]Asource point
[in]Bend point.
[in]Nlength of A and B.
Returns
the squared distance from point A to point B. If N < 1, return 0.
Precondition
A and B are arrays of at least length N.

◆ squared_distance() [2/6]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE double axom::primal::squared_distance ( const Point< T, NDIMS > &  A,
const Point< T, NDIMS > &  B 
)
inline

Computes the squared distance from point A to point B.

Parameters
[in]Asource point
[in]Bend point.
Returns
the squared distance from point A to point B.

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

◆ squared_distance() [3/6]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE double axom::primal::squared_distance ( const Point< T, NDIMS > &  P,
const BoundingBox< T, NDIMS > &  B 
)
inline

Computes the minimum squared distance from a query point, P, to a given axis-aligned bounding box B.

Parameters
[in]Pthe query point.
[in]Bthe axis-aligned bounding box.
Returns
the squared distance from P to the closest point on box B or axom::numerics::floating_point_limits<T>::max() if B is invalid.

References axom::utilities::clampVal(), axom::primal::BoundingBox< T, NDIMS >::contains(), axom::primal::BoundingBox< T, NDIMS >::getMax(), axom::primal::BoundingBox< T, NDIMS >::getMin(), axom::primal::BoundingBox< T, NDIMS >::isValid(), axom::utilities::max(), and squared_distance().

◆ squared_distance() [4/6]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE double axom::primal::squared_distance ( const BoundingBox< T, NDIMS > &  A,
const BoundingBox< T, NDIMS > &  B 
)
inline

Computes the minimum squared distance between 2 axis-aligned boxes.

Parameters
[in]Athe first axis-aligned bounding box.
[in]Bthe second axis-aligned bounding box. If the boxes overlap, the minimum distance is zero.
Returns
the squared distance between the closest points on A and B or axom::numerics::floating_point_limits<T>::max() if either box is invalid.

References axom::primal::BoundingBox< T, NDIMS >::getMax(), axom::primal::BoundingBox< T, NDIMS >::getMin(), axom::primal::BoundingBox< T, NDIMS >::isValid(), axom::utilities::max(), and axom::primal::Vector< T, NDIMS >::squared_norm().

◆ squared_distance() [5/6]

template<typename T , int NDIMS>
double axom::primal::squared_distance ( const Point< T, NDIMS > &  P,
const Segment< T, NDIMS > &  S 
)
inline

Computes the minimum squared distance from a query point, P, to the given segment, S.

Parameters
[in]Pthe query point.
[in]Sthe input segment.
Returns
the minimum squared-distance from P to the segment S.

References squared_distance().

◆ squared_distance() [6/6]

template<typename T , int NDIMS>
double axom::primal::squared_distance ( const Point< T, NDIMS > &  P,
const Triangle< T, NDIMS > &  tri 
)
inline

Computes the minimum squared distance from a query point, P, to the closest point on the given triangle.

Parameters
[in]Pthe query point.
[in]trithe supplied triangle.
Returns
the squared distance from P to the closest point on the triangle T.

References squared_distance().

◆ winding_number() [1/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 2 > &  q,
const Segment< T, 2 > &  s,
double  edge_tol = 1e-8 
)

◆ winding_number() [2/14]

template<typename T >
int axom::primal::winding_number ( const Point< T, 2 > &  q,
const Triangle< T, 2 > &  tri,
bool  includeBoundary = false,
double  edge_tol = 1e-8 
)

References winding_number().

◆ winding_number() [3/14]

template<typename T >
int axom::primal::winding_number ( const Point< T, 2 > &  R,
const Polygon< T, 2 > &  P,
bool &  isOnEdge,
bool  includeBoundary,
double  edge_tol 
)

Computes the winding number for a 2D point wrt a 2D polygon.

Parameters
[in]RThe query point to test
[in]PThe Polygon object to test for containment
[in]includeBoundaryIf true, points on the boundary are considered interior
[in]isOnEdgeAn optional return parameter if the point is on the boundary
[in]edge_tolThe distance at which a point is considered on the boundary

Uses an adapted ray-casting approach that counts quarter-rotation of vertices around the query point. Current policy is to return 1 on edges without strict inclusion, 0 on edges with strict inclusion.

The polygon is assumed to be closed, so the winding number is an integer

Directly uses algorithm in Kai Hormann, Alexander Agathos, "The point in polygon problem for arbitrary polygons" Computational Geometry, Volume 20, Issue 3, 2001,

Returns
The integer winding number

References axom::numerics::determinant(), axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::numVertices(), and squared_distance().

◆ winding_number() [4/14]

template<typename T >
int axom::primal::winding_number ( const Point< T, 2 > &  R,
const Polygon< T, 2 > &  P,
bool  includeBoundary = false,
double  edge_tol = 1e-8 
)

Computes the winding number for a 2D point wrt a 2D polygon.

Parameters
[in]RThe query point to test
[in]PThe Polygon object to test for containment
[in]includeBoundaryIf true, points on the boundary are considered interior
[in]edge_tolThe distance at which a point is considered on the boundary

Computes the integer winding number for a polygon without an additional return parameter for whether the point is on the boundary.

Returns
The integer winding number

References winding_number().

◆ winding_number() [5/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 2 > &  q,
const BezierCurve< T, 2 > &  c,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Computes the GWN for a 2D point wrt a 2D Bezier curve.

Parameters
[in]queryThe query point to test
[in]cThe Bezier curve object
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Computes the GWN using a recursive, bisection algorithm that constructs a polygon with the same integer WN as the curve closed with a linear segment. The generalized WN of the closing line is then subtracted from the integer WN to return the GWN of the original curve.

Nearly-linear Bezier curves are the base case for recursion.

See Algorithm 2 in Jacob Spainhour, David Gunderman, and Kenneth Weiss. 2024. Robust Containment Queries over Collections of Rational Parametric Curves via Generalized Winding Numbers. ACM Trans. Graph. 43, 4, Article 38 (July 2024)

Returns
The GWN.

References axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::addVertex(), axom::primal::BezierCurve< T, NDIMS >::boundingBox(), axom::primal::BoundingBox< T, NDIMS >::contains(), axom::primal::BoundingBox< T, NDIMS >::expand(), axom::primal::BezierCurve< T, NDIMS >::getOrder(), axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::numVertices(), and winding_number().

◆ winding_number() [6/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 2 > &  q,
const NURBSCurve< T, 2 > &  n,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Computes the GWN for a 2D point wrt a 2D NURBS curve.

Parameters
[in]queryThe query point to test
[in]nThe NURBS curve object
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Computes the GWN by decomposing into rational Bezier curves and summing the resulting GWNs. Far-away curves can be evaluated without decomposition using direct formula.

Returns
The GWN.

References axom::primal::NURBSCurve< T, NDIMS >::boundingBox(), axom::primal::BoundingBox< T, NDIMS >::contains(), axom::primal::BoundingBox< T, NDIMS >::expand(), axom::primal::NURBSCurve< T, NDIMS >::extractBezier(), axom::primal::NURBSCurve< T, NDIMS >::getDegree(), axom::primal::NURBSCurve< T, NDIMS >::getNumControlPoints(), and winding_number().

◆ winding_number() [7/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 2 > &  q,
const CurvedPolygon< T, 2 > &  cpoly,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Computes the GWN for a 2D point wrt to a 2D curved polygon.

Parameters
[in]queryThe query point to test
[in]cpolyThe CurvedPolygon object
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Computes the GWN for the curved polygon by summing the GWN for each curved edge

Returns
The GWN.

References axom::primal::CurvedPolygon< T, NDIMS >::numEdges(), and winding_number().

◆ winding_number() [8/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 2 > &  q,
const axom::Array< BezierCurve< T, 2 >> &  carray,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Computes the GWN for a 2D point wrt to a collection of 2D Bezier curves.

Parameters
[in]queryThe query point to test
[in]carrayThe array of Bezier curves
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Sums the GWN at query for each curved edge

Returns
The GWN.

References AXOM_UNUSED_VAR, and winding_number().

◆ winding_number() [9/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 2 > &  q,
const axom::Array< NURBSCurve< T, 2 >> &  narray,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Computes the GWN for a 2D point wrt to a collection of 2D NURBS curves.

Parameters
[in]queryThe query point to test
[in]narrayThe array of NURBS curves
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Sums the GWN at query for each curved edge

Returns
The GWN.

References AXOM_UNUSED_VAR, and winding_number().

◆ winding_number() [10/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 3 > &  q,
const Triangle< T, 3 > &  tri,
bool &  isOnFace,
const double  edge_tol = 1e-8,
const double  EPS = 1e-8 
)

Computes the GWN for a 3D point wrt a 3D triangle.

Parameters
[in]queryThe query point to test
[in]triThe 3D Triangle object
[in]isOnFaceAn optional return parameter if the point is on the triangle
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Computes the GWN as the solid angle modulo 4pi using the formula from Oosterom, Strackee, "The Solid Angle of a Plane Triangle" IEEE Transactions on Biomedical Engineering, Vol BME-30, No. 2, February 1983 with extra adjustments if the triangle takes up a full octant

Returns
The GWN.

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

◆ winding_number() [11/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 3 > &  q,
const Triangle< T, 3 > &  tri,
const double  edge_tol = 1e-8,
const double  EPS = 1e-8 
)

Computes the GWN for a 3D point wrt a 3D triangle.

Parameters
[in]queryThe query point to test
[in]triThe 3D Triangle object
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances

Computes the GWN for the triangle without an additional return parameter

Returns
The GWN.

References winding_number().

◆ winding_number() [12/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 3 > &  q,
const Polygon< T, 3 > &  poly,
bool &  isOnFace,
const double  edge_tol = 1e-8,
const double  EPS = 1e-8 
)

Computes the GWN for a 3D point wrt a 3D planar polygon.

Parameters
[in]queryThe query point to test
[in]polyThe Polygon object
[in]isOnFaceReturn variable to show if the point is on the polygon
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances
Precondition
Assumes the polygon is planar. Otherwise, a meaningless value is returned.

Triangulates the polygon and computes the triangular GWN for each component

Returns
The GWN.

References axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::numVertices(), and winding_number().

◆ winding_number() [13/14]

template<typename T >
double axom::primal::winding_number ( const Point< T, 3 > &  q,
const Polygon< T, 3 > &  poly,
const double  edge_tol = 1e-8,
const double  EPS = 1e-8 
)

Computes the GWN for a 3D point wrt a 3D planar polygon.

Parameters
[in]queryThe query point to test
[in]polyThe Polygon object
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances
Precondition
Assumes the polygon is planar. Otherwise, a meaningless value is returned.

Computes the GWN for the polygon without an additional return parameter

Returns
The GWN.

References winding_number().

◆ winding_number() [14/14]

template<typename T >
int axom::primal::winding_number ( const Point< T, 3 > &  query,
const Polyhedron< T, 3 > &  poly,
bool  includeBoundary = false,
double  edge_tol = 1e-8,
double  EPS = 1e-8 
)

Computes the winding number for a 3D point wrt a 3D convex polyhedron.

Parameters
[in]queryThe query point to test
[in]polyThe Polyhedron object
[in]includeBoundaryIf true, points on the boundary are considered interior.
[in]edge_tolThe physical distance level at which objects are considered indistinguishable
[in]EPSMiscellaneous numerical tolerance level for nonphysical distances
Precondition
Expects the polyhedron to be convex and closed so that the returned value is an integer.

Computes the faces of the polyhedron and computes the GWN for each. The sum is then rounded to the nearest integer, as the shape is assumed to be closed.

Returns
The integer winding number.

References axom::primal::Polygon< T, NDIMS, ARRAY_TYPE, MAX_VERTS >::addVertex(), axom::Array< T, DIM, SPACE >::data(), axom::primal::Polyhedron< T, NDIMS >::getFaces(), axom::primal::Polyhedron< T, NDIMS >::hasNeighbors(), axom::primal::Polyhedron< T, NDIMS >::numVertices(), SLIC_ASSERT, and winding_number().