PhysicsBasedAnimationToolkit 0.0.10
Cross-platform C++20 library of algorithms and data structures commonly used in computer graphics research on physically-based simulation.
Loading...
Searching...
No Matches
pbat::geometry Namespace Reference

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.
 

Detailed Description

Geometric queries, quantities and data structures.

Function Documentation

◆ AllNearestNeighbours()

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

Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FDistanceLowerBoundCallable with signature TScalar(TIndex node)
FDistanceCallable with signature TScalar(TIndex o)
FOnFoundCallable with signature void(TIndex o, TScalar d, TIndex k)
NMax number of children per node
TScalarType of the scalar distance
TIndexType of the index
Parameters
kHeapCapacityMaximum capacity of the min-distance-heap
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fLowerFunction to compute the lower bound of the distance to node
fDistanceFunction to compute the distance to object
fOnFoundFunction to call when a nearest neighbour is found
fUpperUpper bound of the distance to the nearest neighbour
epsEpsilon to consider objects at the same distance
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ ClustersToAabbs() [1/2]

template<auto kDims, auto kClusterNodes, class FCluster, class TDerivedB>
void pbat::geometry::ClustersToAabbs ( FCluster fCluster,
Index nClusters,
Eigen::DenseBase< TDerivedB > & B )
inline

Computes AABBs of nClusters kDims-dimensional point clusters.

Template Parameters
kDimsSpatial dimension
kClusterNodesNumber of nodes in a cluster
FClusterFunction with signature auto (Index) -> std::convertible_to<Matrix<kDims, kClusterNodes>>
TDerivedBEigen dense expression type
Parameters
fClusterFunction to get the cluster at index c
nClustersNumber of clusters
B2*kDims x |# clusters| output AABBs

◆ ClustersToAabbs() [2/2]

template<auto kDims, auto kClusterNodes, class FCluster, class TDerivedL, class TDerivedU>
void pbat::geometry::ClustersToAabbs ( FCluster fCluster,
Index nClusters,
Eigen::DenseBase< TDerivedL > & L,
Eigen::DenseBase< TDerivedU > & U )
inline

Computes AABBs of nClusters kDims-dimensional point clusters.

Template Parameters
kDimsSpatial dimension
kClusterNodesNumber of nodes in a cluster
FClusterFunction with signature auto (Index) -> std::convertible_to<Matrix<kDims, kClusterNodes>>
TDerivedLEigen dense expression type
TDerivedUEigen dense expression type
Parameters
fClusterFunction to get the cluster at index c
nClustersNumber of clusters
LkDims x |# clusters| output AABBs lower bounds
UkDims x |# clusters| output AABBs upper bounds

◆ DfsAllNearestNeighbours() [1/2]

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

Note
Although it's possible to obtain nearest neighbour(s) using the KNearestNeighbours function, we also expose this function which focuses on finding the distance minimizer. Both functions are branch and bound tree traversal algorithms. This function exposes 2 variants of the traversal:
  1. Pre-order traversal: This is the default traversal method. It visits the children of a node in the order they are stored in the branch and bound tree, i.e. natural order, or random order, using depth-first search. This has the advantage of using fewer CPU-cycles per node visit, but might require more visits than a guided traversal.
  2. Best-first search: This is an alternative traversal method. At each node, it sorts the children by their lower bounds and visits them in that order. The sort costs non-trivially more stack memory and CPU cycles, but might require fewer visits than a natural order depth-first traversal.
Another approach for tweaking the performance of the NN search is to provide a tight initial upper bound fUpper, which will help prune most of the search space.
Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FDistanceLowerBoundCallable with signature TScalar(TIndex node)
FDistanceCallable with signature TScalar(TIndex o)
FOnFoundCallable with signature void(TIndex o, TScalar d, TIndex k)
NMax number of children per node
TScalarType of the scalar distance
TIndexType of the index
kStackDepthMaximum depth of the traversal's stack
kQueueSizeMaximum size of the nearest neighbour queue
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fLowerFunction to compute the lower bound of the distance to node
fDistanceFunction to compute the distance to object
fOnFoundFunction to call when a nearest neighbour is found
bUseBestFirstSearchUse best-first search instead of pre-order traversal
fUpperUpper bound of the distance to the nearest neighbour
epsEpsilon to consider objects at the same distance
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ DfsAllNearestNeighbours() [2/2]

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

