forked from arborx/ArborX
-
Notifications
You must be signed in to change notification settings - Fork 0
ArborX::BoundingVolumeHierarchy::BoundingVolumeHierarchy
Andrey Prokopenko edited this page Oct 28, 2024
·
7 revisions
ArborX
/ Spatial indexes / ArborX::BVH
BVH() noexcept; // (1)
template <typename ExecutionSpace, typename Values>
BVH(ExecutionSpace const& space, Values const& values); // (2)
template <typename ExecutionSpace, typename Values>
BVH(ExecutionSpace const& space, Values const& values, IndexableGetter const& indexable_getter); // (3)
- Default constructor. Constructs an empty tree.
- Constructs a bounding volume hierarchy from a given data source.
- Constructs a bounding volume hierarchy from a given data source and a function tranforming a value into indexable.
space
: execution space
values
: values stored in the container
indexable_getter
: a function object to transform a value into an indexable
-
MemorySpace
must be accessible fromExecutionSpace
(i.e.,Kokkos::SpaceAccessibility<ExecutionSpace, MemorySpace>::accessible
must betrue
). - A specialization of
ArborX::AccessTraits
must matchValues
as the template argument. - The return type of
ArborX::AccessTraits<Values, ArborX::PrimitivesTag>::get()
must decay toValue
. ArborX provides specializations for Kokkos views but a user may specialize it for their types. -
IndexableGetter
must provide anoperator()
taking in a type that decays toValue
and returning an indexable (geometric object).
O(N log N), where N is the number of primitives passed to the constructor (ArborX::AccessTraits<Values, ArborX::PrimitivesTag>::size(values)
).
Memory allocation with Kokkos may throw.
#include <ArborX_LinearBVH.hpp>
#include <Kokkos_Core.hpp>
#include <iostream>
std::ostream &operator<<(std::ostream &os, ArborX::Point const &p)
{
os << '(' << p[0] << ',' << p[1] << ',' << p[2] << ')';
return os;
}
int main(int argc, char *argv[])
{
Kokkos::ScopeGuard guard(argc, argv);
Kokkos::View<ArborX::Point<3> *> cloud("point_cloud", 1000);
Kokkos::parallel_for(1000, KOKKOS_LAMBDA(int i) {
cloud[i] = {(float)i, (float)i, (float)i};
});
using memory_space = decltype(cloud)::memory_space; // where to store the tree
using execution_space = decltype(cloud)::execution_space; // where to execute code
ArborX::BVH bvh{execution_space{}, cloud};
auto const box = bvh.bounds();
std::cout << box.minCorner() << " - " << box.maxCorner() << '\n';
return 0;
}
Output
(0,0,0) - (999,999,999)
query
: search for all primitives that meet some predicates.
bounds
: returns the bounding volume that contains all leaves.