Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Use 64-bit Morton indices by default in the BVH construction #637

Merged
merged 3 commits into from
Mar 3, 2022
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 1 addition & 1 deletion src/ArborX_LinearBVH.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -241,7 +241,7 @@ BasicBoundingVolumeHierarchy<MemorySpace, BoundingVolume, Enable>::
Kokkos::Profiling::pushRegion("ArborX::BVH::BVH::assign_morton_codes");

// calculate Morton codes of all objects
Kokkos::View<unsigned int *, MemorySpace> morton_indices(
Kokkos::View<unsigned long long *, MemorySpace> morton_indices(
Kokkos::view_alloc(space, Kokkos::WithoutInitializing,
"ArborX::BVH::BVH::morton"),
size());
Expand Down
2 changes: 2 additions & 0 deletions src/details/ArborX_DetailsBatchedQueries.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -66,6 +66,8 @@ struct BatchedQueries
using Details::returnCentroid;
Point xyz = returnCentroid(getGeometry(Access::get(predicates, i)));
translateAndScale(xyz, xyz, scene_bounding_box);
// Use 32-bit Morton indices instead of 64-bit as in construction. For
// most (all?) situations, 64-bit just adds a penalty with no benefit.
morton_codes(i) = morton32(xyz[0], xyz[1], xyz[2]);
});

Expand Down
38 changes: 23 additions & 15 deletions src/details/ArborX_DetailsTreeConstruction.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -69,7 +69,7 @@ assignMortonCodesImpl(ExecutionSpace const &space, Primitives const &primitives,
Point xyz;
centroid(Access::get(primitives, i), xyz);
translateAndScale(xyz, xyz, scene_bounding_box);
morton_codes(i) = morton32(xyz[0], xyz[1], xyz[2]);
morton_codes(i) = morton64(xyz[0], xyz[1], xyz[2]);
});
}

Expand All @@ -87,15 +87,16 @@ assignMortonCodesImpl(ExecutionSpace const &space, Primitives const &primitives,
Kokkos::RangePolicy<ExecutionSpace>(space, 0, n), KOKKOS_LAMBDA(int i) {
Point xyz;
translateAndScale(Access::get(primitives, i), xyz, scene_bounding_box);
morton_codes(i) = morton32(xyz[0], xyz[1], xyz[2]);
morton_codes(i) = morton64(xyz[0], xyz[1], xyz[2]);
});
}

template <typename ExecutionSpace, typename Primitives,
typename... MortonCodesViewProperties>
inline void assignMortonCodes(
ExecutionSpace const &space, Primitives const &primitives,
Kokkos::View<unsigned int *, MortonCodesViewProperties...> morton_codes,
Kokkos::View<unsigned long long *, MortonCodesViewProperties...>
morton_codes,
Box const &scene_bounding_box)
{
using Access = AccessTraits<Primitives, PrimitivesTag>;
Expand Down Expand Up @@ -152,7 +153,7 @@ class GenerateHierarchy
ExecutionSpace const &space, Primitives const &primitives,
Kokkos::View<unsigned int const *, PermutationIndicesViewProperties...>
permutation_indices,
Kokkos::View<unsigned int const *, MortonCodesViewProperties...>
Kokkos::View<unsigned long long const *, MortonCodesViewProperties...>
sorted_morton_codes,
Kokkos::View<Node *, LeafNodesViewProperties...> leaf_nodes,
Kokkos::View<Node *, InternalNodesViewProperties...> internal_nodes)
Expand All @@ -176,7 +177,7 @@ class GenerateHierarchy
}

