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

Classes

class  BoundingBox
 
class  NeighborCollection
 Represents a collection of neighbor relations between vertices. 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  NumericArray
 A simple statically sized array of data with component-wise operators. 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 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  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
typedef Point< double, 2 > Point2D
 
typedef Point< double, 3 > Point3D
 
Pre-defined Vector types
typedef Vector< double, 2 > Vector2D
 
typedef Vector< double, 3 > Vector3D
 

Enumerations

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

Functions

template<typename T , int NDIMS>
bool operator== (const BoundingBox< T, NDIMS > &lhs, const BoundingBox< T, NDIMS > &rhs)
 Free functions implementing comparison and arithmetic operators. More...
 
template<typename T , int SIZE>
bool operator== (const NumericArray< T, SIZE > &lhs, const NumericArray< T, SIZE > &rhs)
 Free functions implementing comparison and arithmetic operators. More...
 
template<typename T , int SIZE>
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>
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)
 Performs component-wise subtraction of two numeric arrays. 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 >
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 >
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...
 
template<typename T , int NDIMS>
std::ostream & operator<< (std::ostream &os, const Polygon< T, NDIMS > &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 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 , int NDIMS>
Vector< T, NDIMS > operator* (const Vector< T, NDIMS > &vec, const T scalar)
 Free functions involving vectors. More...
 
template<typename T , int NDIMS>
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>
Vector< T, NDIMS > operator+ (const Vector< T, NDIMS > &vec1, const Vector< T, NDIMS > &vec2)
 Adds vectors A, B and stores the result into a new vector C. More...
 
template<typename T , int NDIMS>
Vector< T, NDIMS > operator- (const Vector< T, NDIMS > &vec1, const Vector< T, NDIMS > &vec2)
 Subtracts vectors A, B and stores the result into a new vector C. 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_HOST_DEVICE Polyhedron< T, 3 > clip (const Octahedron< T, 3 > &oct, const Tetrahedron< T, 3 > &tet, double eps=1.e-10)
 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 , int NDIMS>
AXOM_HOST_DEVICE Point< T, NDIMS > closest_point (const Point< T, NDIMS > &P, const Triangle< T, NDIMS > &tri, int *loc=nullptr)
 Computes the closest point from a point, P, to a given triangle. More...
 
template<typename T , int NDIMS>
Point< T, NDIMS > closest_point (const Point< T, NDIMS > &pt, const OrientedBoundingBox< T, NDIMS > &obb)
 Computes the closest point from a point to a given OBB. 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 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 Polyhedron< T, NDIMS > &poly)
 Creates a bounding box around a Polyhedron. 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)
 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)
 Tests whether a query point lies inside a 3D tetrahedron's circumsphere. More...
 
template<typename T >
int orientation (const Point< T, 3 > &p, const Triangle< T, 3 > &tri)
 Computes the orientation of the given point, p, with respect to a supplied oriented triangle. More...
 
template<typename T >
int orientation (const Point< T, 2 > &p, const Segment< T, 2 > &seg)
 Computes the orientation of the given point, p, with respect to a supplied 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>
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>
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>
bool operator== (const Segment< T, NDIMS > &lhs, const Segment< T, NDIMS > &rhs)
 Equality comparison operator for Segment. More...
 
template<typename T , int NDIMS>
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 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>
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>
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>
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 , int NDIMS>
bool intersect (const BezierCurve< T, NDIMS > &c1, const BezierCurve< T, NDIMS > &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 >
AXOM_HOST_DEVICE bool intersect (const Plane< T, 3 > &plane, const Segment< T, 3 > &seg, T &t)
 Determines if a 3D plane intersects a 3D segment. More...
 

Typedef Documentation

◆ Point2D

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

◆ Point3D

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

◆ Vector2D

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

◆ Vector3D

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

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

Function Documentation

◆ operator==() [1/7]

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.

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

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<<() [1/14]

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

Overloaded output operator for bounding boxes.

◆ operator==() [2/7]

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

Free functions implementing comparison and arithmetic operators.

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

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

◆ operator==() [3/7]

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.

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.

Checks if two numeric arrays are component-wise equal.

◆ 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/8]

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/8]

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/8]

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<<() [2/14]

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.

