CGAL 6.0 - Quadtrees, Octrees, and Orthtrees
Loading...
Searching...
No Matches
User Manual

Authors
Jackson Campolattaro, Simon Giraudot, Cédric Portaneri, Tong Zhao, Pierre Alliez

Introduction

Quadtrees are tree data structures in which each node encloses a rectangular section of space, and each internal node has exactly 4 children. Octrees are a similar data structure in 3D in which each node encloses a rectangular cuboid section of space, and each internal node has exactly 8 children.

We call the generalization of such data structure "orthtrees", as orthants are generalizations of quadrants and octants. The term "hyperoctree" can also be found in literature to name such data structures in dimensions 4 and higher.

This package provides a general data structure Orthtree along with aliases for Quadtree and Octree. These trees can be constructed with custom contents and split predicates, and iterated on with various traversal methods. Orthants can be orthotopes and not only hypercubes.


Figure 94.1 Top: an orthtree in 3D (octree) from a point cloud (top); Bottom: an orthtree in 3D (octree) from the triangle faces of a mesh.


Building

A common purpose for an orthtree is to subdivide a collection of points, and the Orthtree package provides a traits class for this purpose. The points are not copied: the provided point range is used directly and rearranged by the orthtree. Altering the point range after creating the orthtree may leave the tree in an invalid state. The constructor returns a tree with a single (root) node that contains all the points.

The method refine() must be called to subdivide space further. This method uses a split predicate which takes a node as input and returns true if this node should be split, false otherwise: this enables users to choose on what criterion should the orthtree be refined. Predefined predicates are provided for the point-set orthtree, including Maximum_depth and Maximum_contained_elements.

The simplest way to create a point-set orthtree is from a vector of points. The constructor generally expects a separate point range and map, but the point map defaults to Identity_property_map if none is provided.

The split predicate is a user-defined functor that determines whether a node needs to be split. Custom predicates can easily be defined if the existing ones do not match users' needs.

Building a Quadtree

The Orthtree class may be templated with Orthtree_traits_point<> while specifying a 2d ambient space and thus behave as a quadtree. For convenience, the alias CGAL::Quadtree is provided.

The following example shows the construction of quadtree from a vector of Point_2 objects. quadtree.refine(10, 5) uses the default split predicate, which enforces a max-depth and a maximum of elements contained per node ("bucket size"). Nodes are split if their depth is less than 10, and they contain more than 5 points.


File Orthtree/quadtree_build_from_point_vector.cpp

#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Quadtree.h>
#include <CGAL/Random.h>
// Type Declarations
using Point_vector = std::vector<Point_2>;
int main()
{
CGAL::Random r;
Point_vector points_2d;
for (std::size_t i = 0; i < 5; ++ i)
points_2d.emplace_back(r.get_double(-1., 1.),
r.get_double(-1., 1.));
Quadtree quadtree(points_2d);
quadtree.refine(10, 5);
return EXIT_SUCCESS;
}
A data structure using an axis-aligned hyperrectangle decomposition of dD space for efficient access ...
Definition: Orthtree.h:117
Orthtree< Orthtree_traits_point< GeomTraits, PointRange, PointMap, squared_nodes, 2 > > Quadtree
Alias that specializes the Orthtree class to a 2D quadtree storing 2D points.
Definition: Quadtree.h:38

Building an Octree

Orthtree_traits_point<> can also be templated with dimension 3 and thus behave as an octree. For convenience, the alias CGAL::Octree is provided.

The following example shows how to create an octree from a vector of Point_3 objects. As with the quadtree example, we use the default split predicate. In this case the maximum depth is 10, and the bucket size is 1.


File Orthtree/octree_build_from_point_vector.cpp

