--- /dev/null
+CMAKE_MINIMUM_REQUIRED(VERSION 2.8.9)
+INCLUDE(../setup_testsubproject.cmake)
+PROJECT(testsuite CXX)
+IF(DEAL_II_WITH_ARBORX)
+ DEAL_II_PICKUP_TESTS()
+ENDIF()
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Check ArborX wrapper: intersection of boundary boxes
+
+
+#include <deal.II/arborx/access_traits.h>
+#include <deal.II/arborx/bvh.h>
+
+#include <deal.II/base/bounding_box.h>
+
+#include "../tests.h"
+
+
+void
+test_2d()
+{
+ std::vector<BoundingBox<2>> bounding_boxes;
+ std::vector<Point<2>> points;
+
+ unsigned int n_points_1d = 5;
+ for (unsigned int i = 0; i < n_points_1d; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d; ++j)
+ {
+ points.emplace_back(j, i);
+ }
+ }
+
+ for (unsigned int i = 0; i < n_points_1d - 1; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d - 1; ++j)
+ {
+ unsigned int point_index = j + i * n_points_1d;
+ bounding_boxes.push_back(
+ std::make_pair(points[point_index],
+ points[point_index + n_points_1d + 1]));
+ }
+ }
+
+ Point<2> point_a(0.5, 0.5);
+ Point<2> point_b(1.5, 1.5);
+ Point<2> point_c(2.2, 2.2);
+ Point<2> point_d(2.6, 2.6);
+
+ std::vector<BoundingBox<2>> query_bounding_boxes;
+ query_bounding_boxes.emplace_back(std::make_pair(point_a, point_b));
+ query_bounding_boxes.emplace_back(std::make_pair(point_c, point_d));
+
+ ArborXWrappers::BVH bvh(bounding_boxes);
+ ArborXWrappers::BoundingBoxIntersect bb_intersect(query_bounding_boxes);
+ auto indices_offset = bvh.query(bb_intersect);
+ std::vector<int> indices = indices_offset.first;
+ std::vector<int> offset = indices_offset.second;
+
+ std::vector<int> indices_ref = {0, 1, 4, 5, 10};
+ std::vector<int> offset_ref = {0, 4, 5};
+
+ AssertThrow(indices.size() == indices_ref.size(), ExcInternalError());
+ AssertThrow(offset.size() == offset_ref.size(), ExcInternalError());
+ for (unsigned int i = 0; i < offset.size() - 1; ++i)
+ {
+ for (unsigned int j = offset[i]; j < offset[i + 1]; ++j)
+ {
+ // The indices associated to each query are not ordered.
+ bool found = false;
+ for (unsigned int k = offset[i]; k < offset[i + 1]; ++k)
+ {
+ if (indices[j] == indices_ref[k])
+ {
+ found = true;
+ break;
+ }
+ }
+ AssertThrow(found, ExcInternalError());
+ }
+ }
+ for (unsigned int i = 0; i < offset.size(); ++i)
+ AssertThrow(offset[i] == offset_ref[i], ExcInternalError());
+
+ deallog << "OK" << std::endl;
+}
+
+void
+test_3d()
+{
+ std::vector<BoundingBox<3>> bounding_boxes;
+ std::vector<Point<3>> points;
+
+ unsigned int n_points_1d = 5;
+ for (unsigned int i = 0; i < n_points_1d; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d; ++j)
+ {
+ for (unsigned int k = 0; k < n_points_1d; ++k)
+ {
+ points.emplace_back(k, j, i);
+ }
+ }
+ }
+
+ for (unsigned int i = 0; i < n_points_1d - 1; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d - 1; ++j)
+ {
+ for (unsigned int k = 0; k < n_points_1d - 1; ++k)
+ {
+ unsigned int point_index =
+ k + j * n_points_1d + i * n_points_1d * n_points_1d;
+ bounding_boxes.push_back(
+ std::make_pair(points[point_index],
+ points[point_index + n_points_1d * n_points_1d +
+ n_points_1d + 1]));
+ }
+ }
+ }
+
+ Point<3> point_a(0.5, 0.5, 0.5);
+ Point<3> point_b(1.5, 1.5, 1.5);
+ Point<3> point_c(2.2, 2.2, 2.2);
+ Point<3> point_d(2.6, 2.6, 2.6);
+
+ std::vector<BoundingBox<3>> query_bounding_boxes;
+ query_bounding_boxes.emplace_back(std::make_pair(point_a, point_b));
+ query_bounding_boxes.emplace_back(std::make_pair(point_c, point_d));
+
+ ArborXWrappers::BVH bvh(bounding_boxes);
+ ArborXWrappers::BoundingBoxIntersect bb_intersect(query_bounding_boxes);
+ auto indices_offset = bvh.query(bb_intersect);
+ std::vector<int> indices = indices_offset.first;
+ std::vector<int> offset = indices_offset.second;
+
+ std::vector<int> indices_ref = {0, 1, 4, 5, 16, 17, 20, 21, 42};
+ std::vector<int> offset_ref = {0, 8, 9};
+
+ AssertThrow(indices.size() == indices_ref.size(), ExcInternalError());
+ AssertThrow(offset.size() == offset_ref.size(), ExcInternalError());
+ for (unsigned int i = 0; i < offset.size() - 1; ++i)
+ {
+ for (unsigned int j = offset[i]; j < offset[i + 1]; ++j)
+ {
+ // The indices associated to each query are not ordered.
+ bool found = false;
+ for (unsigned int k = offset[i]; k < offset[i + 1]; ++k)
+ {
+ if (indices[j] == indices_ref[k])
+ {
+ found = true;
+ break;
+ }
+ }
+ AssertThrow(found, ExcInternalError());
+ }
+ }
+ for (unsigned int i = 0; i < offset.size(); ++i)
+ AssertThrow(offset[i] == offset_ref[i], ExcInternalError());
+ deallog << "OK" << std::endl;
+}
+
+int
+main(int argc, char **argv)
+{
+ // Initialize Kokkos
+ Kokkos::initialize(argc, argv);
+
+ initlog();
+
+ // tests
+ test_2d();
+ test_3d();
+
+ Kokkos::finalize();
+}
--- /dev/null
+
+DEAL::OK
+DEAL::OK
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+
+// Check ArborX wrapper: intersection of points with bounding boxes
+
+
+#include <deal.II/arborx/access_traits.h>
+#include <deal.II/arborx/bvh.h>
+
+#include <deal.II/base/point.h>
+
+#include "../tests.h"
+
+
+void
+test_2d()
+{
+ std::vector<BoundingBox<2>> bounding_boxes;
+ std::vector<Point<2>> points;
+
+ unsigned int n_points_1d = 5;
+ for (unsigned int i = 0; i < n_points_1d; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d; ++j)
+ {
+ points.emplace_back(j, i);
+ }
+ }
+
+ for (unsigned int i = 0; i < n_points_1d - 1; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d - 1; ++j)
+ {
+ unsigned int point_index = j + i * n_points_1d;
+ bounding_boxes.push_back(
+ std::make_pair(points[point_index],
+ points[point_index + n_points_1d + 1]));
+ }
+ }
+
+ std::vector<Point<2>> query_points;
+ query_points.emplace_back(0.5, 0.5);
+ query_points.emplace_back(1.5, 1.5);
+ query_points.emplace_back(2.2, 2.2);
+ query_points.emplace_back(2.6, 2.6);
+
+
+ ArborXWrappers::BVH bvh(bounding_boxes);
+ ArborXWrappers::PointIntersect pt_intersect(query_points);
+ auto indices_offset = bvh.query(pt_intersect);
+ std::vector<int> indices = indices_offset.first;
+ std::vector<int> offset = indices_offset.second;
+
+ std::vector<int> indices_ref = {0, 5, 10, 10};
+ std::vector<int> offset_ref = {0, 1, 2, 3, 4};
+
+ AssertThrow(indices.size() == indices_ref.size(), ExcInternalError());
+ AssertThrow(offset.size() == offset_ref.size(), ExcInternalError());
+ for (unsigned int i = 0; i < offset.size() - 1; ++i)
+ {
+ for (unsigned int j = offset[i]; j < offset[i + 1]; ++j)
+ {
+ // The indices associated to each query are not ordered.
+ bool found = false;
+ for (unsigned int k = offset[i]; k < offset[i + 1]; ++k)
+ {
+ if (indices[j] == indices_ref[k])
+ {
+ found = true;
+ break;
+ }
+ }
+ AssertThrow(found, ExcInternalError());
+ }
+ }
+ for (unsigned int i = 0; i < offset.size(); ++i)
+ AssertThrow(offset[i] == offset_ref[i], ExcInternalError());
+
+ deallog << "OK" << std::endl;
+}
+
+void
+test_3d()
+{
+ std::vector<BoundingBox<3>> bounding_boxes;
+ std::vector<Point<3>> points;
+
+ unsigned int n_points_1d = 5;
+ for (unsigned int i = 0; i < n_points_1d; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d; ++j)
+ {
+ for (unsigned int k = 0; k < n_points_1d; ++k)
+ {
+ points.emplace_back(k, j, i);
+ }
+ }
+ }
+
+ for (unsigned int i = 0; i < n_points_1d - 1; ++i)
+ {
+ for (unsigned int j = 0; j < n_points_1d - 1; ++j)
+ {
+ for (unsigned int k = 0; k < n_points_1d - 1; ++k)
+ {
+ unsigned int point_index =
+ k + j * n_points_1d + i * n_points_1d * n_points_1d;
+ bounding_boxes.push_back(
+ std::make_pair(points[point_index],
+ points[point_index + n_points_1d * n_points_1d +
+ n_points_1d + 1]));
+ }
+ }
+ }
+
+ std::vector<Point<3>> query_points;
+ query_points.emplace_back(0.5, 0.5, 0.5);
+ query_points.emplace_back(1.5, 1.5, 1.5);
+ query_points.emplace_back(2.2, 2.2, 2.2);
+ query_points.emplace_back(2.6, 2.6, 2.6);
+
+
+ ArborXWrappers::BVH bvh(bounding_boxes);
+ ArborXWrappers::PointIntersect pt_intersect(query_points);
+ auto indices_offset = bvh.query(pt_intersect);
+ std::vector<int> indices = indices_offset.first;
+ std::vector<int> offset = indices_offset.second;
+
+
+ std::vector<int> indices_ref = {0, 21, 42, 42};
+ std::vector<int> offset_ref = {0, 1, 2, 3, 4};
+
+ AssertThrow(indices.size() == indices_ref.size(), ExcInternalError());
+ AssertThrow(offset.size() == offset_ref.size(), ExcInternalError());
+ for (unsigned int i = 0; i < offset.size() - 1; ++i)
+ {
+ for (unsigned int j = offset[i]; j < offset[i + 1]; ++j)
+ {
+ // The indices associated to each query are not ordered.
+ bool found = false;
+ for (unsigned int k = offset[i]; k < offset[i + 1]; ++k)
+ {
+ if (indices[j] == indices_ref[k])
+ {
+ found = true;
+ break;
+ }
+ }
+ AssertThrow(found, ExcInternalError());
+ }
+ }
+ for (unsigned int i = 0; i < offset.size(); ++i)
+ AssertThrow(offset[i] == offset_ref[i], ExcInternalError());
+ deallog << "OK" << std::endl;
+}
+
+int
+main(int argc, char **argv)
+{
+ // Initialize Kokkos
+ Kokkos::initialize(argc, argv);
+
+ initlog();
+
+ // Initialize ArborX
+ Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv);
+
+ // tests
+ test_2d();
+ test_3d();
+
+ Kokkos::finalize();
+}
--- /dev/null
+
+DEAL::OK
+DEAL::OK
FOREACH(_var ${_variables})
IF( _var MATCHES "^(TEST|DEAL_II|ALLOW|WITH|FORCE|COMPONENT)_" OR
_var MATCHES "^(DOCUMENTATION|EXAMPLES)" OR
- _var MATCHES "^(ADOLC|ARPACK|BOOST|OPENCASCADE|MUPARSER|HDF5|METIS|MPI)_" OR
+ _var MATCHES "^(ADOLC|ARBORX|ARPACK|BOOST|OPENCASCADE|MUPARSER|HDF5|KOKKOS|METIS|MPI)_" OR
_var MATCHES "^(GINKGO|P4EST|PETSC|SCALAPACK|SLEPC|THREADS|TBB|TRILINOS)_" OR
_var MATCHES "^(UMFPACK|ZLIB|LAPACK|MUPARSER|CUDA)_" OR
_var MATCHES "^(CMAKE|DEAL_II)_(C|CXX|Fortran|BUILD)_(COMPILER|FLAGS)" OR