KOKKOS_FUNCTION
int delta(int const i) const
long long delta(int const i) const
{
// Per Apetrei:
// Because we already know where the highest differing bit is for each
Expand All @@ -191,7 +192,7 @@ class GenerateHierarchy
// This check is here simply to avoid code complications in the main
// operator
if (i < 0 || i >= _num_internal_nodes)
return INT_MAX;
return LLONG_MAX;

// The Apetrei's paper does not mention dealing with duplicate indices. We
// follow the original Karras idea in this situation:
Expand All @@ -202,10 +203,17 @@ class GenerateHierarchy
// concatenation.
// In this case, if the Morton indices are the same, we want to compare is.
// We also want the result in this situation to always be less than any
// Morton comparison. Thus, we add INT_MIN to it.
// We also avoid if/else statement by doing a "x + !x*<blah>" trick.
// Morton comparison. Thus, we add LLONG_MIN to it.
auto x = _sorted_morton_codes(i) ^ _sorted_morton_codes(i + 1);
return x + (!x) * (INT_MIN + (i ^ (i + 1)));
if (x != 0)
{
// When using 63 bits for Morton codes, the LLONG_MAX is actually a valid
// code. As we want the return statement above to return a value always
// greater than anything here, we downshift by 1.
return x - 1;
}

return LLONG_MIN + (i ^ (i + 1));
aprokop marked this conversation as resolved.
Show resolved Hide resolved
}