#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
// Type Declarations
using Point = Kernel::Point_3;
using Point_vector = std::vector<Point>;
int main() {
// Here, our point set is a vector
Point_vector points;
// Add a few points to the vector
points.emplace_back(1, 1, 1);
points.emplace_back(2, 1, -11);
points.emplace_back(2, 1, 1);
points.emplace_back(1, -2, 1);
points.emplace_back(1, 1, 1);
points.emplace_back(-1, 1, 1);
// Create an octree from the points
Octree octree(points);
// Build the octree
octree.refine(10, 1);
// Print out the tree
std::cout << octree;
return EXIT_SUCCESS;
}
Orthtree< Orthtree_traits_point< GeomTraits, PointRange, PointMap, cubic_nodes, 3 > > Octree
Alias that specializes the Orthtree class to a 3D octree storing 3D points.
Definition: Octree.h:38

Building an Octree from a Point_set_3

Some data structures such as Point_set_3 require a non-default point map type and object. This example illustrates how to create an octree from a Point_set_3 loaded from a file. It also shows a more explicit way of setting the split predicate when refining the tree.

An octree is constructed from the point set and its corresponding map, and then written to the standard output.

The split predicate is manually constructed and passed to the refine method. In this case, we use a maximum number of contained elements with no corresponding maximum depth. This means that nodes will continue to be split until none contain more than 10 points.


File Orthtree/octree_build_from_point_set.cpp

#include <fstream>
#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
#include <CGAL/Point_set_3.h>
#include <CGAL/Point_set_3/IO.h>
// Type Declarations
using Point = Kernel::Point_3;
using Point_set = CGAL::Point_set_3<Point>;
using Point_map = Point_set::Point_map;
int main(int argc, char **argv) {
// Point set will be used to hold our points
Point_set points;
// Load points from a file.
std::ifstream stream((argc > 1) ? argv[1] : CGAL::data_file_path("points_3/cube.pwn"));
stream >> points;
if (0 == points.number_of_points()) {
std::cerr << "Error: cannot read file" << std::endl;
return EXIT_FAILURE;
}
std::cout << "loaded " << points.number_of_points() << " points\n" << std::endl;
// Create an octree from the points
Octree octree(points, points.point_map());
// Build the octree with a small bucket size, using a more verbose method
// Print out the tree
std::cout << octree << std::endl;
return EXIT_SUCCESS;
}
A class used to choose when a node should be split depending on the number of contained elements.
Definition: Split_predicates.h:36

Building an Octree with a Custom Split Predicate

The following example illustrates how to refine an octree using a split predicate that isn't provided by default. This particular predicate sets a node's bucket size as a ratio of its depth. For example, for a ratio of 2, a node at depth 2 can hold 4 points, a node at depth 7 can hold 14.


File Orthtree/octree_build_with_custom_split.cpp

#include <fstream>
#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
#include <CGAL/Point_set_3.h>
#include <CGAL/Point_set_3/IO.h>
// Type Declarations
using Point = Kernel::Point_3;
using FT = Kernel::FT;
using Point_set = CGAL::Point_set_3<Point>;
using Point_map = Point_set::Point_map;
// Split Predicate
// The predicate is a functor which returns a Boolean value, whether a node needs to be split or not
struct Split_by_ratio {
std::size_t ratio;
explicit Split_by_ratio(std::size_t ratio) : ratio(ratio) {}
template<typename Node_index, typename Tree>
bool operator()(Node_index i, const Tree &tree) const {
std::size_t num_points = tree.data(i).size();
std::size_t depth = tree.depth(i);
return num_points > (ratio * depth);
}
};
int main(int argc, char **argv) {
// Point set will be used to hold our points
Point_set points;
// Load points from a file.
std::ifstream stream((argc > 1) ? argv[1] : CGAL::data_file_path("points_3/cube.pwn"));
stream >> points;
if (0 == points.number_of_points()) {
std::cerr << "Error: cannot read file" << std::endl;
return EXIT_FAILURE;
}
std::cout << "loaded " << points.number_of_points() << " points" << std::endl;
// Create an octree from the points
Octree octree(points, points.point_map());
// Build the octree using our custom split predicate
octree.refine(Split_by_ratio(2));
// Print out the tree
std::cout << octree;
return EXIT_SUCCESS;
}