Note
Although it's possible to obtain nearest neighbour(s) using the KNearestNeighbours function, we also expose this function which focuses on finding the distance minimizer. Both functions are branch and bound tree traversal algorithms. This function exposes 2 variants of the traversal:
  1. Pre-order traversal: This is the default traversal method. It visits the children of a node in the order they are stored in the branch and bound tree, i.e. natural order, or random order, using depth-first search. This has the advantage of using fewer CPU-cycles per node visit, but might require more visits than a guided traversal.
  2. Best-first search: This is an alternative traversal method. At each node, it sorts the children by their lower bounds and visits them in that order. The sort costs non-trivially more stack memory and CPU cycles, but might require fewer visits than a natural order depth-first traversal.
Another approach for tweaking the performance of the NN search is to provide a tight initial upper bound fUpper, which will help prune most of the search space.
Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FDistanceLowerBoundCallable with signature TScalar(TIndex node)
FDistanceCallable with signature TScalar(TIndex o)
FOnFoundCallable with signature void(TIndex o, TScalar d, TIndex k)
NMax number of children per node
TScalarType of the scalar distance
TIndexType of the index
kStackDepthMaximum depth of the traversal's stack
kQueueSizeMaximum size of the nearest neighbour queue
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fLowerFunction to compute the lower bound of the distance to node
fDistanceFunction to compute the distance to object
fOnFoundFunction to call when a nearest neighbour is found
bUseBestFirstSearchUse best-first search instead of pre-order traversal
fUpperUpper bound of the distance to the nearest neighbour
epsEpsilon to consider objects at the same distance
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ EdgeEdgeCcd()

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

Template Parameters
TP1TType of the input matrix P1T
TQ1TType of the input matrix Q1T
TP2TType of the input matrix P2T
TQ2TType of the input matrix Q2T
TP1Type of the input matrix P1
TQ1Type of the input matrix Q1
TP2Type of the input matrix P2
TQ2Type of the input matrix Q2
TScalarType of the scalar
Parameters
P1TMatrix of the initial positions of the point P on edge 1
Q1TMatrix of the initial positions of the point Q on edge 1
P2TMatrix of the initial positions of the point P on edge 2
Q2TMatrix of the initial positions of the point Q on edge 2
P1Matrix of the final positions of the point P on edge 1
Q1Matrix of the final positions of the point Q on edge 1
P2Matrix of the final positions of the point P on edge 2
Q2Matrix of the final positions of the point Q on edge 2
Returns
3-vector containing the time of impact and the barycentric coordinates
Postcondition
If no intersection is found, 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

◆ ExpandBits()

PBAT_HOST_DEVICE MortonCodeType pbat::geometry::ExpandBits ( MortonCodeType v)
inline

Expands a 10-bit integer into 30 bits by inserting 2 zeros after each bit.

Note
We make this (otherwise non-templated) function inline so that it gets compiled by nvcc whenever this header is included in cuda sources.
Parameters
v10-bit integer
Returns
Expanded 30-bit integer

◆ HashByXorOfPrimeMultiples()

template<common::CIndex THashed>
auto pbat::geometry::HashByXorOfPrimeMultiples ( )
inline

Returns a hash function used in [13] for 2 or 3-vectors suitable as input to HashGrid's API.

Template Parameters
THashedInteger type for the hash function's output
Returns
Hash function suitable for 2 or 3-vectors

◆ KNearestNeighbours() [1/2]

template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N, class TScalar, class TIndex, auto kHeapCapacity>
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.