KOKKOS_FUNCTION Node *getNodePtr(int i) const
Expand Down Expand Up @@ -240,7 +248,7 @@ class GenerateHierarchy
template <typename Tag = typename Node::Tag>
KOKKOS_FUNCTION
std::enable_if_t<std::is_same<Tag, NodeWithLeftChildAndRopeTag>{}>
setRope(Node *node, int range_right, int delta_right) const
setRope(Node *node, int range_right, long long delta_right) const
{
int rope;
if (range_right != _num_internal_nodes)
Expand Down Expand Up @@ -284,8 +292,8 @@ class GenerateHierarchy
int range_left = i - leaf_nodes_shift;
int range_right = range_left;

int delta_left = delta(range_left - 1);
int delta_right = delta(range_right);
auto delta_left = delta(range_left - 1);
auto delta_right = delta(range_right);

setRope(leaf_node, range_right, delta_right);

Expand Down Expand Up @@ -381,7 +389,7 @@ class GenerateHierarchy
private:
Primitives _primitives;
Kokkos::View<unsigned int const *, MemorySpace> _permutation_indices;
Kokkos::View<unsigned int const *, MemorySpace> _sorted_morton_codes;
Kokkos::View<unsigned long long const *, MemorySpace> _sorted_morton_codes;
Kokkos::View<Node *, MemorySpace> _leaf_nodes;
Kokkos::View<Node *, MemorySpace> _internal_nodes;
Kokkos::View<int *, MemorySpace> _ranges;
Expand All @@ -397,15 +405,15 @@ void generateHierarchy(
ExecutionSpace const &space, Primitives const &primitives,
Kokkos::View<unsigned int *, PermutationIndicesViewProperties...>
permutation_indices,
Kokkos::View<unsigned int *, MortonCodesViewProperties...>
Kokkos::View<unsigned long long *, MortonCodesViewProperties...>
sorted_morton_codes,
Kokkos::View<Node *, LeafNodesViewProperties...> leaf_nodes,
Kokkos::View<Node *, InternalNodesViewProperties...> internal_nodes)
{
using ConstPermutationIndices =
Kokkos::View<unsigned int const *, PermutationIndicesViewProperties...>;
using ConstMortonCodes =
Kokkos::View<unsigned int const *, MortonCodesViewProperties...>;
Kokkos::View<unsigned long long const *, MortonCodesViewProperties...>;

using MemorySpace = typename decltype(internal_nodes)::memory_space;

Expand Down
48 changes: 30 additions & 18 deletions test/tstDetailsTreeConstruction.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -30,24 +30,35 @@ namespace tt = boost::test_tools;
BOOST_AUTO_TEST_CASE_TEMPLATE(assign_morton_codes, DeviceType,
ARBORX_DEVICE_TYPES)
{
std::vector<ArborX::Point> points = {
{{0.0, 0.0, 0.0}}, {{0.25, 0.75, 0.25}}, {{0.75, 0.25, 0.25}},
{{0.75, 0.75, 0.25}}, {{1.33, 2.33, 3.33}}, {{1.66, 2.66, 3.66}},
{{1024.0, 1024.0, 1024.0}},
};
// N is the number of Morton grid cells in each dimension for 64-bit Morton
// codes.
constexpr unsigned long long N = 1 << 21;
std::vector<ArborX::Point> points = {{{0.0, 0.0, 0.0}},
{{0.25, 0.75, 0.25}},
{{0.75, 0.25, 0.25}},
{{0.75, 0.75, 0.25}},
{{1.33, 2.33, 3.33}},
{{1.66, 2.66, 3.66}},
{{(float)N, (float)N, (float)N}}};
int const n = points.size();
// lower left front corner corner of the octant the points fall in
std::vector<std::array<unsigned int, 3>> anchors = {
{{0, 0, 0}}, {{0, 0, 0}}, {{0, 0, 0}}, {{0, 0, 0}},
{{1, 2, 3}}, {{1, 2, 3}}, {{1023, 1023, 1023}}};
auto fun = [](std::array<unsigned int, 3> const &anchor) {
std::vector<std::array<unsigned long long, 3>> anchors = {
{{0, 0, 0}},
{{0, 0, 0}},
{{0, 0, 0}},
{{0, 0, 0}},
{{1, 2, 3}},
{{1, 2, 3}},
{{N - 1, N - 1, N - 1}}};
auto fun = [](std::array<unsigned long long, 3> const &anchor) {
using ArborX::Details::expandBitsBy2;
unsigned int i = std::get<0>(anchor);
unsigned int j = std::get<1>(anchor);
unsigned int k = std::get<2>(anchor);
auto i = std::get<0>(anchor);
auto j = std::get<1>(anchor);
auto k = std::get<2>(anchor);
return 4 * expandBitsBy2(i) + 2 * expandBitsBy2(j) + expandBitsBy2(k);
};
std::vector<unsigned int> ref(n, std::numeric_limits<unsigned int>::max());
std::vector<unsigned long long> ref(
n, std::numeric_limits<unsigned long long>::max());
for (int i = 0; i < n; ++i)
ref[i] = fun(anchors[i]);
// using points rather than boxes for convenience here but still have to
Expand All @@ -64,9 +75,10 @@ BOOST_AUTO_TEST_CASE_TEMPLATE(assign_morton_codes, DeviceType,
space, boxes, scene_host);

BOOST_TEST(ArborX::Details::equals(
scene_host, {{{0., 0., 0.}}, {{1024., 1024., 1024.}}}));
scene_host, {{{0.0, 0.0, 0.0}}, {{(float)N, (float)N, (float)N}}}));

Kokkos::View<unsigned int *, DeviceType> morton_codes("morton_codes", n);
Kokkos::View<unsigned long long *, DeviceType> morton_codes("morton_codes",
n);
ArborX::Details::TreeConstruction::assignMortonCodes(
space, boxes, morton_codes, scene_host);
auto morton_codes_host = Kokkos::create_mirror_view(morton_codes);
Expand Down Expand Up @@ -231,16 +243,16 @@ BOOST_AUTO_TEST_CASE_TEMPLATE(example_tree_construction, DeviceType,
// See
// https://devblogs.nvidia.com/parallelforall/thinking-parallel-part-iii-tree-construction-gpu/
int const n = 8;
Kokkos::View<unsigned int *, DeviceType> sorted_morton_codes(
Kokkos::View<unsigned long long *, DeviceType> sorted_morton_codes(
"sorted_morton_codes", n);
std::vector<std::string> s{
"00001", "00010", "00100", "00101", "10011", "11000", "11001", "11110",
};
for (int i = 0; i < n; ++i)
{
std::bitset<6> b(s[i]);
BOOST_TEST_MESSAGE(b << " " << b.to_ulong());
sorted_morton_codes(i) = b.to_ulong();
BOOST_TEST_MESSAGE(b << " " << b.to_ullong());
sorted_morton_codes(i) = b.to_ullong();
}

Kokkos::View<Test::FakePrimitive *, DeviceType> primitives(
Expand Down