Building an Orthtree

An orthtree can also be used with an arbitrary number of dimensions. The Orthtree_traits_point template can infer the arbitrary dimension count from the d-dimensional kernel.

The following example shows how to build a generalized orthtree in dimension 4. As std::vector<Point_d> is manually filled with 4-dimensional points. The vector is used as the point set, and an Identity_property_map is automatically set as the orthtree's map type, so a map does not need to be provided.


File Orthtree/orthtree_build.cpp

#include <iostream>
#include <CGAL/Epick_d.h>
#include <CGAL/Cartesian_d.h>
#include <CGAL/Orthtree.h>
#include <CGAL/Orthtree_traits_point.h>
#include <CGAL/Random.h>
// Type Declarations
const int dimension = 4;
using Point_d = Kernel::Point_d;
using Point_vector = std::vector<Point_d>;
using Orthtree = CGAL::Orthtree<Traits>;
int main()
{
CGAL::Random r;
Point_vector points_dd;
for (std::size_t i = 0; i < 20; ++ i)
{
std::array<double, dimension> init{};
for (double& v : init)
v = r.get_double(-1., 1.);
points_dd.emplace_back (init.begin(), init.end());
}
Orthtree orthtree(points_dd);
orthtree.refine(10, 5);
std::cout << orthtree << std::endl;
return EXIT_SUCCESS;
}
Traits class for defining an orthtree of points using the class CGAL::Orthtree.
Definition: Orthtree_traits_point.h:85

Properties

The Orthtree uses a mechanism to attach properties to nodes at run-time which follows Surface mesh properties. Each property is identified by a string and its value type and stored in a consecutive block of memory.

Traversal

Traversal is the act of navigating among the nodes of the tree. The Orthtree class provides a number of different solutions for traversing the tree.

Manual Traversal

Because our orthtree is a form of connected acyclic undirected graph, it is possible to navigate between any two nodes. What that means in practice is that given a node on the tree, it is possible to access any other node using the right set of operations. The Node_index type provides a handle on a node, and the orthtree class provides methods that enable the user to retrieve the indices of each of its children as well as its parent (if they exist).

Given any node index nid, the nth child of that node can be found with orthtree.child(nid, n). For an octree, values of n from 0-7 provide access to the different children. For non-root nodes, it is possible to access parent nodes using the orthtree.parent() accessor.

To access grandchildren, it isn't necessary to use nested orthtree.child() calls. Instead, the orthtree.descendant(node, a, b, ...) convenience method is provided. This retrieves the bth child of the ath child, to any depth.

In most cases, we want to find the descendants of the root node. For this case, there is another convenience method orthtree.node(a, b, c, ...). This retrieves the node specified by the descent a, b, c. It is equivalent to orthtree.descendant(orthtree.root(), a, b, c, ...).

The following example demonstrates the use of several of these accessors.


File Orthtree/octree_traversal_manual.cpp