◆ operator==() [4/7]

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

Free functions implementing comparison and arithmetic operators.

Checks if two numeric arrays are component-wise equal.

◆ operator+() [2/4]

template<typename T , int SIZE>
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*() [4/8]

template<typename T , int SIZE>
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-() [3/6]

template<typename T , int SIZE>
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<<() [3/14]

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::Point< T, NDIMS >::print().

◆ operator==() [5/7]

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<<() [4/14]

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<<() [5/14]

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

◆ make_plane() [1/4]

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/4]

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.

◆ make_plane() [3/4]

template<typename T >
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() [4/4]

template<typename T >
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==() [6/7]

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>
bool axom::primal::operator!= ( const Point< T, NDIMS > &  lhs,
const Point< T, NDIMS > &  rhs 
)

◆ operator<<() [6/14]

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.

◆ operator<<() [7/14]

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

Overloaded output operator for polygons.

Free functions implementing Polygon's operators.

◆ operator<<() [8/14]

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.

◆ operator<<() [9/14]

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::Vector< T, NDIMS >::print().

◆ operator==() [7/7]

template<typename T , int NDIMS>
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>
bool axom::primal::operator!= ( const Segment< T, NDIMS > &  lhs,
const Segment< T, NDIMS > &  rhs 
)

Inequality comparison operator for Segment.

Inequality operator for segments.

◆ operator<<() [10/14]

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.

◆ operator<<() [11/14]

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

◆ operator<<() [12/14]

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::Point< T, NDIMS >::print().

◆ operator<<() [13/14]

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.

◆ operator+() [3/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-() [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>
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*() [5/8]

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.

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

Scalar multiplication of vector; Scalar on rhs.

◆ operator*() [6/8]

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>
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]

◆ operator<<() [14/14]

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

◆ operator*() [7/8]

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

Free functions involving vectors.

Scalar multiplication of vector; Scalar on rhs.

◆ operator*() [8/8]

template<typename T , int NDIMS>
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+() [4/4]

template<typename T , int NDIMS>
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-() [6/6]

template<typename T , int NDIMS>
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 \)

◆ clip() [1/2]

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::Polygon< T, NDIMS >::addVertex(), 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().

Referenced by axom::quest::InOutOctree< DIM >::generateIndex(), and axom::quest::IntersectionShaper::setExecPolicy().

◆ clip() [2/2]

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 
)

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
Returns
A polyhedron of the octahedron clipped against the tetrahedron.
Note
Function is based off clipPolyhedron() in Mike Owen's PolyClipper.

◆ closest_point() [1/2]

template<typename T , int NDIMS>
AXOM_HOST_DEVICE Point<T, NDIMS> axom::primal::closest_point ( const Point< T, NDIMS > &  P,
const Triangle< T, NDIMS > &  tri,
int *  loc = nullptr 
)
inline

Computes the closest point from a point, P, to a given triangle.

Parameters
[in]Pthe query point
[in]triuser-supplied triangle.
[out]locint pointer to store location of closest point (optional).
Returns
cp the closest point from a point P and a triangle.
Note
If the optional int pointer is supplied for loc, the method returns the location of the closest point, which is illustrated in the schematic diagram below and encoded as follows:
  • loc \( \in [0,2] \), loc corresponds to the triangle node index
  • loc \( \in [-3,-1] \), abs(loc) corresponds to an edge
  • loc >= 3, loc is on a triangle face
*
*            2
*           /\
*    (-3)--/  \--(-2)
*         /    \
*        /_ _ _ \
*       0   |    1
*           |
*         (-1)
*
* 
Precondition
NDIMS==2 || NDIMS==3
Note
Implementation is based on "Real Time Collision Detection, Chapter 5.1.5 Closest Point on Triangle to Point".

References A, B, C, and axom::primal::Vector< T, NDIMS >::dot_product().

Referenced by axom::quest::SignedDistance< NDIMS, ExecSpace >::computeDistances(), and squared_distance().

◆ closest_point() [2/2]

template<typename T , int NDIMS>
Point<T, NDIMS> axom::primal::closest_point ( const Point< T, NDIMS > &  pt,
const OrientedBoundingBox< T, NDIMS > &  obb 
)
inline