Note
This function is a branch and bound tree traversal algorithm, similar to NearestNeighbour. However, we maintain a min-distance-heap of all visited nodes so far, whose space complexities hould scale linearly w.r.t. the height of the tree, i.e. O(log(n)), where n is the number of nodes in the tree. By always visiting nodes closest to the query first, we quickly restrict the upper bound, pruning most of the search space. Each node visit costs approximately O(log(log(n)), due to heap insertions and deletions, but the total number of visits is expected to be much smaller than a generic depth-first traversal. However, the stack memory requirements are higher, due to storing node indices, distances, and types (i.e. node or leaf object).
Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FDistanceLowerBoundCallable with signature TScalar(TIndex node)
FDistanceCallable with signature TScalar(TIndex o)
FOnFoundCallable with signature void(TIndex o, TScalar d, TIndex k)
NMax number of children per node
TScalarType of the scalar distance
TIndexType of the index
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fLowerFunction to compute the lower bound of the distance to node
fDistanceFunction to compute the distance to object
fOnFoundFunction to call when a nearest neighbour is found
KNumber of nearest neighbours to find
fUpperUpper bound of the distance to the nearest neighbour
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ KNearestNeighbours() [2/2]

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

Note
This function is a branch and bound tree traversal algorithm, similar to NearestNeighbour. However, we maintain a min-distance-heap of all visited nodes so far, whose space complexities hould scale linearly w.r.t. the height of the tree, i.e. O(log(n)), where n is the number of nodes in the tree. By always visiting nodes closest to the query first, we quickly restrict the upper bound, pruning most of the search space. Each node visit costs approximately O(log(log(n)), due to heap insertions and deletions, but the total number of visits is expected to be much smaller than a generic depth-first traversal. However, the stack memory requirements are higher, due to storing node indices, distances, and types (i.e. node or leaf object).
Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FDistanceLowerBoundCallable with signature TScalar(TIndex node)
FDistanceCallable with signature TScalar(TIndex o)
FOnFoundCallable with signature void(TIndex o, TScalar d, TIndex k)
NMax number of children per node
TScalarType of the scalar distance
TIndexType of the index
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fLowerFunction to compute the lower bound of the distance to node
fDistanceFunction to compute the distance to object
fOnFoundFunction to call when a nearest neighbour is found
KNumber of nearest neighbours to find
fUpperUpper bound of the distance to the nearest neighbour
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ MeshToAabbs() [1/2]

template<auto kDims, auto kElemNodes, class TDerivedX, class TDerivedE, class TDerivedB>
void pbat::geometry::MeshToAabbs ( Eigen::DenseBase< TDerivedX > const & X,
Eigen::DenseBase< TDerivedE > const & E,
Eigen::DenseBase< TDerivedB > & B )
inline

Computes AABBs of nElemNodes simplex mesh elements in kDims dimensions.

Template Parameters
kDimsSpatial dimension
kElemNodesNumber of nodes in an element
TDerivedXEigen dense expression type
TDerivedEEigen dense expression type
TDerivedBEigen dense expression type
Parameters
XkDims x |# nodes| matrix of node positions
EkElemNodes x |# elements| matrix of element node indices
B2*kDims x |# elements| output AABBs

◆ MeshToAabbs() [2/2]

template<auto kDims, auto kElemNodes, class TDerivedX, class TDerivedE, class TDerivedL, class TDerivedU>
void pbat::geometry::MeshToAabbs ( Eigen::DenseBase< TDerivedX > const & X,
Eigen::DenseBase< TDerivedE > const & E,
Eigen::DenseBase< TDerivedL > & L,
Eigen::DenseBase< TDerivedU > & U )
inline

Computes AABBs of nElemNodes simplex mesh elements in kDims dimensions.

Template Parameters
kDimsSpatial dimension
kElemNodesNumber of nodes in an element
TDerivedXEigen dense expression type
TDerivedEEigen dense expression type
TDerivedLEigen dense expression type for lower bounds
TDerivedUEigen dense expression type for upper bounds
Parameters
XkDims x |# nodes| matrix of node positions
EkElemNodes x |# elements| matrix of element node indices
LkDims x |# elements| output AABBs lower bounds
UkDims x |# elements| output AABBs upper bounds

◆ Morton3D()

template<detail::CMorton3dPoint Point>
PBAT_HOST_DEVICE MortonCodeType pbat::geometry::Morton3D ( Point x)
inline

Calculates a 30-bit Morton code for the given 3D point located within the unit cube [0,1].

Template Parameters
PointType of the point
Parameters
x3D point located within the unit cube [0,1]
Returns
Morton code of x

◆ NearestToFurthestNeighbours() [1/2]

template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FDistanceLowerBound, class FDistance, class FOnFound, auto N, class TScalar, class TIndex, auto kHeapCapacity>
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.

Note
This function is a branch and bound tree traversal algorithm, similar to NearestNeighbour. However, we maintain a min-distance-heap of all visited nodes so far, whose space complexities hould scale linearly w.r.t. the height of the tree, i.e. O(log(n)), where n is the number of nodes in the tree. By always visiting nodes closest to the query first, we quickly restrict the upper bound, pruning most of the search space. Each node visit costs approximately O(log(log(n)), due to heap insertions and deletions, but the total number of visits is expected to be much smaller than a generic depth-first traversal. However, the stack memory requirements are much higher, due to storing node indices, distances, and types (i.e. node or leaf object), and due to the search order not necessarily obeying a depth-first traversal.
Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FDistanceLowerBoundCallable with signature TScalar(TIndex node)
FDistanceCallable with signature TScalar(TIndex o)
FOnFoundCallable with signature bool(TIndex o, TScalar d, TIndex k)
NMax number of children per node
TScalarType of the scalar distance
TIndexType of the index
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fLowerFunction to compute the lower bound of the distance to node
fDistanceFunction to compute the distance to object
fOnFoundFunction to call when a nearest neighbour is found, returns true if the search should continue
fUpperUpper bound of the distance to the nearest neighbour
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ NearestToFurthestNeighbours() [2/2]

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

Note
This function is a branch and bound tree traversal algorithm, similar to NearestNeighbour. However, we maintain a min-distance-heap of all visited nodes so far, whose space complexities hould scale linearly w.r.t. the height of the tree, i.e. O(log(n)), where n is the number of nodes in the tree. By always visiting nodes closest to the query first, we quickly restrict the upper bound, pruning most of the search space. Each node visit costs approximately O(log(log(n)), due to heap insertions and deletions, but the total number of visits is expected to be much smaller than a generic depth-first traversal. However, the stack memory requirements are much higher, due to storing node indices, distances, and types (i.e. node or leaf object), and due to the search order not necessarily obeying a depth-first traversal.
Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FDistanceLowerBoundCallable with signature TScalar(TIndex node)
FDistanceCallable with signature TScalar(TIndex o)
FOnFoundCallable with signature bool(TIndex o, TScalar d, TIndex k)
NMax number of children per node
TScalarType of the scalar distance
TIndexType of the index
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fLowerFunction to compute the lower bound of the distance to node
fDistanceFunction to compute the distance to object
fOnFoundFunction to call when a nearest neighbour is found, returns true if the search should continue
fUpperUpper bound of the distance to the nearest neighbour
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ Overlaps() [1/3]

template<class FChild, class FIsLeaf, class FLeafSize, class FLeafObject, class FNodeOverlap, class FObjectOverlap, class FOnFound, auto N, class TIndex, auto kStackDepth>
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.

Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FNodeOverlapsCallable with signature bool(TIndex node)
FObjectOverlapsCallable with signature bool(TIndex o)
FOnFoundCallable with signature void(TIndex o, TIndex k)
NMax number of children per node
TIndexType of the index
kStackDepthMaximum depth of the traversal's stack
kQueueSizeMaximum size of the nearest neighbour queue
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fNodeOverlapsFunction to determine if a node overlaps with the query
fObjectOverlapsFunction to determine if an object overlaps with the query
fOnFoundFunction to call when a nearest neighbour is found
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ Overlaps() [2/3]

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

Template Parameters
FChildCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FNodeOverlapsCallable with signature bool(TIndex node)
FObjectOverlapsCallable with signature bool(TIndex o)
FOnFoundCallable with signature void(TIndex o, TIndex k)
NMax number of children per node
TIndexType of the index
kStackDepthMaximum depth of the traversal's stack
kQueueSizeMaximum size of the nearest neighbour queue
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fNodeOverlapsFunction to determine if a node overlaps with the query
fObjectOverlapsFunction to determine if an object overlaps with the query
fOnFoundFunction to call when a nearest neighbour is found
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ Overlaps() [3/3]

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

Template Parameters
FChildLhsCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafLhsCallable with signature bool(TIndex node)
FLeafSizeLhsCallable with signature TIndex(TIndex node)
FLeafObjectLhsCallable with signature TIndex(TIndex node, TIndex i)
FChildRhsCallable with signature template <auto c> TIndex(TIndex node)
FIsLeafRhsCallable with signature bool(TIndex node)
FLeafSizeRhsCallable with signature TIndex(TIndex node)
FLeafObjectRhsCallable with signature TIndex(TIndex node, TIndex i)
FNodesOverlapCallable with signature bool(TIndex nlhs, TIndex nrhs)
FObjectsOverlapCallable with signature bool(TIndex olhs, TIndex orhs)
FOnFoundCallable with signature void(TIndex olhs, TIndex orhs, TIndex k)
NLhsNumber of children per node in the left tree
NRhsNumber of children per node in the right tree
TIndexType of the index
kStackDepthMaximum depth of the traversal's stack
Parameters
fChildLhsFunction to get child c of a node in the left tree. Returns the child index or -1 if no child.
fIsLeafLhsFunction to determine if a node is a leaf node in the left tree
fLeafSizeLhsFunction to get the number of leaf objects in a node in the left tree
fLeafObjectLhsFunction to get the i-th leaf object in a node in the left tree
fChildRhsFunction to get child c of a node in the right tree. Returns the child index or -1 if no child.
fIsLeafRhsFunction to determine if a node is a leaf node in the right tree
fLeafSizeRhsFunction to get the number of leaf objects in a node in the right tree
fLeafObjectRhsFunction to get the i-th leaf object in a node in the right tree
fNodesOverlapFunction to determine if two nodes overlap
fObjectsOverlapFunction to determine if two objects overlap
fOnFoundFunction to call when an overlap is found
rootLhsIndex of the root node in the left tree to start the search from
rootRhsIndex of the root node in the right tree to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ PointTriangleCcd()

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

Template Parameters
TXTType of the input matrix XT
TATType of the input matrix AT
TBTType of the input matrix BT
TCTType of the input matrix CT
TXType of the input matrix X
TAType of the input matrix A
TBType of the input matrix B
TCType of the input matrix C
TScalarType of the scalar
Parameters
XTMatrix of the initial positions of the point
ATMatrix of the initial positions of the triangle vertex A
BTMatrix of the initial positions of the triangle vertex B
CTMatrix of the initial positions of the triangle vertex C
XMatrix of the final positions of the point
AMatrix of the final positions of the triangle vertex A
BMatrix of the final positions of the triangle vertex B
CMatrix of the final positions of the triangle vertex C
Returns
A 4-vector r containing the time of impact and the barycentric coordinates.
Postcondition
If no intersection is found, r[0] < 0, otherwise r[0] >= 0 is the earliest time of impact (subject to floating point error).

◆ SelfOverlaps()

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

Template Parameters
FChildCallable with signature template <auto c> Index(TIndex node)
FIsLeafCallable with signature bool(TIndex node)
FLeafSizeCallable with signature TIndex(TIndex node)
FLeafObjectCallable with signature TIndex(TIndex node, TIndex i)
FNodesOverlapCallable with signature bool(TIndex n1, TIndex n2)
FObjectsOverlapCallable with signature bool(TIndex o1, TIndex o2)
FOnFoundCallable with signature void(TIndex o1, TIndex o2, TIndex k)
NMax number of children per node
TIndexType of the index
kStackDepthMaximum depth of the traversal's stack
kQueueSizeMaximum size of the nearest neighbour queue
Parameters
fChildFunction to get child c of a node. Returns the child index or -1 if no child.
fIsLeafFunction to determine if a node is a leaf node
fLeafSizeFunction to get the number of leaf objects in a node
fLeafObjectFunction to get the i-th leaf object in a node
fNodesOverlapFunction to determine if 2 nodes overlap
fObjectsOverlapFunction to determine if 2 objects overlap
fOnFoundFunction to call when a nearest neighbour is found
rootIndex of the root node to start the search from
Returns
true if the traversal completed, false if it was stopped early due to insufficient stack capacity

◆ SimplexMeshBoundary()

template<common::CIndex TIndex = Index>
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.

Note
Only works for triangle (C.rows()==3) and tetrahedral (C.rows()==4) meshes.
Template Parameters
TIndexThe index type used in the mesh (default: Index)
Parameters
CThe connectivity matrix of the mesh (i.e. the simplices)
nThe number of vertices in the mesh. If -1, the number of vertices is computed from C.
Returns
A tuple containing the boundary vertices and the boundary facets