#include <fstream>
#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
#include <CGAL/Point_set_3.h>
#include <CGAL/Point_set_3/IO.h>
// Type Declarations
using Point = Kernel::Point_3;
using Point_set = CGAL::Point_set_3<Point>;
using Point_map = Point_set::Point_map;
int main(int argc, char** argv) {
// Point set will be used to hold our points
Point_set points;
// Load points from a file.
std::ifstream stream((argc > 1) ? argv[1] : CGAL::data_file_path("points_3/cube.pwn"));
stream >> points;
if (0 == points.number_of_points()) {
std::cerr << "Error: cannot read file" << std::endl;
return EXIT_FAILURE;
}
std::cout << "loaded " << points.number_of_points() << " points\n" << std::endl;
// Create an octree from the points
Octree octree(points, points.point_map());
// Build the octree using the default arguments
octree.refine();
// Print out a few nodes
std::cout << "Navigation relative to the root node" << std::endl;
std::cout << "~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~" << std::endl;
std::cout << "the root node: " << std::endl;
std::cout << octree.to_string(octree.root()) << std::endl;
std::cout << "the first child of the root node: " << std::endl;
std::cout << octree.to_string(octree.child(octree.root(), 0)) << std::endl;
std::cout << "the fifth child: " << std::endl;
std::cout << octree.to_string(octree.child(octree.root(), 4)) << std::endl;
std::cout << "the fifth child, accessed without going through root: " << std::endl;
std::cout << octree.to_string(octree.node(4)) << std::endl;
std::cout << "the second child of the fourth child: " << std::endl;
std::cout << octree.to_string(octree.child(octree.child(octree.root(), 4), 1)) << std::endl;
std::cout << "the second child of the fourth child, accessed without going through root: " << std::endl;
std::cout << octree.to_string(octree.node(4, 1)) << std::endl;
std::cout << std::endl;
// Retrieve one of the deeper children
Octree::Node_index cur = octree.node(4, 3);
std::cout << "Navigation relative to a child node" << std::endl;
std::cout << "~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~" << std::endl;
std::cout << "the third child of the fourth child: " << std::endl;
std::cout << octree.to_string(cur) << std::endl;
std::cout << "the third child: " << std::endl;
std::cout << octree.to_string(octree.parent(cur)) << std::endl;
std::cout << "the next sibling of the third child of the fourth child: " << std::endl;
std::cout << octree.to_string(octree.child(octree.parent(cur), octree.local_coordinates(cur).to_ulong() + 1))
<< std::endl;
return EXIT_SUCCESS;
}
typename Traits::Node_index Node_index
Index of a given node in the tree; the root always has index 0.
Definition: Orthtree.h:141

Preorder Traversal

It is often useful to be able to iterate over the nodes of the tree in a particular order. For example, the stream operator << uses a traversal to print out each node. A few traversals are provided, among them Preorder_traversal and Postorder_traversal. Traversing a tree in preorder means to visit each parent immediately followed by its children, whereas in postorder traversal the children are visited first.

The following example illustrates how to use the provided traversals.

A tree is constructed, and a traversal is used to create a range that can be iterated over using a for-each loop. The default output operator for the orthtree uses the preorder traversal to do a pretty-print of the tree structure. In this case, we print out the nodes of the tree without indentation instead.


File Orthtree/octree_traversal_preorder.cpp

#include <fstream>
#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
#include <CGAL/Point_set_3.h>
#include <CGAL/Point_set_3/IO.h>
// Type Declarations
using Point = Kernel::Point_3;
using Point_set = CGAL::Point_set_3<Point>;
using Point_map = Point_set::Point_map;
int main(int argc, char **argv) {
// Point set will be used to hold our points
Point_set points;
// Load points from a file.
std::ifstream stream((argc > 1) ? argv[1] : CGAL::data_file_path("points_3/cube.pwn"));
stream >> points;
if (0 == points.number_of_points()) {
std::cerr << "Error: cannot read file" << std::endl;
return EXIT_FAILURE;
}
std::cout << "loaded " << points.number_of_points() << " points\n" << std::endl;
// Create an octree from the points
Octree octree(points, points.point_map());
// Build the octree
octree.refine();
// Print out the octree using preorder traversal
for (auto node : octree.traverse<Preorder_traversal>()) {
std::cout << octree.to_string(node) << std::endl;
}
return EXIT_SUCCESS;
}
A class used for performing a preorder traversal.
Definition: Traversals.h:38

Custom Traversal

Users can define their own traversal methods by creating models of the OrthtreeTraversal concept. The custom traversal must provide a method which returns the starting point of the traversal (first_index()) and another method which returns the next node in the sequence (next_index()). next_index() returns an empty optional when the end of the traversal is reached. The following example shows how to define a custom traversal that only traverses the first branch an octree:


