CGAL 5.5.5 - 3D Fast Intersection and Distance Computation (AABB Tree)
CGAL::AABB_tree< AABBTraits > Class Template Reference

#include <CGAL/AABB_tree.h>

Inherited by CGAL::AABB_trees::internal::AABB_covered_triangle_tree< Kernel, Point >.

Definition

Public Member Functions

AABBTraits::Primitive::Datum_reference datum (Primitive &p) const
 Returns the datum (geometric object) represented p.
 

Types

typedef AABBTraits::FT FT
 Number type returned by the distance queries.
 
typedef AABBTraits::Point_3 Point
 Type of 3D point.
 
typedef AABBTraits::Primitive Primitive
 Type of input primitive.
 
typedef Primitive::Id Primitive_id
 Identifier for a primitive in the tree.
 
typedef Primitives::size_type size_type
 Unsigned integral size type.
 
typedef AABBTraits::Bounding_box Bounding_box
 Type of bounding box.
 
typedef AABBTraits::Point_and_primitive_id Point_and_primitive_id
 3D Point and Primitive Id type
 
typedef AABBTraits::Object_and_primitive_id Object_and_primitive_id
 
template<typename Query >
using Intersection_and_primitive_id = AABBTraits::Intersection_and_primitive_id< Query >
 An alias to AABBTraits::Intersection_and_primitive_id<Query>
 

Creation

 AABB_tree (const AABBTraits &traits=AABBTraits())
 constructs an empty tree, and initializes the internally stored traits class using traits. More...
 
 AABB_tree (Self &&) noexcept
 move constructor
 
Selfoperator= (Self &&) noexcept
 assignment operator
 
 AABB_tree (const Self &)=delete
 
Selfoperator= (const Self &)=delete
 
template<typename InputIterator , typename ... T>
 AABB_tree (InputIterator first, InputIterator beyond, T &&...)
 Builds the data structure from a sequence of primitives. More...
 
template<typename ... T>
void build (T &&...)
 triggers the (re)construction of the internal tree structure. More...
 

Operations

template<typename ConstPrimitiveIterator , typename ... T>
void rebuild (ConstPrimitiveIterator first, ConstPrimitiveIterator beyond, T &&...)
 is equivalent to calling clear(), insert(first,last,t...), and build()
 
template<typename InputIterator , typename ... T>
void insert (InputIterator first, InputIterator beyond, T &&...)
 adds a sequence of primitives to the set of primitives of the AABB tree. More...
 
void insert (const Primitive &p)
 adds a primitive to the set of primitives of the tree.
 
 ~AABB_tree ()
 clears and destroys the tree.
 
const AABBTraitstraits () const
 returns a const reference to the internally stored traits class.
 
void clear ()
 clears the tree and the search tree if it was constructed, and switches on the usage of the search tree to find the hint for the distance queries
 
const Bounding_box bbox () const
 returns the axis-aligned bounding box of the whole tree. More...
 
size_type size () const
 returns the number of primitives in the tree.
 
bool empty () const
 returns true, iff the tree contains no primitive.
 

Intersection Tests

template<typename Query >
bool do_intersect (const Query &query) const
 returns true, iff the query intersects at least one of the input primitives. More...
 
template<typename Query >
size_type number_of_intersected_primitives (const Query &query) const
 returns the number of primitives intersected by the query. More...
 
template<typename Query , typename OutputIterator >
OutputIterator all_intersected_primitives (const Query &query, OutputIterator out) const
 puts in out the ids of all intersected primitives. More...
 
template<typename Query >
boost::optional< Primitive_idany_intersected_primitive (const Query &query) const
 returns the id of the intersected primitive that is encountered first in the tree traversal, iff the query intersects at least one of the input primitives. More...
 

Intersections

template<typename Query , typename OutputIterator >
OutputIterator all_intersections (const Query &query, OutputIterator out) const
 puts in out all intersections, as objects of Intersection_and_primitive_id<Query>::Type, between the query and the input data to the iterator. More...
 
template<typename Query >
boost::optional< typename Intersection_and_primitive_id< Query >::Type > any_intersection (const Query &query) const
 returns if any the intersection that is encountered first in the tree traversal. More...
 
