Update branch from master after trailing whitespaces and tabs removal

This commit is contained in:
Sébastien Loriot 2020-03-26 19:28:07 +01:00
commit c67d093ebf
7 changed files with 799 additions and 532 deletions

View File

@ -5,27 +5,20 @@ cmake_minimum_required(VERSION 3.1...3.15)
project( AABB_traits_benchmark)
# CGAL and its components
find_package( CGAL QUIET)
if ( CGAL_FOUND )
include( ${CGAL_USE_FILE} )
else ()
message(STATUS "This project requires the CGAL library, and will not be compiled.")
return()
endif()
find_package( CGAL REQUIRED )
include( ${CGAL_USE_FILE} )
# Boost and its components
find_package( Boost REQUIRED )
# include for local directory
if ( NOT Boost_FOUND )
message(STATUS "This project requires the Boost library, and will not be compiled.")
return()
endif()
# google benchmark
find_package(benchmark REQUIRED)
# include for local package
include_directories( BEFORE ../../include )
add_executable (test_ test.cpp)
# add_executable (test_ test.cpp) // TODO: fix this benchmark
add_executable(tree_creation tree_creation.cpp)
target_link_libraries(tree_creation benchmark::benchmark)

View File

@ -0,0 +1,69 @@
#include <benchmark/benchmark.h>
#include <CGAL/AABB_tree.h>
#include <CGAL/AABB_traits.h>
#include <CGAL/AABB_face_graph_triangle_primitive.h>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Surface_mesh.h>
#include <fstream>
typedef CGAL::Simple_cartesian<double> K;
typedef CGAL::Surface_mesh<K::Point_3> Surface_mesh;
typedef CGAL::AABB_face_graph_triangle_primitive<Surface_mesh> Primitive;
typedef CGAL::AABB_traits<K, Primitive> Traits;
typedef CGAL::AABB_tree<Traits> Tree;
typedef K::Segment_3 Segment;
typedef K::Point_3 Point_3;
Surface_mesh mesh;
static void BM_TreeCreation(benchmark::State& state)
{
for (auto _ : state)
{
benchmark::DoNotOptimize([]() {
Tree tree{mesh.faces_begin(), mesh.faces_end(), mesh};
tree.build();
return 0;
}());
}
}
BENCHMARK(BM_TreeCreation);
static void BM_Intersections(benchmark::State& state)
{
Point_3 p(-0.5, 0.03, 0.04);
Point_3 q(-0.5, 0.04, 0.06);
Tree tree{mesh.faces_begin(), mesh.faces_end(), mesh};
tree.accelerate_distance_queries();
Segment segment_query(p, q);
for (auto _ : state)
{
benchmark::DoNotOptimize([&]() {
tree.number_of_intersected_primitives(segment_query);
Point_3 point_query(2.0, 2.0, 2.0);
Point_3 closest = tree.closest_point(point_query);
return 0;
}());
}
}
BENCHMARK(BM_Intersections);
int main(int argc, char** argv)
{
const char* default_file = "data/handle.off";
const char* filename = argc > 2? argv[2] : default_file;
{
std::ifstream input(filename);
input >> mesh;
}
::benchmark::Initialize(&argc, argv);
::benchmark::RunSpecifiedBenchmarks();
return EXIT_SUCCESS;
}

File diff suppressed because it is too large Load Diff

View File