File Orthtree/octree_traversal_custom.cpp

#include <fstream>
#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
#include <CGAL/Orthtree/IO.h>
#include <CGAL/Point_set_3.h>
#include <CGAL/Point_set_3/IO.h>
// Type Declarations
using Point = Kernel::Point_3;
using Point_set = CGAL::Point_set_3<Point>;
using Point_map = Point_set::Point_map;
template <typename Tree>
struct First_branch_traversal {
using Node_index = typename Tree::Node_index;
const Tree& m_orthtree;
explicit First_branch_traversal(const Tree& orthtree) : m_orthtree(orthtree) {}
Node_index first_index() const {
return m_orthtree.root();
}
std::optional<Node_index> next_index(Node_index n) const {
// Stop when we reach the base of the tree
if (m_orthtree.is_leaf(n))
return {};
// Always descend the first child
return m_orthtree.child(n, 0);
}
};
int main(int argc, char** argv) {
// Point set will be used to hold our points
Point_set points;
// Load points from a file.
std::ifstream stream((argc > 1) ? argv[1] : CGAL::data_file_path("points_3/cube.pwn"));
stream >> points;
if (0 == points.number_of_points()) {
std::cerr << "Error: cannot read file" << std::endl;
return EXIT_FAILURE;
}
std::cout << "loaded " << points.number_of_points() << " points\n" << std::endl;
// Create an octree from the points
Octree octree(points, points.point_map());
// Build the octree
octree.refine();
// Print out the first branch using custom traversal
for (auto node: octree.traverse<First_branch_traversal<Octree>>()) {
std::cout << octree.to_string(node) << std::endl;
}
return EXIT_SUCCESS;
}

Comparison of Traversals

Figure Figure 94.2 shows in which order nodes are visited depending on the traversal method used.

Figure 94.2 Quadtree visualized as a graph. Each node is labelled according to the order in which it is visited by the traversal. When using leaves and level traversals, the quadtree is only partially traversed.


Acceleration of Common Tasks

Once an orthtree is built, its structure can be used to accelerate different tasks.

Finding the Nearest Neighbor of a Point

The naive way of finding the nearest neighbor of a point requires finding the distance to all elements. An orthtree can be used to perform the same task in significantly less time. For large numbers of elements, this can be a large enough difference to outweigh the time spent building the tree.

Note that a kd-tree is expected to outperform the orthtree for this task on points, it should be preferred unless features specific to the orthtree are needed.

The following example illustrates how to use an octree to accelerate the search for points close to a location.

Points are loaded from a file and an octree is built. The nearest neighbor method is invoked for several input points. A k value of 1 is used to find the single closest point. Results are put in a vector, and then printed.


File Orthtree/octree_find_nearest_neighbor.cpp