template<typename Ray , typename SkipFunctor >
boost::optional< typename Intersection_and_primitive_id< Ray >::Type > first_intersection (const Ray &query, const SkipFunctor &skip) const
 returns the intersection and primitive id closest to the source point of the ray query. More...
 
template<typename Ray , typename SkipFunctor >
boost::optional< Primitive_idfirst_intersected_primitive (const Ray &query, const SkipFunctor &skip) const
 returns the primitive id closest to the source point of the ray query. More...
 

Distance Queries

FT squared_distance (const Point &query) const
 returns the minimum squared distance between the query point and all input primitives. More...
 
Point closest_point (const Point &query) const
 returns the point in the union of all input primitives which is closest to the query. More...
 
Point_and_primitive_id closest_point_and_primitive (const Point &query) const
 returns a Point_and_primitive_id which realizes the smallest distance between the query point and all input primitives. More...
 

Accelerating the Distance Queries

In the following paragraphs, we discuss details of the implementation of the distance queries.

We explain the internal use of hints, how the user can pass his own hints to the tree, and how the user can influence the construction of the secondary data structure used for accelerating distance queries. Internally, the distance queries algorithms are initialized with some hint, which has the same type as the return type of the query, and this value is refined along a traversal of the tree, until it is optimal, that is to say until it realizes the shortest distance to the primitives. In particular, the exact specification of these internal algorithms is that they minimize the distance to the object composed of the union of the primitives and the hint. It follows that

  • in order to return the exact distance to the set of primitives, the algorithms need the hint to be exactly on the primitives;
  • if this is not the case, and if the hint happens to be closer to the query point than any of the primitives, then the hint is returned.

This second observation is reasonable, in the sense that providing a hint to the algorithm means claiming that this hint belongs to the union of the primitives. These considerations about the hints being exactly on the primitives or not are important: in the case where the set of primitives is a triangle soup, and if some of the primitives are large, one may want to provide a much better hint than a vertex of the triangle soup could be. It could be, for example, the barycenter of one of the triangles. But, except with the use of a kernel with exact constructions, one cannot easily construct points other than the vertices, that lie exactly on a triangle soup. Hence, providing a good hint sometimes means not being able to provide it exactly on the primitives. In rare occasions, this hint can be returned as the closest point. In order to accelerate distance queries significantly, the AABB tree builds an internal KD-tree containing a set of potential hints. This KD-tree provides very good hints that allow the algorithms to run much faster than when do_not_accelerate_distance_queries() that makes the hint to always be the reference_point of the first primitive. The set of potential hints is a sampling of the union of the primitives, which is obtained, by default, by calling the method reference_point of each of the primitives. However, such a sampling with one point per primitive may not be the most relevant one: if some primitives are very large, it helps inserting more than one sample on them. Conversely, a sparser sampling with less than one point per input primitive is relevant in some cases. The internal KD-tree is always used if no call to do_not_accelerate_distance_queries() was made since object creation or the last call to clear(). It will be built by the first distance query or by a call to accelerate_distance_queries().

bool accelerate_distance_queries ()
 constructs the internal search tree from a point set taken on the internal primitives returns true iff successful memory allocation
 
void do_not_accelerate_distance_queries ()
 turns off the usage of the internal search tree and clears it if it was already constructed.
 
template<typename ConstPointIterator >
bool accelerate_distance_queries (ConstPointIterator first, ConstPointIterator beyond)
 constructs an internal KD-tree containing the specified point set, to be used as the set of potential hints for accelerating the distance queries. More...
 
FT squared_distance (const Point &query, const Point &hint) const
 returns the minimum squared distance between the query point and all input primitives. More...
 
Point closest_point (const Point &query, const Point &hint) const
 returns the point in the union of all input primitives which is closest to the query. More...
 
Point_and_primitive_id closest_point_and_primitive (const Point &query, const Point_and_primitive_id &hint) const
 returns a Point_and_primitive_id which realizes the smallest distance between the query point and all input primitives. More...
 

Member Typedef Documentation

◆ Object_and_primitive_id

Constructor & Destructor Documentation

◆ AABB_tree() [1/2]

template<typename AABBTraits >
CGAL::AABB_tree< AABBTraits >::AABB_tree ( const AABBTraits traits = AABBTraits())

constructs an empty tree, and initializes the internally stored traits class using traits.

◆ AABB_tree() [2/2]

