PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
|
Geometric queries, quantities and data structures. More...
Namespaces | |
namespace | ClosestPointQueries |
This namespace contains functions to answer closest point queries. | |
namespace | DistanceQueries |
This namespace contains functions to answer distance queries. | |
namespace | IntersectionQueries |
This namespace contains functions to answer intersection queries. | |
namespace | OverlapQueries |
This namespace contains functions to answer overlap queries. | |
namespace | sdf |
Namespace for signed distance functions (SDFs) and related operations. | |
Classes | |
class | AabbKdTreeHierarchy |
Axis-aligned k-D tree hierarchy of axis-aligned bounding boxes. More... | |
class | AabbRadixTreeHierarchy |
Axis-aligned radix tree hierarchy of axis-aligned bounding boxes. More... | |
class | AxisAlignedBoundingBox |
Axis-aligned bounding box class. More... | |
class | BoundingVolumeHierarchy |
CRTP base class for BVHs. More... | |
class | HashGrid |
HashGrid is a spatial partitioning data structure that divides 3D space into a sparse grid of cells, allowing for efficient querying of point neighbours within a certain region. More... | |
class | HierarchicalHashGrid |
Spatial partitioning data structure that divides 3D space into a set of sparse grids. Allowing for efficient querying of point neighbours within a certain region. Implements [4]. More... | |
class | KdTree |
KDTree class. More... | |
struct | KdTreeNode |
Node of a KDTree. More... | |
class | TetrahedralAabbHierarchy |
Tetrahedral AABB hierarchy class. More... | |
class | TriangleAabbHierarchy2D |
Bounding volume hierarchy for triangles in 2D. More... | |
class | TriangleAabbHierarchy3D |
Bounding volume hierarchy for triangles in 3D. More... | |
Typedefs | |
using | MortonCodeType = std::uint32_t |
Type used to represent Morton codes. | |
Functions | |
template<auto kDims, auto kClusterNodes, class FCluster, class TDerivedL, class TDerivedU> | |
void | ClustersToAabbs (FCluster fCluster, Index nClusters, Eigen::DenseBase< TDerivedL > &L, Eigen::DenseBase< TDerivedU > &U) |
Computes AABBs of nClusters kDims-dimensional point clusters. | |
template<auto kDims, auto kClusterNodes, class FCluster, class TDerivedB> | |
void | ClustersToAabbs (FCluster fCluster, Index nClusters, Eigen::DenseBase< TDerivedB > &B) |
Computes AABBs of nClusters kDims-dimensional point clusters. | |
template<auto kDims, auto kElemNodes, class TDerivedX, class TDerivedE, class TDerivedL, class TDerivedU> | |
void | MeshToAabbs (Eigen::DenseBase< TDerivedX > const &X, Eigen::DenseBase< TDerivedE > const &E, Eigen::DenseBase< TDerivedL > &L, Eigen::DenseBase< TDerivedU > &U) |
Computes AABBs of nElemNodes simplex mesh elements in kDims dimensions. | |
template<auto kDims, auto kElemNodes, class TDerivedX, class TDerivedE, class TDerivedB> | |
void | MeshToAabbs (Eigen::DenseBase< TDerivedX > const &X, Eigen::DenseBase< TDerivedE > const &E, Eigen::DenseBase< TDerivedB > &B) |
Computes AABBs of nElemNodes simplex mesh elements in kDims dimensions. | |
template<math::linalg::mini::CReadableVectorizedMatrix TP1T, math::linalg::mini::CReadableVectorizedMatrix TQ1T, math::linalg::mini::CReadableVectorizedMatrix TP2T, math::linalg::mini::CReadableVectorizedMatrix TQ2T, math::linalg::mini::CReadableVectorizedMatrix TP1, math::linalg::mini::CReadableVectorizedMatrix TQ1, math::linalg::mini::CReadableVectorizedMatrix TP2, math::linalg::mini::CReadableVectorizedMatrix TQ2, class TScalar = typename TP1T::ScalarType> | |
PBAT_HOST_DEVICE auto | EdgeEdgeCcd (TP1T const &P1T, TQ1T const &Q1T, TP2T const &P2T, TQ2T const &Q2T, TP1 const &P1, TQ1 const &Q1, TP2 const &P2, TQ2 const &Q2) -> math::linalg::mini::SVector< TScalar, 3 > |
Computes the time of impact \( t^* \) and barycentric coordinates \( \mathbf{\beta} \) of the intersection point between an edge \( (P1,Q1) \) and another edge \( (P2,Q2) \) moving along a linear trajectory. | |
template<common::CIndex THashed> | |
auto | HashByXorOfPrimeMultiples () |
Returns a hash function used in [13] for 2 or 3-vectors suitable as input to HashGrid's API. | |
template<common::CIndex TIndex = Index> | |
auto | SimplexMeshBoundary (Eigen::Ref< Eigen::Matrix< TIndex, Eigen::Dynamic, Eigen::Dynamic > const > const &C, TIndex n) -> std::tuple< Eigen::Vector< TIndex, Eigen::Dynamic >, Eigen::Matrix< TIndex, Eigen::Dynamic, Eigen::Dynamic > > |
Obtains the boundary mesh of a simplex mesh. | |
PBAT_HOST_DEVICE MortonCodeType | ExpandBits (MortonCodeType v) |
Expands a 10-bit integer into 30 bits by inserting 2 zeros after each bit. | |
template<detail::CMorton3dPoint Point> | |
PBAT_HOST_DEVICE MortonCodeType | Morton3D (Point x) |
Calculates a 30-bit Morton code for the given 3D point located within the unit cube [0,1]. | |
template<math::linalg::mini::CReadableVectorizedMatrix TXT, math::linalg::mini::CReadableVectorizedMatrix TAT, math::linalg::mini::CReadableVectorizedMatrix TBT, math::linalg::mini::CReadableVectorizedMatrix TCT, math::linalg::mini::CReadableVectorizedMatrix TX, math::linalg::mini::CReadableVectorizedMatrix TA, math::linalg::mini::CReadableVectorizedMatrix TB, math::linalg::mini::CReadableVectorizedMatrix TC, class TScalar = typename TXT::ScalarType> | |
PBAT_HOST_DEVICE auto | PointTriangleCcd (TXT const &XT, TAT const &AT, TBT const &BT, TCT const &CT, TX const &X, TA const &A, TB const &B, TC const &C) -> math::linalg::mini::SVector< TScalar, 4 > |
Computes the time of impact \( t^* \) and barycentric coordinates \( \mathbf{\beta} \) of the intersection point between a point and a triangle moving along a linear trajectory. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FNodeOverlap, class FObjectOverlap, class FOnFound, auto N = 2, class TIndex = Index, auto kStackDepth = 64> | |
PBAT_HOST_DEVICE bool | Overlaps (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FNodeOverlap fNodeOverlaps, FObjectOverlap fObjectOverlaps, FOnFound fOnFound, TIndex root=0) |
Find all nodes in a branch and bound tree that overlap with a given object. | |
template<class FChildLhs, class FIsLeafLhs, class FLeafSizeLhs, class FLeafObjectLhs, class FChildRhs, class FIsLeafRhs, class FLeafSizeRhs, class FLeafObjectRhs, class FNodesOverlap, class FObjectsOverlap, class FOnFound, auto NLhs = 2, auto NRhs = 2, class TIndex = Index, auto kStackDepth = 128> | |
PBAT_HOST_DEVICE bool | Overlaps (FChildLhs fChildLhs, FIsLeafLhs fIsLeafLhs, FLeafSizeLhs fLeafSizeLhs, FLeafObjectLhs fLeafObjectLhs, FChildRhs fChildRhs, FIsLeafRhs fIsLeafRhs, FLeafSizeRhs fLeafSizeRhs, FLeafObjectRhs fLeafObjectRhs, FNodesOverlap fNodesOverlap, FObjectsOverlap fObjectsOverlap, FOnFound fOnFound, TIndex rootLhs=0, TIndex rootRhs=0) |
Computes overlaps between two branch and bound trees. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FNodesOverlap, class FObjectsOverlap, class FOnFound, auto N = 2, class TIndex = Index, auto kStackDepth = 128> | |
PBAT_HOST_DEVICE bool | SelfOverlaps (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FNodesOverlap fNodesOverlap, FObjectsOverlap fObjectsOverlap, FOnFound fOnFound, TIndex root=0) |
Compute overlapping nodes from a branch and bound tree rooted in root. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N = 2, class TScalar = Scalar, class TIndex = Index, auto kStackDepth = 64, auto kQueueSize = 8> | |
PBAT_HOST_DEVICE bool | DfsAllNearestNeighbours (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FDistanceLowerBound fLower, FDistance fDistance, FOnFound fOnFound, bool bUseBestFirstSearch=false, TScalar fUpper=std::numeric_limits< TScalar >::max(), TScalar eps=TScalar(0), TIndex root=0) |
Find distance minimizing objects in branch and bound tree rooted in root. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N = 2, class TScalar = Scalar, class TIndex = Index, auto kHeapCapacity = N * 1024> | |
PBAT_HOST_DEVICE bool | NearestToFurthestNeighbours (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FDistanceLowerBound fLower, FDistance fDistance, FOnFound fOnFound, TScalar fUpper=std::numeric_limits< TScalar >::max(), TIndex root=0) |
Traverse objects in a branch and bound tree rooted in root in increasing order of distance. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N = 2, class TScalar = Scalar, class TIndex = Index, auto kHeapCapacity = N * 1024> | |
PBAT_HOST_DEVICE bool | AllNearestNeighbours (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FDistanceLowerBound fLower, FDistance fDistance, FOnFound fOnFound, TScalar fUpper=std::numeric_limits< TScalar >::max(), TScalar eps=TScalar(0), TIndex root=0) |
Find all global distance minimizers in a branch and bound tree rooted in root. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N = 2, class TScalar = Scalar, class TIndex = Index, auto kHeapCapacity = N * 1024> | |
PBAT_HOST_DEVICE bool | KNearestNeighbours (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FDistanceLowerBound fLower, FDistance fDistance, FOnFound fOnFound, TIndex K, TScalar fUpper=std::numeric_limits< TScalar >::max(), TIndex root=0) |
Find the (sorted) K objects with the smallest distances in a branch and bound tree rooted in root. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FNodeOverlap, class FObjectOverlap, class FOnFound, auto N, class TIndex, auto kStackDepth> | |
bool | Overlaps (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FNodeOverlap fNodeOverlaps, FObjectOverlap fObjectOverlaps, FOnFound fOnFound, TIndex root) |
Find all nodes in a branch and bound tree that overlap with a given object. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N, class TScalar, class TIndex, auto kStackDepth, auto kQueueSize> | |
bool | DfsAllNearestNeighbours (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FDistanceLowerBound fLower, FDistance fDistance, FOnFound fOnFound, bool bUseBestFirstSearch, TScalar fUpper, TScalar eps, TIndex root) |
Find distance minimizing objects in branch and bound tree rooted in root. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N, class TScalar, class TIndex, auto kHeapCapacity> | |
bool | NearestToFurthestNeighbours (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FDistanceLowerBound fLower, FDistance fDistance, FOnFound fOnFound, TScalar fUpper, TIndex root) |
Traverse objects in a branch and bound tree rooted in root in increasing order of distance. | |
template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N, class TScalar, class TIndex, auto kHeapCapacity> | |
bool | KNearestNeighbours (FChild fChild, FIsLeaf fIsLeaf, FLeafSize fLeafSize, FLeafObject fLeafObject, FDistanceLowerBound fLower, FDistance fDistance, FOnFound fOnFound, TIndex K, TScalar fUpper, TIndex root) |
Find the (sorted) K objects with the smallest distances in a branch and bound tree rooted in root. | |
Geometric queries, quantities and data structures.
PBAT_HOST_DEVICE bool pbat::geometry::AllNearestNeighbours | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FDistanceLowerBound | fLower, | ||
FDistance | fDistance, | ||
FOnFound | fOnFound, | ||
TScalar | fUpper = std::numeric_limits<TScalar>::max(), | ||
TScalar | eps = TScalar(0), | ||
TIndex | root = 0 ) |
Find all global distance minimizers in a branch and bound tree rooted in root.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FDistanceLowerBound | Callable with signature TScalar(TIndex node) |
FDistance | Callable with signature TScalar(TIndex o) |
FOnFound | Callable with signature void(TIndex o, TScalar d, TIndex k) |
N | Max number of children per node |
TScalar | Type of the scalar distance |
TIndex | Type of the index |
kHeapCapacity | Maximum capacity of the min-distance-heap |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fLower | Function to compute the lower bound of the distance to node |
fDistance | Function to compute the distance to object |
fOnFound | Function to call when a nearest neighbour is found |
fUpper | Upper bound of the distance to the nearest neighbour |
eps | Epsilon to consider objects at the same distance |
root | Index of the root node to start the search from |
|
inline |
Computes AABBs of nClusters kDims-dimensional point clusters.
kDims | Spatial dimension |
kClusterNodes | Number of nodes in a cluster |
FCluster | Function with signature auto (Index) -> std::convertible_to<Matrix<kDims, kClusterNodes>> |
TDerivedB | Eigen dense expression type |
fCluster | Function to get the cluster at index c |
nClusters | Number of clusters |
B | 2*kDims x |# clusters| output AABBs |
|
inline |
Computes AABBs of nClusters kDims-dimensional point clusters.
kDims | Spatial dimension |
kClusterNodes | Number of nodes in a cluster |
FCluster | Function with signature auto (Index) -> std::convertible_to<Matrix<kDims, kClusterNodes>> |
TDerivedL | Eigen dense expression type |
TDerivedU | Eigen dense expression type |
fCluster | Function to get the cluster at index c |
nClusters | Number of clusters |
L | kDims x |# clusters| output AABBs lower bounds |
U | kDims x |# clusters| output AABBs upper bounds |
bool pbat::geometry::DfsAllNearestNeighbours | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FDistanceLowerBound | fLower, | ||
FDistance | fDistance, | ||
FOnFound | fOnFound, | ||
bool | bUseBestFirstSearch = false, | ||
TScalar | fUpper = std::numeric_limits< TScalar >::max(), | ||
TScalar | eps = TScalar(0), | ||
TIndex | root = 0 ) |
Find distance minimizing objects in branch and bound tree rooted in root.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FDistanceLowerBound | Callable with signature TScalar(TIndex node) |
FDistance | Callable with signature TScalar(TIndex o) |
FOnFound | Callable with signature void(TIndex o, TScalar d, TIndex k) |
N | Max number of children per node |
TScalar | Type of the scalar distance |
TIndex | Type of the index |
kStackDepth | Maximum depth of the traversal's stack |
kQueueSize | Maximum size of the nearest neighbour queue |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fLower | Function to compute the lower bound of the distance to node |
fDistance | Function to compute the distance to object |
fOnFound | Function to call when a nearest neighbour is found |
bUseBestFirstSearch | Use best-first search instead of pre-order traversal |
fUpper | Upper bound of the distance to the nearest neighbour |
eps | Epsilon to consider objects at the same distance |
root | Index of the root node to start the search from |
PBAT_HOST_DEVICE bool pbat::geometry::DfsAllNearestNeighbours | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FDistanceLowerBound | fLower, | ||
FDistance | fDistance, | ||
FOnFound | fOnFound, | ||
bool | bUseBestFirstSearch = false, | ||
TScalar | fUpper = std::numeric_limits< TScalar >::max(), | ||
TScalar | eps = TScalar(0), | ||
TIndex | root = 0 ) |
Find distance minimizing objects in branch and bound tree rooted in root.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FDistanceLowerBound | Callable with signature TScalar(TIndex node) |
FDistance | Callable with signature TScalar(TIndex o) |
FOnFound | Callable with signature void(TIndex o, TScalar d, TIndex k) |
N | Max number of children per node |
TScalar | Type of the scalar distance |
TIndex | Type of the index |
kStackDepth | Maximum depth of the traversal's stack |
kQueueSize | Maximum size of the nearest neighbour queue |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fLower | Function to compute the lower bound of the distance to node |
fDistance | Function to compute the distance to object |
fOnFound | Function to call when a nearest neighbour is found |
bUseBestFirstSearch | Use best-first search instead of pre-order traversal |
fUpper | Upper bound of the distance to the nearest neighbour |
eps | Epsilon to consider objects at the same distance |
root | Index of the root node to start the search from |
PBAT_HOST_DEVICE auto pbat::geometry::EdgeEdgeCcd | ( | TP1T const & | P1T, |
TQ1T const & | Q1T, | ||
TP2T const & | P2T, | ||
TQ2T const & | Q2T, | ||
TP1 const & | P1, | ||
TQ1 const & | Q1, | ||
TP2 const & | P2, | ||
TQ2 const & | Q2 ) -> math::linalg::mini::SVector<TScalar, 3> |
Computes the time of impact \( t^* \) and barycentric coordinates \( \mathbf{\beta} \) of the intersection point between an edge \( (P1,Q1) \) and another edge \( (P2,Q2) \) moving along a linear trajectory.
Solves for inexact roots (if any) in the range [0,1] of the polynomial
\[\langle \mathbf{n}(t), \mathbf{q}(t) \rangle = 0 , \]
where \( \mathbf{n}(t) = (\mathbf{q}_1(t) - \mathbf{p}_1(t)) \times (\mathbf{q}_2(t) - \mathbf{p}_2(t)) \) and \( \mathbf{q}(t) = \mathbf{p}_2(t) - \mathbf{p}_1(t) \) using polynomial root finder from [16].
See [11] and [15] for more details.
TP1T | Type of the input matrix P1T |
TQ1T | Type of the input matrix Q1T |
TP2T | Type of the input matrix P2T |
TQ2T | Type of the input matrix Q2T |
TP1 | Type of the input matrix P1 |
TQ1 | Type of the input matrix Q1 |
TP2 | Type of the input matrix P2 |
TQ2 | Type of the input matrix Q2 |
TScalar | Type of the scalar |
P1T | Matrix of the initial positions of the point P on edge 1 |
Q1T | Matrix of the initial positions of the point Q on edge 1 |
P2T | Matrix of the initial positions of the point P on edge 2 |
Q2T | Matrix of the initial positions of the point Q on edge 2 |
P1 | Matrix of the final positions of the point P on edge 1 |
Q1 | Matrix of the final positions of the point Q on edge 1 |
P2 | Matrix of the final positions of the point P on edge 2 |
Q2 | Matrix of the final positions of the point Q on edge 2 |
r[0] < 0
, otherwise r[0]
is the earliest time of impact r[1]
and r[2]
are the barycentric coordinates of the intersection point along \((P1,Q1) \) and \( (P2,Q2) \) respectively
|
inline |
Expands a 10-bit integer into 30 bits by inserting 2 zeros after each bit.
v | 10-bit integer |
|
inline |
bool pbat::geometry::KNearestNeighbours | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FDistanceLowerBound | fLower, | ||
FDistance | fDistance, | ||
FOnFound | fOnFound, | ||
TIndex | K, | ||
TScalar | fUpper = std::numeric_limits< TScalar >::max(), | ||
TIndex | root = 0 ) |
Find the (sorted) K objects with the smallest distances in a branch and bound tree rooted in root.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FDistanceLowerBound | Callable with signature TScalar(TIndex node) |
FDistance | Callable with signature TScalar(TIndex o) |
FOnFound | Callable with signature void(TIndex o, TScalar d, TIndex k) |
N | Max number of children per node |
TScalar | Type of the scalar distance |
TIndex | Type of the index |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fLower | Function to compute the lower bound of the distance to node |
fDistance | Function to compute the distance to object |
fOnFound | Function to call when a nearest neighbour is found |
K | Number of nearest neighbours to find |
fUpper | Upper bound of the distance to the nearest neighbour |
root | Index of the root node to start the search from |
PBAT_HOST_DEVICE bool pbat::geometry::KNearestNeighbours | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FDistanceLowerBound | fLower, | ||
FDistance | fDistance, | ||
FOnFound | fOnFound, | ||
TIndex | K, | ||
TScalar | fUpper = std::numeric_limits< TScalar >::max(), | ||
TIndex | root = 0 ) |
Find the (sorted) K objects with the smallest distances in a branch and bound tree rooted in root.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FDistanceLowerBound | Callable with signature TScalar(TIndex node) |
FDistance | Callable with signature TScalar(TIndex o) |
FOnFound | Callable with signature void(TIndex o, TScalar d, TIndex k) |
N | Max number of children per node |
TScalar | Type of the scalar distance |
TIndex | Type of the index |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fLower | Function to compute the lower bound of the distance to node |
fDistance | Function to compute the distance to object |
fOnFound | Function to call when a nearest neighbour is found |
K | Number of nearest neighbours to find |
fUpper | Upper bound of the distance to the nearest neighbour |
root | Index of the root node to start the search from |
|
inline |
Computes AABBs of nElemNodes simplex mesh elements in kDims dimensions.
kDims | Spatial dimension |
kElemNodes | Number of nodes in an element |
TDerivedX | Eigen dense expression type |
TDerivedE | Eigen dense expression type |
TDerivedB | Eigen dense expression type |
X | kDims x |# nodes| matrix of node positions |
E | kElemNodes x |# elements| matrix of element node indices |
B | 2*kDims x |# elements| output AABBs |
|
inline |
Computes AABBs of nElemNodes simplex mesh elements in kDims dimensions.
kDims | Spatial dimension |
kElemNodes | Number of nodes in an element |
TDerivedX | Eigen dense expression type |
TDerivedE | Eigen dense expression type |
TDerivedL | Eigen dense expression type for lower bounds |
TDerivedU | Eigen dense expression type for upper bounds |
X | kDims x |# nodes| matrix of node positions |
E | kElemNodes x |# elements| matrix of element node indices |
L | kDims x |# elements| output AABBs lower bounds |
U | kDims x |# elements| output AABBs upper bounds |
|
inline |
Calculates a 30-bit Morton code for the given 3D point located within the unit cube [0,1].
Point | Type of the point |
x | 3D point located within the unit cube [0,1] |
bool pbat::geometry::NearestToFurthestNeighbours | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FDistanceLowerBound | fLower, | ||
FDistance | fDistance, | ||
FOnFound | fOnFound, | ||
TScalar | fUpper = std::numeric_limits< TScalar >::max(), | ||
TIndex | root = 0 ) |
Traverse objects in a branch and bound tree rooted in root in increasing order of distance.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FDistanceLowerBound | Callable with signature TScalar(TIndex node) |
FDistance | Callable with signature TScalar(TIndex o) |
FOnFound | Callable with signature bool(TIndex o, TScalar d, TIndex k) |
N | Max number of children per node |
TScalar | Type of the scalar distance |
TIndex | Type of the index |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fLower | Function to compute the lower bound of the distance to node |
fDistance | Function to compute the distance to object |
fOnFound | Function to call when a nearest neighbour is found, returns true if the search should continue |
fUpper | Upper bound of the distance to the nearest neighbour |
root | Index of the root node to start the search from |
PBAT_HOST_DEVICE bool pbat::geometry::NearestToFurthestNeighbours | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FDistanceLowerBound | fLower, | ||
FDistance | fDistance, | ||
FOnFound | fOnFound, | ||
TScalar | fUpper = std::numeric_limits< TScalar >::max(), | ||
TIndex | root = 0 ) |
Traverse objects in a branch and bound tree rooted in root in increasing order of distance.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FDistanceLowerBound | Callable with signature TScalar(TIndex node) |
FDistance | Callable with signature TScalar(TIndex o) |
FOnFound | Callable with signature bool(TIndex o, TScalar d, TIndex k) |
N | Max number of children per node |
TScalar | Type of the scalar distance |
TIndex | Type of the index |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fLower | Function to compute the lower bound of the distance to node |
fDistance | Function to compute the distance to object |
fOnFound | Function to call when a nearest neighbour is found, returns true if the search should continue |
fUpper | Upper bound of the distance to the nearest neighbour |
root | Index of the root node to start the search from |
bool pbat::geometry::Overlaps | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FNodeOverlap | fNodeOverlaps, | ||
FObjectOverlap | fObjectOverlaps, | ||
FOnFound | fOnFound, | ||
TIndex | root = 0 ) |
Find all nodes in a branch and bound tree that overlap with a given object.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FNodeOverlaps | Callable with signature bool(TIndex node) |
FObjectOverlaps | Callable with signature bool(TIndex o) |
FOnFound | Callable with signature void(TIndex o, TIndex k) |
N | Max number of children per node |
TIndex | Type of the index |
kStackDepth | Maximum depth of the traversal's stack |
kQueueSize | Maximum size of the nearest neighbour queue |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fNodeOverlaps | Function to determine if a node overlaps with the query |
fObjectOverlaps | Function to determine if an object overlaps with the query |
fOnFound | Function to call when a nearest neighbour is found |
root | Index of the root node to start the search from |
PBAT_HOST_DEVICE bool pbat::geometry::Overlaps | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FNodeOverlap | fNodeOverlaps, | ||
FObjectOverlap | fObjectOverlaps, | ||
FOnFound | fOnFound, | ||
TIndex | root = 0 ) |
Find all nodes in a branch and bound tree that overlap with a given object.
FChild | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FNodeOverlaps | Callable with signature bool(TIndex node) |
FObjectOverlaps | Callable with signature bool(TIndex o) |
FOnFound | Callable with signature void(TIndex o, TIndex k) |
N | Max number of children per node |
TIndex | Type of the index |
kStackDepth | Maximum depth of the traversal's stack |
kQueueSize | Maximum size of the nearest neighbour queue |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fNodeOverlaps | Function to determine if a node overlaps with the query |
fObjectOverlaps | Function to determine if an object overlaps with the query |
fOnFound | Function to call when a nearest neighbour is found |
root | Index of the root node to start the search from |
PBAT_HOST_DEVICE bool pbat::geometry::Overlaps | ( | FChildLhs | fChildLhs, |
FIsLeafLhs | fIsLeafLhs, | ||
FLeafSizeLhs | fLeafSizeLhs, | ||
FLeafObjectLhs | fLeafObjectLhs, | ||
FChildRhs | fChildRhs, | ||
FIsLeafRhs | fIsLeafRhs, | ||
FLeafSizeRhs | fLeafSizeRhs, | ||
FLeafObjectRhs | fLeafObjectRhs, | ||
FNodesOverlap | fNodesOverlap, | ||
FObjectsOverlap | fObjectsOverlap, | ||
FOnFound | fOnFound, | ||
TIndex | rootLhs = 0, | ||
TIndex | rootRhs = 0 ) |
Computes overlaps between two branch and bound trees.
FChildLhs | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeafLhs | Callable with signature bool(TIndex node) |
FLeafSizeLhs | Callable with signature TIndex(TIndex node) |
FLeafObjectLhs | Callable with signature TIndex(TIndex node, TIndex i) |
FChildRhs | Callable with signature template <auto c> TIndex(TIndex node) |
FIsLeafRhs | Callable with signature bool(TIndex node) |
FLeafSizeRhs | Callable with signature TIndex(TIndex node) |
FLeafObjectRhs | Callable with signature TIndex(TIndex node, TIndex i) |
FNodesOverlap | Callable with signature bool(TIndex nlhs, TIndex nrhs) |
FObjectsOverlap | Callable with signature bool(TIndex olhs, TIndex orhs) |
FOnFound | Callable with signature void(TIndex olhs, TIndex orhs, TIndex k) |
NLhs | Number of children per node in the left tree |
NRhs | Number of children per node in the right tree |
TIndex | Type of the index |
kStackDepth | Maximum depth of the traversal's stack |
fChildLhs | Function to get child c of a node in the left tree. Returns the child index or -1 if no child. |
fIsLeafLhs | Function to determine if a node is a leaf node in the left tree |
fLeafSizeLhs | Function to get the number of leaf objects in a node in the left tree |
fLeafObjectLhs | Function to get the i-th leaf object in a node in the left tree |
fChildRhs | Function to get child c of a node in the right tree. Returns the child index or -1 if no child. |
fIsLeafRhs | Function to determine if a node is a leaf node in the right tree |
fLeafSizeRhs | Function to get the number of leaf objects in a node in the right tree |
fLeafObjectRhs | Function to get the i-th leaf object in a node in the right tree |
fNodesOverlap | Function to determine if two nodes overlap |
fObjectsOverlap | Function to determine if two objects overlap |
fOnFound | Function to call when an overlap is found |
rootLhs | Index of the root node in the left tree to start the search from |
rootRhs | Index of the root node in the right tree to start the search from |
PBAT_HOST_DEVICE auto pbat::geometry::PointTriangleCcd | ( | TXT const & | XT, |
TAT const & | AT, | ||
TBT const & | BT, | ||
TCT const & | CT, | ||
TX const & | X, | ||
TA const & | A, | ||
TB const & | B, | ||
TC const & | C ) -> math::linalg::mini::SVector<TScalar, 4> |
Computes the time of impact \( t^* \) and barycentric coordinates \( \mathbf{\beta} \) of the intersection point between a point and a triangle moving along a linear trajectory.
Solves for roots (if any) in the range [0,1] of the polynomial
\[\langle \mathbf{n}(t), \mathbf{q}(t) \rangle = 0 , \]
where \( \mathbf{n}(t) = (\mathbf{b}(t) - \mathbf{a}(t)) \times (\mathbf{c}(t) - \mathbf{a}(t)) \) and \( \mathbf{q}(t) = \mathbf{x}(t) - \mathbf{a}(t) \) using polynomial root finder from [16].
See [11] and [15] for more details.
TXT | Type of the input matrix XT |
TAT | Type of the input matrix AT |
TBT | Type of the input matrix BT |
TCT | Type of the input matrix CT |
TX | Type of the input matrix X |
TA | Type of the input matrix A |
TB | Type of the input matrix B |
TC | Type of the input matrix C |
TScalar | Type of the scalar |
XT | Matrix of the initial positions of the point |
AT | Matrix of the initial positions of the triangle vertex A |
BT | Matrix of the initial positions of the triangle vertex B |
CT | Matrix of the initial positions of the triangle vertex C |
X | Matrix of the final positions of the point |
A | Matrix of the final positions of the triangle vertex A |
B | Matrix of the final positions of the triangle vertex B |
C | Matrix of the final positions of the triangle vertex C |
r
containing the time of impact and the barycentric coordinates. r[0] < 0
, otherwise r[0] >= 0
is the earliest time of impact (subject to floating point error). PBAT_HOST_DEVICE bool pbat::geometry::SelfOverlaps | ( | FChild | fChild, |
FIsLeaf | fIsLeaf, | ||
FLeafSize | fLeafSize, | ||
FLeafObject | fLeafObject, | ||
FNodesOverlap | fNodesOverlap, | ||
FObjectsOverlap | fObjectsOverlap, | ||
FOnFound | fOnFound, | ||
TIndex | root = 0 ) |
Compute overlapping nodes from a branch and bound tree rooted in root.
FChild | Callable with signature template <auto c> Index(TIndex node) |
FIsLeaf | Callable with signature bool(TIndex node) |
FLeafSize | Callable with signature TIndex(TIndex node) |
FLeafObject | Callable with signature TIndex(TIndex node, TIndex i) |
FNodesOverlap | Callable with signature bool(TIndex n1, TIndex n2) |
FObjectsOverlap | Callable with signature bool(TIndex o1, TIndex o2) |
FOnFound | Callable with signature void(TIndex o1, TIndex o2, TIndex k) |
N | Max number of children per node |
TIndex | Type of the index |
kStackDepth | Maximum depth of the traversal's stack |
kQueueSize | Maximum size of the nearest neighbour queue |
fChild | Function to get child c of a node. Returns the child index or -1 if no child. |
fIsLeaf | Function to determine if a node is a leaf node |
fLeafSize | Function to get the number of leaf objects in a node |
fLeafObject | Function to get the i-th leaf object in a node |
fNodesOverlap | Function to determine if 2 nodes overlap |
fObjectsOverlap | Function to determine if 2 objects overlap |
fOnFound | Function to call when a nearest neighbour is found |
root | Index of the root node to start the search from |
auto pbat::geometry::SimplexMeshBoundary | ( | Eigen::Ref< Eigen::Matrix< TIndex, Eigen::Dynamic, Eigen::Dynamic > const > const & | C, |
TIndex | n ) -> std::tuple< Eigen::Vector<TIndex, Eigen::Dynamic>, Eigen::Matrix<TIndex, Eigen::Dynamic, Eigen::Dynamic>> |
Obtains the boundary mesh of a simplex mesh.
C.rows()==3
) and tetrahedral (C.rows()==4
) meshes.TIndex | The index type used in the mesh (default: Index ) |
C | The connectivity matrix of the mesh (i.e. the simplices) |
n | The number of vertices in the mesh. If -1, the number of vertices is computed from C. |