#include <fstream>
#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
#include <CGAL/Point_set_3.h>
#include <CGAL/Point_set_3/IO.h>
#include <boost/iterator/function_output_iterator.hpp>
// Type Declarations
using Point = Kernel::Point_3;
using Point_set = CGAL::Point_set_3<Point>;
using Point_map = Point_set::Point_map;
int main(int argc, char** argv) {
// Point set will be used to hold our points
Point_set points;
// Load points from a file.
std::ifstream stream((argc > 1) ? argv[1] : CGAL::data_file_path("points_3/cube.pwn"));
stream >> points;
if (0 == points.number_of_points()) {
std::cerr << "Error: cannot read file" << std::endl;
return EXIT_FAILURE;
}
std::cout << "loaded " << points.number_of_points() << " points" << std::endl;
// Create an octree from the points
Octree octree(points, points.point_map());
// Build the octree
octree.refine(10, 20);
// Find the nearest points to a few locations
std::vector<Point> points_to_find = {
{0, 0, 0},
{1, 1, 1},
{-1, -1, -1},
{-0.46026, -0.25353, 0.32051},
{-0.460261, -0.253533, 0.320513}
};
for (const Point& p : points_to_find)
octree.nearest_k_neighbors
(p, 1, // k=1 to find the single closest point
boost::make_function_output_iterator
([&](const Point_set::Index& nearest)
{
std::cout << "the nearest point to (" << p <<
") is (" << points.point(nearest) << ")" << std::endl;
}));
typename Octree::Sphere s(points_to_find[0], 0.0375);
std::cout << std::endl << "Closest points within the sphere around " << s.center() << " with squared radius of " << s.squared_radius() << ":" << std::endl;
octree.neighbors_within_radius
(s,
boost::make_function_output_iterator
([&](const Point_set::Index& nearest)
{
std::cout << points.point(nearest) << " dist: " << (Point(0, 0, 0) - points.point(nearest)).squared_length() << std::endl;
}));
std::size_t k = 3;
std::cout << std::endl << "The up to " << k << " closest points to(" << s.center() << ") within a squared radius of " << s.squared_radius() << " are:" << std::endl;
octree.nearest_k_neighbors_within_radius
(s, k,
boost::make_function_output_iterator
([&](const Point_set::Index& nearest)
{
std::cout << "(" << points.point(nearest) << " dist: " << (Point(0, 0, 0) - points.point(nearest)).squared_length() << std::endl;
}));
return EXIT_SUCCESS;
}
typename Traits::Sphere_d Sphere
Sphere type.
Definition: Orthtree.h:138

Not all octrees are compatible with nearest neighbor functionality, as the idea of a nearest neighbor may not make sense for some tree contents. For the nearest neighbor methods to work, the traits class must implement the CollectionPartitioningOrthtreeTraits concept.

Grading

An orthtree is considered "graded" if the difference of depth between two adjacent leaves is at most 1 for every pair of leaves.

Figure 94.3 Quadtree before and after being graded.


The following example demonstrates how to use the grade method to eliminate large jumps in depth within the orthtree.

A tree is created such that one node is split many more times than those it borders. grade() splits the octree's nodes so that adjacent nodes never have a difference in depth greater than one. The tree is printed before and after grading, so that the differences are visible.


File Orthtree/octree_grade.cpp

#include <iostream>
#include <CGAL/Simple_cartesian.h>
#include <CGAL/Octree.h>
// Type Declarations
using Point = Kernel::Point_3;
using Point_vector = std::vector<Point>;
int main() {
// Here, our point set is a vector
Point_vector points;
// Add a few points to the vector, most of which are in one region
points.emplace_back(1, 1, 1);
points.emplace_back(2, 1, -11);
points.emplace_back(2, 1, 1);
points.emplace_back(1, -2, 1);
points.emplace_back(1, 1, 1);
points.emplace_back(-1, 1, 1);
points.emplace_back(-1.1, 1, 1);
points.emplace_back(-1.01, 1, 1);
points.emplace_back(-1.001, 1, 1);
points.emplace_back(-1.0001, 1, 1);
points.emplace_back(-1.0001, 1, 1);
// Create an octree from the points
Octree octree(points);
// Build the octree with a small bucket size, so we get a deep node
octree.refine(10, 2);
// Print out the tree
std::cout << "\nUn-graded tree" << std::endl;
std::cout << "~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~" << std::endl;
std::cout << octree << std::endl;
// Grade the tree to eliminate large jumps in depth
octree.grade();
// Print out the tree again
std::cout << "\nGraded tree" << std::endl;
std::cout << "~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~" << std::endl;
std::cout << octree << std::endl;
return EXIT_SUCCESS;
}

Performance

Tree Construction

Tree construction benchmarks were conducted by randomly generating a collection of points, and then timing the process of creating a fully refined tree which contains them.

Because of its simplicity, an octree can be constructed faster than a kd-tree.

Figure 94.4 Plot of the time to construct a tree.


Nearest Neighbors