Computes the closest point from a point to a given OBB.

Parameters
[in]ptthe query pt.
[in]obbuser-supplied oriented bounding box.
Returns
cp the closest point from a point pt and an OBB.

References axom::primal::OrientedBoundingBox< T, NDIMS >::getAxes(), axom::primal::OrientedBoundingBox< T, NDIMS >::getCentroid(), axom::primal::OrientedBoundingBox< T, NDIMS >::getExtents(), and axom::primal::OrientedBoundingBox< T, NDIMS >::toLocal().

◆ 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/3]

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

◆ compute_bounding_box() [2/3]

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

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

◆ compute_bounding_box() [3/3]

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

◆ in_sphere() [1/2]

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 
)
inline

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

A triangle's circumscircle 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
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().

◆ in_sphere() [2/2]

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

◆ intersect() [1/24]

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.

Referenced by axom::quest::InOutOctree< DIM >::generateIndex(), axom::quest::getMeshTriangle(), and intersect().

◆ intersect() [2/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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/24]

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

References intersect().

◆ intersect() [17/24]

template<typename T , int DIM>
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/24]

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/24]

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

◆ intersect() [20/24]

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/24]

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/24]

template<typename T , int NDIMS>
bool axom::primal::intersect ( const BezierCurve< T, NDIMS > &  c1,
const BezierCurve< T, NDIMS > &  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 \).

◆ intersect() [23/24]

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/24]

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

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

Parameters
[in]planeA 3D plane
[in]segA 3D 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.

◆ orientation() [1/2]

template<typename T >
int axom::primal::orientation ( const Point< T, 3 > &  p,
const Triangle< T, 3 > &  tri 
)
inline

Computes the orientation of the given point, p, with respect to a supplied oriented triangle.

Parameters
[in]pthe query point.
[in]trioriented triangle.
Returns
The orientation of the point with respect to the given triangle.
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 into the positive half-space. The return value of this routine can be one of the following:
  • ON_BOUNDARY, when the point is coplanar with the triangle
  • ON_POSITIVE_SIDE, when the point lies in the positive half-space
  • ON_NEGATIVE_SIDE, when the point lies in the negative half-space

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 
)
inline

Computes the orientation of the given point, p, with respect to a supplied oriented segment.

Parameters
[in]pthe query point.
[in]segthe user-supplied segment.
Returns
The orientation of the point with respect to the given segment.
Note
The return value can be one of the following:
  • ON_BOUNDARY, when the point is collinear with the points that define the segment.
  • ON_POSITIVE_SIDE, when the point is clockwise, i.e., to the right of the directed segment.
  • ON_NEGATIVE_SIDE, when the point is counter-clockwise, i.e., to the left of the directed segment.

References axom::numerics::determinant(), axom::utilities::isNearlyEqual(), ON_BOUNDARY, ON_NEGATIVE_SIDE, ON_POSITIVE_SIDE, axom::primal::Segment< T, NDIMS >::source(), and axom::primal::Segment< T, NDIMS >::target().

◆ 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 C, and SLIC_ASSERT.

◆ squared_distance() [1/5]

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
d the distance from point A to point B. If N < 1, return 0.
Precondition
A and B are arrays of at least length N.

Referenced by axom::quest::SignedDistance< NDIMS, ExecSpace >::computeDistances(), axom::primal::Octahedron< T, NDIMS >::equals(), axom::quest::InOutOctree< DIM >::generateIndex(), axom::quest::IntersectionShaper::setExecPolicy(), and squared_distance().

◆ squared_distance() [2/5]

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
d the distance from point A to point B.

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

◆ squared_distance() [3/5]

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
d the signed distance from P to the closest point on B.

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

◆ squared_distance() [4/5]

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
d the minimum distance from P on the

References axom::primal::Vector< T, NDIMS >::dot_product(), axom::primal::Segment< T, NDIMS >::source(), axom::primal::Vector< T, NDIMS >::squared_norm(), and axom::primal::Segment< T, NDIMS >::target().

◆ squared_distance() [5/5]

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
d the distance from Q to the closest point on the triangle T.

References closest_point(), and squared_distance().