@ -32,6 +32,9 @@ namespace CGAL {
template<typename AABBTraits>
class AABB_node
{
private:
typedef AABB_node<AABBTraits> Self;
public:
typedef typename AABBTraits::Bounding_box Bounding_box;
@ -41,27 +44,15 @@ public:
, m_p_left_child(nullptr)
, m_p_right_child(nullptr) { };
/// Non virtual Destructor
/// Do not delete children because the tree hosts and delete them
~AABB_node() { };
AABB_node(Self&& node) = default;
// Disabled copy constructor & assignment operator
AABB_node(const Self& src) = delete;
Self& operator=(const Self& src) = delete;
/// Returns the bounding box of the node
const Bounding_box& bbox() const { return m_bbox; }
/**
* @brief Builds the tree by recursive expansion.
* @param first the first primitive to insert
* @param last the last primitive to insert
* @param range the number of primitive of the range
*
* [first,last[ is the range of primitives to be added to the tree.
*/
template<typename ConstPrimitiveIterator>
void expand(ConstPrimitiveIterator first,
ConstPrimitiveIterator beyond,
const std::size_t range,
const AABBTraits&);
/**
* @brief General traversal query
* @param query the query
@ -93,8 +84,17 @@ public:
{ return *static_cast<Primitive*>(m_p_left_child); }
const Primitive& right_data() const
{ return *static_cast<Primitive*>(m_p_right_child); }
template <class Left, class Right>
void set_children(Left& l, Right& r)
{
m_p_left_child = static_cast<void*>(std::addressof(l));
m_p_right_child = static_cast<void*>(std::addressof(r));
}
void set_bbox(const Bounding_box& bbox)
{
m_bbox = bbox;
}
private:
Node& left_child() { return *static_cast<Node*>(m_p_left_child); }
Node& right_child() { return *static_cast<Node*>(m_p_right_child); }
Primitive& left_data() { return *static_cast<Primitive*>(m_p_left_child); }
@ -109,49 +109,8 @@ private:
void *m_p_left_child;
void *m_p_right_child;
private:
// Disabled copy constructor & assignment operator
typedef AABB_node<AABBTraits> Self;
AABB_node(const Self& src);
Self& operator=(const Self& src);
}; // end class AABB_node
template<typename Tr>
template<typename ConstPrimitiveIterator>
void
AABB_node<Tr>::expand(ConstPrimitiveIterator first,
ConstPrimitiveIterator beyond,
const std::size_t range,
const Tr& traits)
{
m_bbox = traits.compute_bbox_object()(first, beyond);
// sort primitives along longest axis aabb
traits.split_primitives_object()(first, beyond, m_bbox);
switch(range)
{
case 2:
m_p_left_child = &(*first);
m_p_right_child = &(*(++first));
break;
case 3:
m_p_left_child = &(*first);
m_p_right_child = static_cast<Node*>(this)+1;
right_child().expand(first+1, beyond, 2,traits);
break;
default:
const std::size_t new_range = range/2;
m_p_left_child = static_cast<Node*>(this) + 1;
m_p_right_child = static_cast<Node*>(this) + new_range;
left_child().expand(first, first + new_range, new_range,traits);
right_child().expand(first + new_range, beyond, range - new_range,traits);
}
}
template<typename Tr>
template<class Traversal_traits, class Query>
void

View File

@ -83,7 +83,7 @@ namespace CGAL
typedef typename CGAL::Orthogonal_k_neighbor_search<TreeTraits> Neighbor_search;
typedef typename Neighbor_search::Tree Tree;
private:
Tree* m_p_tree;
Tree m_tree;
Point_and_primitive_id get_p_and_p(const Point_and_primitive_id& p)
@ -98,30 +98,22 @@ namespace CGAL
public:
template <class ConstPointIterator>
AABB_search_tree(ConstPointIterator begin, ConstPointIterator beyond)
: m_p_tree(nullptr)
: m_tree{}
{
typedef typename Add_decorated_point<Traits, typename Traits::Primitive::Id>::Point_3 Decorated_point;
typedef typename Add_decorated_point<Traits,typename Traits::Primitive::Id>::Point_3 Decorated_point;
std::vector<Decorated_point> points;
while(begin != beyond) {
Point_and_primitive_id pp = get_p_and_p(*begin);
points.push_back(Decorated_point(pp.first,pp.second));
points.emplace_back(pp.first, pp.second);
++begin;
}
m_p_tree = new Tree(points.begin(), points.end());
if(m_p_tree != nullptr)
m_p_tree->build();
else
std::cerr << "unable to build the search tree!" << std::endl;
m_tree.insert(points.begin(), points.end());
m_tree.build();
}
~AABB_search_tree() {
delete m_p_tree;
}
Point_and_primitive_id closest_point(const Point& query) const
{
Neighbor_search search(*m_p_tree, query, 1);
Neighbor_search search(m_tree, query, 1);
return Point_and_primitive_id(static_cast<Point>(search.begin()->first), search.begin()->first.id());
}
};

View File

@ -0,0 +1,193 @@
#include <iostream>
#include <list>
#include <utility>
#include <CGAL/AABB_face_graph_triangle_primitive.h>
#include <CGAL/AABB_halfedge_graph_segment_primitive.h>
#include <CGAL/AABB_segment_primitive.h>
#include <CGAL/AABB_traits.h>
#include <CGAL/AABB_tree.h>
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
#include <CGAL/Polyhedron_3.h>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Surface_mesh.h>
#include <CGAL/assertions.h>
#include <fstream>
template <int test_number> auto create_tree();
template <int test_number, typename T> void init_tree(T &tree) {}
template <int test_number, typename T> bool test_tree(T &tree);
template <int test_number> class TestUtils;
// test 0 is from "aabb_test_singleton_tree"
template <> struct TestUtils<0> {
typedef CGAL::Simple_cartesian<double> K;
typedef K::FT FT;
typedef K::Point_3 Point;
typedef K::Plane_3 Plane;
typedef K::Segment_3 Segment;
typedef K::Triangle_3 Triangle;
typedef std::vector<Segment>::iterator Iterator;
typedef CGAL::AABB_segment_primitive<K, Iterator> Primitive;
typedef CGAL::AABB_traits<K, Primitive> Traits;
typedef CGAL::AABB_tree<Traits> Tree;
Point a = {1.0, 0.0, 0.0};
Point b = {0.0, 1.0, 0.0};
Point c = {0.0, 0.0, 1.0};
Point d = {0.0, 0.0, 0.0};
std::vector<Segment> segments = {Segment(Point(0, 0, 0), Point(2, 2, 2))};
};
template <> auto create_tree<0>() {
using T = TestUtils<0>;
T utils;
return T::Tree(utils.segments.begin(), utils.segments.end());
}
using Test0Param = decltype(create_tree<0>());
template <> bool test_tree<0, Test0Param>(Test0Param &tree) {
using T = TestUtils<0>;
T utils;
T::Plane plane_query(utils.a, utils.b, utils.d);
T::Triangle triangle_query(utils.a, utils.b, utils.c);
// Test calls to all functions
CGAL::Emptyset_iterator devnull;
tree.all_intersections(triangle_query, devnull);
tree.all_intersected_primitives(triangle_query, devnull);
assert(tree.any_intersected_primitive(triangle_query));
assert(tree.any_intersection(triangle_query));
const CGAL::Bbox_3 bbox = tree.bbox();
assert(bbox == CGAL::Bbox_3(0, 0, 0, 2, 2, 2));
tree.clear();
tree.insert(utils.segments.begin(), utils.segments.end());
tree.build();
assert(tree.closest_point(T::Point(-0.1, -0.1, -0.1)) == T::Point(0, 0, 0));
assert(tree.closest_point(T::Point(-0.1, -0.1, -0.1), T::Point(0, 0, 0)) ==
T::Point(0, 0, 0));
assert(tree.closest_point_and_primitive(T::Point(-0.1, -0.1, -0.1)).second ==
utils.segments.begin());
assert(tree.do_intersect(plane_query) == true);
assert(tree.do_intersect(triangle_query) == true);
assert(!tree.empty());
assert(tree.size() == 1);
tree.clear();
assert(tree.size() == 0);
tree.insert(utils.segments.begin(), utils.segments.end());
assert(tree.size() == 1);
assert(tree.number_of_intersected_primitives(plane_query) == 1);
tree.rebuild(utils.segments.begin(), utils.segments.end());
assert(tree.size() == 1);
assert(tree.number_of_intersected_primitives(triangle_query) == 1);
assert(tree.squared_distance(T::Point(0, 0, 0)) == 0);
return true;
}
// test 1 is from "aabb_test_all_intersected_primitives"
template <> struct TestUtils<1> {
typedef CGAL::Epick K;
typedef K::FT FT;
typedef K::Point_3 Point;
typedef K::Vector_3 Vector;
typedef K::Segment_3 Segment;
typedef K::Ray_3 Ray;
typedef CGAL::Surface_mesh<CGAL::Point_3<CGAL::Epick>> Mesh;
typedef CGAL::AABB_halfedge_graph_segment_primitive<Mesh, CGAL::Default,
CGAL::Tag_false>
S_Primitive;
typedef CGAL::AABB_face_graph_triangle_primitive<Mesh, CGAL::Default,
CGAL::Tag_false>
T_Primitive;
typedef CGAL::AABB_traits<K, T_Primitive> T_Traits;
typedef CGAL::AABB_traits<K, S_Primitive> S_Traits;
typedef CGAL::AABB_tree<T_Traits> T_Tree;
typedef CGAL::AABB_tree<S_Traits> S_Tree;
typedef T_Tree::Primitive_id T_Primitive_id;
typedef S_Tree::Primitive_id S_Primitive_id;
};
template <> auto create_tree<1>() {
using T = TestUtils<1>;
static CGAL::Surface_mesh<CGAL::Point_3<CGAL::Epick>> m1 = {};
static CGAL::Surface_mesh<CGAL::Point_3<CGAL::Epick>> m2 = {};
static bool mesh_loaded = false;
if (!mesh_loaded) {
std::ifstream in("data/cube.off");
assert(in);
in >> m1;
in.close();
in.open("data/tetrahedron.off");
assert(in);
in >> m2;
in.close();
mesh_loaded = true;
}
return std::make_pair(T::T_Tree{faces(m1).first, faces(m1).second, m1},
T::S_Tree{edges(m2).first, edges(m2).second, m2});
}
using Test1Param = decltype(create_tree<1>());
template <> void init_tree<1, Test1Param>(Test1Param &trees) {
trees.first.build();
trees.second.build();
}
template <> bool test_tree<1, Test1Param>(Test1Param &trees) {
using T = TestUtils<1>;
auto &cube_tree = trees.first;
auto &tet_tree = trees.second;
std::list<T::T_Tree::Primitive::Id> t_primitives;
std::list<T::S_Tree::Primitive::Id> s_primitives;
cube_tree.all_intersected_primitives(tet_tree,
std::back_inserter(t_primitives));
CGAL_assertion(t_primitives.size() == 6);
tet_tree.all_intersected_primitives(cube_tree,
std::back_inserter(s_primitives));
CGAL_assertion(s_primitives.size() == 6);
CGAL_assertion(tet_tree.do_intersect(cube_tree));
CGAL_assertion(cube_tree.do_intersect(tet_tree));
std::vector<T::T_Tree::Primitive::Id> all_primitives;
cube_tree.all_intersected_primitives(tet_tree,
std::back_inserter(all_primitives));
bool found_f5 = false;
for (auto prim : all_primitives) {
if ((int)prim.first == 5)
found_f5 = true;
}
CGAL_assertion(found_f5);
CGAL_USE(found_f5);
return true;
}
template <int test_number> bool run_test() {
// create_tree should return prvalue for guaranteed copy elision
auto tree_1 = create_tree<test_number>();
init_tree<test_number>(tree_1);
auto tree_2 = create_tree<test_number>();
init_tree<test_number>(tree_2);
auto tree_3 = create_tree<test_number>();
init_tree<test_number>(tree_3);
decltype(tree_1) tree_ctor{std::move(tree_2)};
decltype(tree_1) tree_assig{};
tree_assig = std::move(tree_3);
bool normal = test_tree<test_number>(tree_1);
bool move_ctor = test_tree<test_number>(tree_ctor);
bool move_ass =
test_tree<test_number>(tree_assig); // test move assignment operator
if (!normal)
std::cout << "Test " << test_number << "failed on the original tree\n";
if (!move_ctor)
std::cout << "Test " << test_number << "failed on move constructed tree\n";
if (!move_ass)
std::cout << "Test " << test_number << "failed on move assigned tree\n";
return normal && move_ctor && move_ass;
}
int main() { return (run_test<0>() && run_test<1>()) ? 0 : 1; }

View File

@ -6,6 +6,9 @@ Release History
Release date: June 2020
### Add Move Semantics to AABB Trees
- Added move constructor and assignment operator to AABB Tree class.
### Surface Mesh Topology (new package)
- This package allows to compute some topological invariants of