Orthtree nodes are uniform, so orthtrees will tend to have deeper hierarchies than equivalent kd-trees. As a result, orthtrees will generally perform worse for nearest neighbor searches. Both nearest neighbor algorithms have a theoretical complexity of \(O(log(n))\), but the orthtree can generally be expected to have a higher coefficient.

Figure 94.5 Plot of the time to find the 10 nearest neighbors of a random point using a pre-constructed tree.


The performance difference between the two trees is large, but both algorithms compare very favorably to the linear complexity of the naive approach, which involves comparing every point to the search point.

Using the orthtree for nearest neighbor computations instead of the kd-tree can be justified either when few queries are needed (as the construction is faster) or when the orthtree is also needed for other purposes.

Figure 94.6 Plot of the time to find nearest neighbors using tree methods and a naive approach.


For nontrivial point counts, the naive approach's calculation time dwarfs that of either the orthtree or kd-tree.

Migrating Code Written Before Release 6.0

The orthtree package changed to allow for custom data stored per node in the orthtree. To migrate existing code written before CGAL 6.0 Orthtree_traits_point can be used for a point-based orthtrees. The aliases CGAL::Quadtree and CGAL::Octree have been extended by a boolean template parameter to allow for non-cubic cells, which is the default. The data is passed via the traits class and no longer directly to the orthtree.

Former code to declare and define an Octree with cubic cells was as follows:

...
Octree octree(points, points.point_map());

CGAL 6.0 code with identical behavior is now:

...
Octree octree(points, points.point_map());

The node class does not exist anymore and has been replaced by the lightweight type Node_index. All information formerly contained in the node class is now accessible via the Orthtree interface.

Former node access was as follows:

Orthtree::Node root = orthtree.root();
Orthtree::Node child = root[0];
bool is_leaf = child.is_leaf();
for (Orthtree::Node::const_iterator it = child.begin(); it != child.end(); it++) {
const Orthtree::Point &p = *it;
...
}
typename Traits::Point_d Point
Point type.
Definition: Orthtree.h:136

CGAL 6.0 node access is now:

using Point = Kernel::Point_3;
using Point_set = CGAL::Point_set_3<Point>;
Point_set points;
...
Orthtree::Node_index root = orthtree.root();
Orthtree::Node_index child = orthtree.child(root, 0);
bool is_leaf = orthtree.is_leaf(child);
for (Octree::Traits::Node_data_element& e : octree.data(child)) {
// Using a pointmap is necessary when using a Point_set_3.
Point& p = points.point(e);
// If the orthtree is used with a std::vector<Point>, Node_data_element is identical to Point.
// Point& p = e;
...
}

The nearest neighbor search behaves as before, however, the output iterator will be required to store CollectionPartitioningOrthtreeTraits::Node_data_element. This may be the actual point or, e.g., in case a Point_set_3 is used, an index which requires the use of a point map to retrieve the point.

The provided traversals, i.e., CGAL::Orthtrees::Leaves_traversal, CGAL::Orthtrees::Preorder_traversal, CGAL::Orthtrees::Postorder_traversal, require the orthtree as template parameter now.

Former traversal use was as follows:

for (Orthtree::Node node : orthtree.traverse<CGAL::Orthtrees::Leaves_traversal>())
...
A class used for performing a traversal on leaves only.
Definition: Traversals.h:113

CGAL 6.0 traversal use is now:

History

A prototype code was implemented by Pierre Alliez and improved by Tong Zhao and Cédric Portaneri. From this prototype code, the package was developed by Jackson Campolattaro as part of the Google Summer of Code 2020. Simon Giraudot, supervisor of the GSoC internship, completed and finalized the package for integration in CGAL 5.3. Pierre Alliez provided kind help and advice all the way through. Starting with CGAL 6.0 an API improvement of the Orthtree package was released. Most notably, the orthtree nodes do not need to store anything. Data to be stored by the node are managed through a mechanism of dynamic properties. This improvement was done by Jackson Campolattaro thanks to a funding provided by INRIA, together with GeometryFactory.