Pierre Alliez,
Stéphane Tayeb,
Camille Wormser
Figure 62.1: AABB tree. Left: surface triangle mesh of a mechanical part. Right: AABB tree constructed.
File: examples/AABB_tree/AABB_triangle_3_example.cpp
#include <iostream> #include <list> #include <CGAL/AABB_tree.h> // must be inserted before kernel #include <CGAL/AABB_traits.h> #include <CGAL/AABB_triangle_primitive.h> #include <CGAL/Simple_cartesian.h> typedef CGAL::Simple_cartesian<double> K; typedef K::FT FT; typedef K::Ray_3 Ray; typedef K::Line_3 Line; typedef K::Point_3 Point; typedef K::Triangle_3 Triangle; typedef std::list<Triangle>::iterator Iterator; typedef CGAL::AABB_triangle_primitive<K,Iterator> Primitive; typedef CGAL::AABB_traits<K, Primitive> AABB_triangle_traits; typedef CGAL::AABB_tree<AABB_triangle_traits> Tree; int main() { 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::list<Triangle> triangles; triangles.push_back(Triangle(a,b,c)); triangles.push_back(Triangle(a,b,d)); triangles.push_back(Triangle(a,d,c)); // constructs AABB tree Tree tree(triangles.begin(),triangles.end()); // counts #intersections Ray ray_query(a,b); std::cout << tree.number_of_intersected_primitives(ray_query) << " intersections(s) with ray query" << std::endl; // compute closest point and squared distance Point point_query(2.0, 2.0, 2.0); Point closest_point = tree.closest_point(point_query); FT sqd = tree.squared_distance(point_query); std::cout << "squared distance: " << sqd << std::endl; return EXIT_SUCCESS; }
File: examples/AABB_tree/AABB_polyhedron_facet_intersection_example.cpp
#include <iostream> #include <CGAL/AABB_tree.h> // must be inserted before kernel #include <CGAL/AABB_traits.h> #include <CGAL/Polyhedron_3.h> #include <CGAL/AABB_polyhedron_triangle_primitive.h> #include <CGAL/Simple_cartesian.h> typedef CGAL::Simple_cartesian<double> K; typedef K::Point_3 Point; typedef K::Plane_3 Plane; typedef K::Vector_3 Vector; typedef K::Segment_3 Segment; typedef CGAL::Polyhedron_3<K> Polyhedron; typedef CGAL::AABB_polyhedron_triangle_primitive<K,Polyhedron> Primitive; typedef CGAL::AABB_traits<K, Primitive> Traits; typedef CGAL::AABB_tree<Traits> Tree; typedef Tree::Object_and_primitive_id Object_and_primitive_id; typedef Tree::Primitive_id Primitive_id; int main() { Point p(1.0, 0.0, 0.0); Point q(0.0, 1.0, 0.0); Point r(0.0, 0.0, 1.0); Point s(0.0, 0.0, 0.0); Polyhedron polyhedron; polyhedron.make_tetrahedron(p, q, r, s); // constructs AABB tree Tree tree(polyhedron.facets_begin(),polyhedron.facets_end()); // constructs segment query Point a(-0.2, 0.2, -0.2); Point b(1.3, 0.2, 1.3); Segment segment_query(a,b); // tests intersections with segment query if(tree.do_intersect(segment_query)) std::cout << "intersection(s)" << std::endl; else std::cout << "no intersection" << std::endl; // computes #intersections with segment query std::cout << tree.number_of_intersected_primitives(segment_query) << " intersection(s)" << std::endl; // computes first encountered intersection with segment query // (generally a point) boost::optional<Object_and_primitive_id> intersection = tree.any_intersection(segment_query); if(intersection) { // gets intersection object Object_and_primitive_id op = *intersection; CGAL::Object object = op.first; Point point; if(CGAL::assign(point,object)) std::cout << "intersection object is a point" << std::endl; } // computes all intersections with segment query (as pairs object - primitive_id) std::list<Object_and_primitive_id> intersections; tree.all_intersections(segment_query, std::back_inserter(intersections)); // computes all intersected primitives with segment query as primitive ids std::list<Primitive_id> primitives; tree.all_intersected_primitives(segment_query, std::back_inserter(primitives)); // constructs plane query Vector vec(0.0,0.0,1.0); Plane plane_query(a,vec); // computes first encountered intersection with plane query // (generally a segment) intersection = tree.any_intersection(plane_query); if(intersection) { // gets intersection object Object_and_primitive_id op = *intersection; CGAL::Object object = op.first; Segment segment; if(CGAL::assign(segment,object)) std::cout << "intersection object is a segment" << std::endl; } return EXIT_SUCCESS; }
File: examples/AABB_tree/AABB_polyhedron_facet_distance_example.cpp
#include <iostream> #include <CGAL/AABB_tree.h> // must be inserted before kernel #include <CGAL/AABB_traits.h> #include <CGAL/Polyhedron_3.h> #include <CGAL/AABB_polyhedron_triangle_primitive.h> #include <CGAL/Simple_cartesian.h> typedef CGAL::Simple_cartesian<double> K; typedef K::FT FT; typedef K::Point_3 Point; typedef K::Segment_3 Segment; typedef CGAL::Polyhedron_3<K> Polyhedron; typedef CGAL::AABB_polyhedron_triangle_primitive<K,Polyhedron> Primitive; typedef CGAL::AABB_traits<K, Primitive> Traits; typedef CGAL::AABB_tree<Traits> Tree; typedef Tree::Object_and_primitive_id Object_and_primitive_id; typedef Tree::Point_and_primitive_id Point_and_primitive_id; int main() { Point p(1.0, 0.0, 0.0); Point q(0.0, 1.0, 0.0); Point r(0.0, 0.0, 1.0); Point s(0.0, 0.0, 0.0); Polyhedron polyhedron; polyhedron.make_tetrahedron(p, q, r, s); // constructs AABB tree and computes internal KD-tree // data structure to accelerate distance queries Tree tree(polyhedron.facets_begin(),polyhedron.facets_end()); tree.accelerate_distance_queries(); // query point Point query(0.0, 0.0, 3.0); // computes squared distance from query FT sqd = tree.squared_distance(query); std::cout << "squared distance: " << sqd << std::endl; // computes closest point Point closest = tree.closest_point(query); std::cout << "closest point: " << closest << std::endl; // computes closest point and primitive id Point_and_primitive_id pp = tree.closest_point_and_primitive(query); std::cout << "closest point: " << pp.first << std::endl; Polyhedron::Face_handle f = pp.second; // closest primitive id return EXIT_SUCCESS; }
File: examples/AABB_tree/AABB_segment_3_example.cpp
#include <iostream> #include <list> #include <CGAL/AABB_tree.h> // must be inserted before kernel #include <CGAL/AABB_traits.h> #include <CGAL/AABB_segment_primitive.h> #include <CGAL/Simple_cartesian.h> 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::list<Segment>::iterator Iterator; typedef CGAL::AABB_segment_primitive<K,Iterator> Primitive; typedef CGAL::AABB_traits<K, Primitive> Traits; typedef CGAL::AABB_tree<Traits> Tree; int main() { 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::list<Segment> segments; segments.push_back(Segment(a,b)); segments.push_back(Segment(a,c)); segments.push_back(Segment(c,d)); // constructs the AABB tree and the internal search tree for // efficient distance computations. Tree tree(segments.begin(),segments.end()); tree.accelerate_distance_queries(); // counts #intersections with a plane query Plane plane_query(a,b,d); std::cout << tree.number_of_intersected_primitives(plane_query) << " intersections(s) with plane" << std::endl; // counts #intersections with a triangle query Triangle triangle_query(a,b,c); std::cout << tree.number_of_intersected_primitives(triangle_query) << " intersections(s) with triangle" << std::endl; // computes the closest point from a point query Point point_query(2.0, 2.0, 2.0); Point closest = tree.closest_point(point_query); return EXIT_SUCCESS; }
File: examples/AABB_tree/AABB_polyhedron_edge_example.cpp
#include <iostream> #include <CGAL/AABB_tree.h> // must be inserted before kernel #include <CGAL/AABB_traits.h> #include <CGAL/Polyhedron_3.h> #include <CGAL/AABB_polyhedron_segment_primitive.h> #include <CGAL/Simple_cartesian.h> typedef CGAL::Simple_cartesian<double> K; typedef K::FT FT; typedef K::Point_3 Point; typedef K::Triangle_3 Triangle; typedef CGAL::Polyhedron_3<K> Polyhedron; typedef CGAL::AABB_polyhedron_segment_primitive<K,Polyhedron> Primitive; typedef CGAL::AABB_traits<K, Primitive> Traits; typedef CGAL::AABB_tree<Traits> Tree; int main() { Point p(1.0, 0.0, 0.0); Point q(0.0, 1.0, 0.0); Point r(0.0, 0.0, 1.0); Point s(0.0, 0.0, 0.0); Polyhedron polyhedron; polyhedron.make_tetrahedron(p, q, r, s); // constructs the AABB tree and the internal search tree for // efficient distance queries. Tree tree(polyhedron.edges_begin(),polyhedron.edges_end()); tree.accelerate_distance_queries(); // counts #intersections with a triangle query Triangle triangle_query(p,q,r); std::cout << tree.number_of_intersected_primitives(triangle_query) << " intersections(s) with triangle" << std::endl; // computes the closest point from a query point Point point_query(2.0, 2.0, 2.0); Point closest = tree.closest_point(point_query); return EXIT_SUCCESS; }
We provide some performance numbers for the case where the AABB tree contains a set of polyhedron triangle facets. We measure the tree construction time, the memory occupancy and the number of queries per second for a variety of intersection and distance queries. The machine used is a PC running Windows XP64 with an Intel CPU Core2 Extreme clocked at 3.06 GHz with 4GB of RAM. By default the kernel used is Simple_cartesian<double> (the fastest in our experiments). The program has been compiled with Visual C++ 2005 compiler with the O2 option which maximizes speed.
The surface triangle mesh chosen for benchmarking the tree construction is the knot model (14,400 triangles) depicted by Figure 62.4.3. We measure the tree construction time (both AABB tree alone and AABB tree with internal KD-tree) for this model as well as for three denser versions subdivided through the Loop subdivision scheme which multiplies the number of triangles by four.
Number of triangles | Construction (in ms) | Construction with internal KD-tree (in ms) |
14,400 | 156 | 157 |
57,600 | 328 | 328 |
230,400 | 1,141 | 1,437 |
921,600 | 4,813 | 5,953 |
When using the polyhedron triangle facet primitive (defined in AABB_polyhedron_triangle_primitive.h) the AABB tree occupies approximately 61 bytes per primitive (without constructing the internal KD-tree). It increases to approximately 150 bytes per primitive when constructing the internal KD-tree with one reference point per primitive (the default mode when calling the function tree.accelerate_distance_queries()). Note that the polyhedron facet primitive primitive stores only one facet handle as primitive id and computes on the fly a 3D triangle from the facet handle stored internally. When explicitly storing a 3D triangle in the primitive the tree occupies approximately 140 bytes per primitive instead of 60 (without constructing the internal KD-tree).
The following table provides orders of memory occupancy in MBytes for an increasing number of triangles. As the internal KD-tree used to accelerate the distance queries dominates the memory occupancy, we recommend to specify for large models a lower number of reference point (evenly distributed) to construct the internal KD-tree through the function tree.accelerate_distance_queries(begin,end) which takes an iterator range as input.
Number of triangles | AABB tree (in MBytes) | AABB tree with internal KD-tree (in MBytes) |
18,400 | 1.10 | 2.76 |
102,400 | 6.33 | 14.73 |
1,022,400 | 59.56 | 151.31 |
1,822,400 | 108.34 | 291.84 |
The following table measures the number of intersection queries per second on the 14,400 triangle version of the knot mesh model for ray, line, segment and plane queries. Each ray query is generated by choosing a random source point within the mesh bounding box and a random vector. A line or segment query is generated by choosing two random points inside the bounding box. A plane query is generated by picking a random point inside the bounding box and a random normal vector. Note that a plane query generally intersects many triangles of the input surface mesh. This explains the low performance numbers for the intersection functions which enumerate all intersections.
Function | Segment | Ray | Line | Plane |
do_intersect() | 187,868 | 185,649 | 206,096 | 377,969 |
any_intersected_primitive() | 190,684 | 190,027 | 208,941 | 360,337 |
any_intersection() | 147,468 | 143,230 | 148,235 | 229,336 |
number_of_intersected_primitives() | 64,389 | 52,943 | 54,559 | 7,906 |
all_intersected_primitives() | 65,553 | 54,838 | 53,183 | 5,693 |
all_intersections() | 46,507 | 38,471 | 36,374 | 2,644 |
Curve 62.4.3 plots the number of queries per second (here the all_intersections function with random segment queries) against the number of input triangles for the knot triangle surface mesh.
Figure 62.2: Number of queries per second against number of triangles for the knot model with 14K (shown), 57K, 230K and 921K triangles. We call the all_intersections function with segment queries randomly chosen within the bounding box.
The following table measures the number of all_intersections() queries per second against several kernels. We use the 14,400 triangle version of the knot mesh model for random segment queries. Note how the Simple_cartesian kernel is substantially faster than the Cartesian kernel.
Kernel | Queries/s (all_intersections() with segment queries) |
Simple_cartesian<double> | 46,507 |
Simple_cartesian<float> | 43,187 |
Cartesian<double> | 5,335 |
Cartesian<float> | 5,522 |
Exact_predicates_inexact_constructions | 18,411 |
The surface triangle mesh chosen for benchmarking distances is again the knot model in four increasing resolutions obtained through Loop subdivision. In the following table we first measure the tree construction time (which includes the construction of the internal KD-tree data structure used to accelerate the distance queries by up to one order of magnitude in our experiments). We then measure the number of queries per second for the three types distance queries (closest_point, squared_distance and closest_point_and_primitive) from point queries randomly chosen inside the bounding box.
Nb triangles | Construction (ms) | Closest_point() | squared_distance() | closest_point_and_primitive() |
14,400 | 157 | 45,132 | 45,626 | 45,770 |
57,600 | 328 | 21,589 | 21,312 | 21,137 |
230,400 | 1,437 | 11,063 | 10,962 | 11,086 |
921,600 | 5,953 | 5,636 | 5,722 | 5,703 |
The experiments described above are neither exhaustive nor conclusive as we have chosen one specific case where the input primitives are the facets of a triangle surface polyhedron. Nevertheless we now provide some general observations and advices about how to put the AABB tree to use with satisfactory performances. While the tree construction times and memory occupancy do not fluctuate much in our experiments depending on the input surface triangle mesh, the performance expressed in number of queries varies greatly depending on a complex combination of criteria: type of kernel, number of input primitives, distribution of primitives in space, type of function queried, type of query and location of query in space.
Camille Wormser and Pierre Alliez started working on a data structure for efficient collision detection in 2007. The generic design for implementing both intersection and distance queries, and for generic queries and primitives was developed by Camille Wormser. In 2009, Pierre Alliez, Stéphane Tayeb and Camille Wormser made the implementation CGAL-compliant, with the help of Laurent Rineau for optimizing the tree construction. The authors wish to thank Andreas Fabri, Jane Tournois, Mariette Yvinec and Sylvain Lefèbvre for helpful comments and discussions.