\( \newcommand{\E}{\mathrm{E}} \) \( \newcommand{\A}{\mathrm{A}} \) \( \newcommand{\R}{\mathrm{R}} \) \( \newcommand{\N}{\mathrm{N}} \) \( \newcommand{\Q}{\mathrm{Q}} \) \( \newcommand{\Z}{\mathrm{Z}} \) \( \def\ccSum #1#2#3{ \sum_{#1}^{#2}{#3} } \def\ccProd #1#2#3{ \sum_{#1}^{#2}{#3} }\)
CGAL 4.12.2 - 3D Fast Intersection and Distance Computation (AABB Tree)
AABB_tree/AABB_triangle_3_example.cpp
// Author(s) : Camille Wormser, Pierre Alliez
#include <iostream>
#include <list>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/AABB_tree.h>
#include <CGAL/AABB_traits.h>
#include <CGAL/AABB_triangle_primitive.h>
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_traits<K, Primitive> AABB_triangle_traits;
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);
std::cerr << "closest point is: " << closest_point << std::endl;
FT sqd = tree.squared_distance(point_query);
std::cout << "squared distance: " << sqd << std::endl;
return EXIT_SUCCESS;
}