AXOM
Axom provides a robust, flexible software infrastructure for the development of multi-physics applications and computational tools.
|
Classes | |
class | BoundingBox |
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 | 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 | 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 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 Polygon< T, NDIMS > &poly) |
Overloaded output operator for polygons. 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 Segment< T, NDIMS > &seg) |
Overloaded output operator for Segment. 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 , int NDIMS> | |
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 > | |
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... | |
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> | |
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> | |
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> | |
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> | |
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> | |
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> | |
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> | |
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> | |
std::ostream & | operator<< (std::ostream &os, const Sphere< T, NDIMS > &s) |
template<typename T , int NDIMS> | |
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> | |
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> | |
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> | |
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... | |
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, Point< T, 2 > &ip) |
Computes the intersection of the given ray, R, with the segment, S. More... | |
template<typename T , int DIM> | |
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 Intersection Routines | |
template<typename T , int DIM> | |
bool | intersect (const Segment< T, DIM > &S, const BoundingBox< T, DIM > &bb, Point< T, DIM > &ip) |
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... | |
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... | |
typedef Point<double, 2> axom::primal::Point2D |
typedef Point<double, 3> axom::primal::Point3D |
typedef Vector<double, 2> axom::primal::Vector2D |
typedef Vector<double, 3> axom::primal::Vector3D |
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().
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.
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const BoundingBox< T, NDIMS > & | bb | ||
) |
Overloaded output operator for bounding boxes.
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.
[in] | lhs | numeric array instance on the left-hand side. |
[in] | rhs | numeric array instance on the right-hand side. |
bool axom::primal::operator!= | ( | const NumericArray< T, SIZE > & | lhs, |
const NumericArray< T, SIZE > & | rhs | ||
) |
Checks if two numeric arrays are not component-wise equal.
[in] | lhs | numeric array instance on the left-hand side. |
[in] | rhs | numeric array instance on the right-hand side. |
|
inline |
Performs component-wise addition of two numeric arrays.
[in] | lhs | numeric array instance on the left-hand side. |
[in] | rhs | numeric array instance on the right-hand side. |
|
inline |
Performs component-wise subtraction of two numeric arrays.
[in] | lhs | numeric array instance on the left-hand side. |
[in] | rhs | numeric array instance on the right-hand side. |
|
inline |
Unary negation of a numeric array instance.
[in] | arr | numeric array instance on the left-hand side. |
|
inline |
Scalar multiplication a numeric array; Scalar on rhs.
[in] | arr | numeric array instance. |
[in] | scalar | user-supplied scalar. |
|
inline |
Scalar multiplication a numeric array; Scalar on lhs.
[in] | scalar | user-supplied scalar. |
[in] | arr | numeric array instance. |
|
inline |
Component-wise multiplication of NumericArrays.
[in] | lhs | numeric array instance on the left-hand side. |
[in] | rhs | numeric array instance on the right-hand side. |
|
inline |
Component-wise division of NumericArrays.
[in] | lhs | numeric array instance on the left-hand side. |
[in] | rhs | numeric array instance on the right-hand side. |
|
inline |
Scalar division of NumericArray; Scalar on rhs.
[in] | arr | numeric array instance |
[in] | scalar | user-supplied scalar |
|
inline |
Coordinate-wise absolute value on the NumericArray.
[in] | arr | numeric array instance |
References axom::utilities::abs().
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const NumericArray< T, SIZE > & | arr | ||
) |
Overloaded output operator for numeric arrays.
[in] | os | C++ output stream |
[in] | arr | numeric array instance. |
|
inline |
Performs component-wise subtraction of two numeric arrays.
[in] | lhs | numeric array instance on the left-hand side. |
[in] | rhs | numeric array instance on the right-hand side. |
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().
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
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().
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.
[in] | pts | array of points |
[in] | n | number of points |
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const Plane< T, NDIMS > & | p | ||
) |
bool axom::primal::operator== | ( | const Point< T, NDIMS > & | lhs, |
const Point< T, NDIMS > & | rhs | ||
) |
Equality comparison operator for points.
bool axom::primal::operator!= | ( | const Point< T, NDIMS > & | lhs, |
const Point< T, NDIMS > & | rhs | ||
) |
Inequality comparison operator for points.
Inequality operator for points.
References A, AXOM_HOST_DEVICE, B, axom::primal::Point< T, NDIMS >::lerp(), axom::primal::Point< T, NDIMS >::make_point(), axom::primal::Point< T, NDIMS >::midpoint(), and axom::primal::Point< T, NDIMS >::print().
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const Point< T, NDIMS > & | pt | ||
) |
Overloaded output operator for points.
Free functions implementing Point's operators.
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const Polygon< T, NDIMS > & | poly | ||
) |
Overloaded output operator for polygons.
Free functions implementing Polygon's operators.
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().
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const Segment< T, NDIMS > & | seg | ||
) |
Overloaded output operator for Segment.
Free functions implementing Segments's operators.
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const Sphere< T, NDIMS > & | s | ||
) |
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().
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const Triangle< T, NDIMS > & | tri | ||
) |
Overloaded output operator for triangles.
Free functions implementing Triangle's operators.
|
inline |
Adds vectors A, B and stores the result into a new vector C.
[in] | A | vector on the left-hand side. |
[in] | B | vector on the right-hand side. |
|
inline |
Subtracts vectors A, B and stores the result into a new vector C.
[in] | A | vector on the left-hand side. |
[in] | B | vector on the right-hand side. |
|
inline |
Unary negation of a vector instance.
[in] | vec1 | vector instance to negate. |
References axom::primal::Vector< T, NDIMS >::negate().
|
inline |
Scalar multiplication of vector; Scalar on rhs.
Free functions involving vectors.
[in] | vec | vector instance. |
[in] | scalar | user-supplied scalar. |
|
inline |
Scalar multiplication of vector; Scalar on lhs.
[in] | scalar | user-supplied scalar. |
[in] | vec | vector instance. |
|
inline |
Scalar division of vector; Scalar on rhs.
[in] | vec | vector instance |
[in] |
std::ostream & axom::primal::operator<< | ( | std::ostream & | os, |
const Vector< T, NDIMS > & | vec | ||
) |
Overloaded output operator for vectors.
[in] | os | C++ output stream |
[in] | vec | vector instance |
References axom::primal::Vector< T, NDIMS >::print().
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.
[in] | tri | The triangle to clip |
[in] | bbox | The bounding box to clip against |
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().
|
inline |
Computes the closest point from a point, P, to a given triangle.
[in] | P | the query point |
[in] | tri | user-supplied triangle. |
[out] | loc | int pointer to store location of closest point (optional). |
loc
, the method returns the location of the closest point, which is illustrated in the schematic diagram below and encoded as follows: * * 2 * /\ * (-3)--/ \--(-2) * / \ * /_ _ _ \ * 0 | 1 * | * (-1) * *
References A, B, C, and axom::primal::Vector< T, NDIMS >::dot_product().
Referenced by axom::quest::SignedDistance< NDIMS >::computeDistance(), and squared_distance().
|
inline |
Computes the closest point from a point to a given OBB.
[in] | pt | the query pt. |
[in] | obb | user-supplied oriented bounding box. |
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().
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.
[in] | l | left obb |
[in] | r | right obb |
References axom::primal::OrientedBoundingBox< T, NDIMS >::contains(), and axom::primal::OrientedBoundingBox< T, NDIMS >::vertices().
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.
[in] | l | left bb |
[in] | r | right bb |
References axom::primal::BoundingBox< T, NDIMS >::addBox().
AXOM_HOST_DEVICE BoundingBox<T, NDIMS> axom::primal::compute_bounding_box | ( | const Triangle< T, NDIMS > & | tri | ) |
Creates a bounding box around a Triangle.
[in] | tri | The Triangle |
References axom::primal::BoundingBox< T, NDIMS >::addPoint().
Referenced by axom::quest::getMeshTriangle().
|
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.
[in] | q | the query point |
[in] | p0 | the first vertex of the triangle |
[in] | p1 | the second vertex of the triangle |
[in] | p2 | the third vertex of the triangle |
References axom::numerics::determinant().
|
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.
[in] | q | the query point |
[in] | p0 | the first vertex of the tetrahedron |
[in] | p1 | the second vertex of the tetrahedron |
[in] | p2 | the third vertex of the tetrahedron |
[in] | p3 | the fourth vertex of the tetrahedron |
References axom::numerics::determinant().
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.
[in] | t1 | The first triangle |
[in] | t2 | The second triangle |
[in] | includeBoundary | Indicates if boundaries should be considered when detecting intersections (default: false) |
[in] | EPS | Tolerance for determining intersections (default: 1E-8) |
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::detail::InOutOctreeValidator< DIM >::checkBlockIndexingConsistency(), axom::quest::InOutOctree< DIM >::generateIndex(), and axom::quest::getMeshTriangle().
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.
[in] | t1 | The first triangle |
[in] | t2 | The second triangle |
[in] | includeBoundary | Indicates if boundaries should be considered when detecting intersections (default: false) |
[in] | EPS | Tolerance for determining intersections (default: 1E-8) |
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.
bool axom::primal::intersect | ( | const Triangle< T, 3 > & | tri, |
const BoundingBox< T, 3 > & | bb | ||
) |
Determines if a triangle and a bounding box intersect.
[in] | tri | user-supplied triangle (with three vertices). |
[in] | bb | user-supplied axis aligned bounding box. |
bool axom::primal::intersect | ( | const Triangle< T, 3 > & | tri, |
const Ray< T, 3 > & | ray | ||
) |
Determines if a 3D triangle intersects a 3D ray.
[in] | tri | A 3D triangle |
[in] | ray | A 3D ray |
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.
[in] | tri | A 3D triangle |
[in] | ray | A 3D ray |
[out] | t | Intersection point of tri and R, w.r.t. parametrization of R |
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.
[in] | tri | A 3D triangle |
[in] | ray | A 3D ray |
[out] | t | Intersection point of tri and R, w.r.t. parametrization of R |
[out] | p | Intersection point of tri and R, in barycentric coordinates relative to tri. |
References axom::primal::Point< T, NDIMS >::array().
bool axom::primal::intersect | ( | const Triangle< T, 3 > & | tri, |
const Segment< T, 3 > & | seg | ||
) |
Determines if a 3D triangle intersects a 3D segment.
[in] | tri | A 3D triangle |
[in] | seg | A 3D line segment |
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.
[in] | tri | A 3D triangle |
[in] | seg | A 3D line segment |
[out] | t | Intersection point of tri and seg, w.r.t. seg's parametrization |
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.
[in] | tri | A 3D triangle |
[in] | seg | A 3D line segment |
[out] | t | Intersection point of tri and seg, w.r.t. seg's parametrization |
[out] | p | Intersection point of tri and R, in barycentric coordinates relative to tri. |
References axom::primal::Point< T, NDIMS >::array().
bool axom::primal::intersect | ( | const Ray< T, 2 > & | R, |
const Segment< T, 2 > & | S, | ||
Point< T, 2 > & | ip | ||
) |
Computes the intersection of the given ray, R, with the segment, S.
[in] | R | the specified ray |
[in] | S | the segment to check |
[out] | ip | the intersection point on S, valid only if status=true. |
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.
[in] | R | the specified ray |
[in] | bb | the user-supplied axis-aligned bounding box |
[out] | ip | the intersection point where R intersects bb. |
bool axom::primal::intersect | ( | const Segment< T, DIM > & | S, |
const BoundingBox< T, DIM > & | bb, | ||
Point< T, DIM > & | ip | ||
) |
Computes the intersection of the given segment, S, with the Box, bb. If an intersection is found, output parameter ip contains an intersection point.
Computes Segment-Box intersection using the slab method from pg 180 of Real Time Collision Detection by Christer Ericson.
bool axom::primal::intersect | ( | const BoundingBox< T, DIM > & | bb1, |
const BoundingBox< T, DIM > & | bb2 | ||
) |
Determines if two axis aligned bounding boxes intersect.
[in] | bb1 | user-supplied axis aligned bounding box. |
[in] | bb2 | user-supplied axis aligned bounding box. |
References axom::primal::BoundingBox< T, NDIMS >::intersectsWith().
bool axom::primal::intersect | ( | const Sphere< T, DIM > & | s1, |
const Sphere< T, DIM > & | s2, | ||
double | TOL = 1.e-9 |
||
) |
Determines if two spheres intersect.
[in] | s1 | user-supplied sphere object to check for intersection. |
[in] | s2 | user-supplied sphere object to check for intersection. |
[in] | TOL | tolerance used for intersection check (optional) |
References axom::primal::Sphere< T, NDIMS >::intersectsWith().
bool axom::primal::intersect | ( | const OrientedBoundingBox< T, 1 > & | b1, |
const OrientedBoundingBox< T, 1 > & | b2 | ||
) |
bool axom::primal::intersect | ( | const OrientedBoundingBox< T, 2 > & | b1, |
const OrientedBoundingBox< T, 2 > & | b2 | ||
) |
Determines if a 2D OBB intersects a 2D OBB.
[in] | b1 | A 2D OrientedBoundingBox |
[in] | b2 | A 2D OrientedBoundingBox |
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.
[in] | b1 | A 3D OrientedBoundingBox |
[in] | b2 | A 3D OrientedBoundingBox |
[in] | EPS | error tolerance for intersection |
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.
[in] | c1 | the first BezierCurve, parametrized in [0,1) |
[in] | c2 | the second BezierCurve, parametrized in [0,1) |
[out] | sp | vector of parameter space intersection points for c1 |
[out] | tp | vector of parameter space intersection points for c2 |
[in] | tol | tolerance parameter for determining if a curve can be approximated by a line segment. |
Finds all intersection points between the two curves.
|
inline |
Computes the orientation of the given point, p, with respect to a supplied oriented triangle.
[in] | p | the query point. |
[in] | tri | oriented triangle. |
References axom::numerics::determinant(), axom::utilities::isNearlyEqual(), ON_BOUNDARY, ON_NEGATIVE_SIDE, and ON_POSITIVE_SIDE.
|
inline |
Computes the orientation of the given point, p, with respect to a supplied oriented segment.
[in] | p | the query point. |
[in] | seg | the user-supplied 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().
|
inline |
Computes the squared distance from point A to point B, represented by arrays of length N.
[in] | A | source point |
[in] | B | end point. |
[in] | N | length of A and B. |
Referenced by axom::quest::SignedDistance< NDIMS >::computeDistance(), axom::spin::BVHTree< T, NDIMS >::contains(), axom::quest::InOutOctree< DIM >::generateIndex(), and squared_distance().
|
inline |
Computes the squared distance from point A to point B.
[in] | A | source point |
[in] | B | end point. |
References axom::primal::Vector< T, NDIMS >::squared_norm().
|
inline |
Computes the minimum squared distance from a query point, P, to a given axis-aligned bounding box B.
[in] | P | the query point. |
[in] | B | the axis-aligned bounding box. |
References axom::primal::BoundingBox< T, NDIMS >::contains(), axom::primal::BoundingBox< T, NDIMS >::getMax(), axom::primal::BoundingBox< T, NDIMS >::getMin(), and squared_distance().
|
inline |
Computes the minimum squared distance from a query point, P, to the given segment, S.
[in] | P | the query point. |
[in] | S | the input segment. |
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().
|
inline |
Computes the minimum squared distance from a query point, P, to the closest point on the given triangle.
[in] | P | the query point. |
[in] | tri | the supplied triangle. |
References closest_point(), and squared_distance().