template<typename AABBTraits >
template<typename InputIterator , typename ... T>
CGAL::AABB_tree< AABBTraits >::AABB_tree ( InputIterator  first,
InputIterator  beyond,
T &&  ... 
)

Builds the data structure from a sequence of primitives.

Parameters
firstiterator over first primitive to insert
beyondpast-the-end iterator

constructs an empty tree followed by a call to insert(first,last,t...). The tree stays empty if the memory allocation is not successful.

Member Function Documentation

◆ accelerate_distance_queries()

template<typename AABBTraits >
template<typename ConstPointIterator >
bool CGAL::AABB_tree< AABBTraits >::accelerate_distance_queries ( ConstPointIterator  first,
ConstPointIterator  beyond 
)

constructs an internal KD-tree containing the specified point set, to be used as the set of potential hints for accelerating the distance queries.

Note that the search tree built in this function will not be invalidated by the insertion of a new primitive, and an explicit call to accelerate_distance_queries() is needed to update the search tree.

Template Parameters
ConstPointIteratoris an iterator with value type Point_and_primitive_id.

◆ all_intersected_primitives()

template<typename Tr >
template<typename Query , typename OutputIterator >
OutputIterator CGAL::AABB_tree< Tr >::all_intersected_primitives ( const Query &  query,
OutputIterator  out 
) const

puts in out the ids of all intersected primitives.

This function does not compute the intersection points and is hence faster than the function all_intersections() function below.

Template Parameters
Querymust be a type for which Do_intersect operators are defined in the traits class AABBTraits.

◆ all_intersections()

template<typename Tr >
template<typename Query , typename OutputIterator >
OutputIterator CGAL::AABB_tree< Tr >::all_intersections ( const Query &  query,
OutputIterator  out 
) const

puts in out all intersections, as objects of Intersection_and_primitive_id<Query>::Type, between the query and the input data to the iterator.

Template Parameters
Querymust be a type for which Do_intersect and Intersection operators are defined in the traits class AABBTraits.

◆ any_intersected_primitive()

template<typename AABBTraits >
template<typename Query >
boost::optional<Primitive_id> CGAL::AABB_tree< AABBTraits >::any_intersected_primitive ( const Query &  query) const

returns the id of the intersected primitive that is encountered first in the tree traversal, iff the query intersects at least one of the input primitives.

No particular order is guaranteed over the tree traversal, such that, e.g, the primitive returned is not necessarily the closest from the source point of a ray query.

Template Parameters
Querymust be a type for which Do_intersect operators are defined in the traits class AABBTraits.

◆ any_intersection()

template<typename AABBTraits >
template<typename Query >
boost::optional< typename Intersection_and_primitive_id<Query>::Type > CGAL::AABB_tree< AABBTraits >::any_intersection ( const Query &  query) const

returns if any the intersection that is encountered first in the tree traversal.

No particular order is guaranteed over the tree traversal, e.g, the primitive returned is not necessarily the closest from the query.

Template Parameters
Querymust be a type for which Do_intersect and Intersection operators are defined in the traits class AABBTraits.

◆ bbox()

template<typename AABBTraits >
const Bounding_box CGAL::AABB_tree< AABBTraits >::bbox ( ) const

returns the axis-aligned bounding box of the whole tree.

Precondition
!empty()

◆ build()

template<typename AABBTraits >
template<typename ... T>
void CGAL::AABB_tree< AABBTraits >::build ( T &&  ...)

triggers the (re)construction of the internal tree structure.

The internal tree structure is automatically invalidated by the insertion of any primitives after one or more calls to insert(). This procedure is called implicitly at the first call to a query member function. An explicit call to build() must be made to ensure that the next call to a query function will not trigger the construction of the data structure. A call to AABBTraits::set_shared_data(t...) is made using the internally stored traits. This procedure has a complexity of \(O(n log(n))\), where \(n\) is the number of primitives of the tree.

◆ closest_point() [1/2]

template<typename Tr >
AABB_tree< Tr >::Point CGAL::AABB_tree< Tr >::closest_point ( const Point query) const

returns the point in the union of all input primitives which is closest to the query.

In case there are several closest points, one arbitrarily chosen closest point is returned.

Precondition
!empty()

◆ closest_point() [2/2]

template<typename Tr >
AABB_tree< Tr >::Point CGAL::AABB_tree< Tr >::closest_point ( const Point query,
const Point hint 
) const

returns the point in the union of all input primitives which is closest to the query.

In case there are several closest points, one arbitrarily chosen closest point is returned. The internal KD-tree is not used.

Precondition
!empty()

◆ closest_point_and_primitive() [1/2]

template<typename Tr >
AABB_tree< Tr >::Point_and_primitive_id CGAL::AABB_tree< Tr >::closest_point_and_primitive ( const Point query) const

returns a Point_and_primitive_id which realizes the smallest distance between the query point and all input primitives.

Precondition
!empty()

◆ closest_point_and_primitive() [2/2]

template<typename Tr >
AABB_tree< Tr >::Point_and_primitive_id CGAL::AABB_tree< Tr >::closest_point_and_primitive ( const Point query,
const Point_and_primitive_id hint 
) const

returns a Point_and_primitive_id which realizes the smallest distance between the query point and all input primitives.

The internal KD-tree is not used.

Precondition
!empty()

◆ do_intersect()

template<typename Tr >
template<typename Query >
bool CGAL::AABB_tree< Tr >::do_intersect ( const Query &  query) const

returns true, iff the query intersects at least one of the input primitives.

Template Parameters
Querymust be a type for which Do_intersect operators are defined in the traits class AABBTraits.

◆ first_intersected_primitive()

template<typename AABBTraits >
template<typename Ray , typename SkipFunctor >
boost::optional<Primitive_id> CGAL::AABB_tree< AABBTraits >::first_intersected_primitive ( const Ray &  query,
const SkipFunctor &  skip 
) const

returns the primitive id closest to the source point of the ray query.

Template Parameters
Raymust be the same as AABBTraits::Ray_3 and do_intersect predicates and intersections for it must be defined.
Skipa functor with an operator bool operator()(const Primitive_id& id) const that returns true in order to skip the primitive. Defaults to a functor that always returns false.

AABBTraits must be a model of AABBRayIntersectionTraits to call this member function.

◆ first_intersection()

template<typename AABBTraits >
template<typename Ray , typename SkipFunctor >
boost::optional< typename Intersection_and_primitive_id<Ray>::Type > CGAL::AABB_tree< AABBTraits >::first_intersection ( const Ray &  query,
const SkipFunctor &  skip 
) const

returns the intersection and primitive id closest to the source point of the ray query.

Template Parameters
Raymust be the same as AABBTraits::Ray_3 and do_intersect predicates and intersections for it must be defined.
Skipa functor with an operator bool operator()(const Primitive_id& id) const that returns true in order to skip the primitive. Defaults to a functor that always returns false.
Note
skip might be given some primitives that are not intersected by query because the intersection test is done after the skip test. Also note that the order the primitives are given to skip is not necessarily the intersection order with query.

AABBTraits must be a model of AABBRayIntersectionTraits to call this member function.

◆ insert()

template<typename AABBTraits >
template<typename InputIterator , typename ... T>
void CGAL::AABB_tree< AABBTraits >::insert ( InputIterator  first,
InputIterator  beyond,
T &&  ... 
)

adds a sequence of primitives to the set of primitives of the AABB tree.

InputIterator is any iterator and the parameter pack T contains any types such that Primitive has a constructor with the following signature: Primitive(InputIterator, T...). If Primitive is a model of the concept AABBPrimitiveWithSharedData, a call to AABBTraits::set_shared_data(t...) is made using the internally stored traits.

◆ number_of_intersected_primitives()

template<typename AABBTraits >
template<typename Query >
size_type CGAL::AABB_tree< AABBTraits >::number_of_intersected_primitives ( const Query &  query) const

returns the number of primitives intersected by the query.

Template Parameters
Querymust be a type for which Do_intersect operators are defined in the traits class AABBTraits.

◆ squared_distance() [1/2]

template<typename Tr >
AABB_tree< Tr >::FT CGAL::AABB_tree< Tr >::squared_distance ( const Point query) const

returns the minimum squared distance between the query point and all input primitives.

Precondition
!empty()

◆ squared_distance() [2/2]

template<typename Tr >
AABB_tree< Tr >::FT CGAL::AABB_tree< Tr >::squared_distance ( const Point query,
const Point hint 
) const

returns the minimum squared distance between the query point and all input primitives.

The internal KD-tree is not used.

Precondition
!empty()