Octree
The Octree class implements the generalized octree, a hierarchical tree
structure.  When built on data in two dimensions, it is also called a
‘quadtree’.  The generalized octree is most useful on data in only two or three
dimensions, as its number of children is exponential in the data dimension:
so, e.g., an octree node in two dimensions has up to four children; in three
dimensions has up to eight; and so on.
mlpack’s Octree implementation differs from many textbook descriptions of
quadtrees and octrees in that mlpack allows the bounding box for an Octree
node to be a general hyperrectangle instead of a hypercube (e.g. each dimension
may have a different width).  This allows the Octree to fit the data points
with tighter bounds for each node—very important for pruning in distance-based
tasks.
For higher-dimensional use cases (e.g. data dimension of 5 or more), the
KDTree, other BinarySpaceTree variants,
or virtually any other type of tree will perform better.
mlpack’s Octree implementation supports three template parameters for
configurable behavior, and implements all the functionality required by the
TreeType API, plus some
additional functionality specific to octrees.
- Template parameters
 - Constructors
 - Basic tree properties
 - Bounding distances with the tree
 - Tree traversals
 - Example usage
 
🔗 See also
- Octree on Wikipedia
 - Quadtree on Wikipedia
 - Binary space partitioning on Wikipedia
 - Tree-Independent Dual-Tree Algorithms (pdf)
 
🔗 Template parameters
In accordance with the TreeType
API
(see also this more detailed section),
the Octree class takes three template parameters:
Octree<DistanceType, StatisticType, MatType>
DistanceType: the distance metric to use for distance computations. For theOctree, this must be anLMetric. By default, this isEuclideanDistance.StatisticType: this holds auxiliary information in each tree node. By default,EmptyStatisticis used, which holds no information.MatType: the type of matrix used to represent points. Must be a type matching the Armadillo API. By default,arma::matis used, but other types such asarma::fmator similar will work just fine.
If no template parameters are explicitly specified, then defaults are used:
Octree<> = Octree<EuclideanDistance, EmptyStatistic, arma::mat>
🔗 Constructors
Octrees are efficiently constructed by permuting points in a dataset in a
quicksort-like algorithm.  However, this means that the ordering of points in
the tree’s dataset (accessed with node.Dataset()) after construction may be
different.
node = Octree(data, maxLeafSize=20)node = Octree(data, oldFromNew, maxLeafSize=20)node = Octree(data, oldFromNew, newFromOld, maxLeafSize=20)- Construct an 
Octreeon the givendata, usingmaxLeafSizeas the maximum number of points held in a leaf. - By default, 
datais copied. Avoid a copy by usingstd::move()(e.g.std::move(data)); when doing this,datawill be set to an empty matrix. - Optionally, construct mappings from old points to new points.  
oldFromNewandnewFromOldwill have lengthdata.n_cols, and:oldFromNew[i]indicates that pointiin the tree’s dataset was originally pointoldFromNew[i]indata; that is,node.Dataset().col(i)is the pointdata.col(oldFromNew[i]).newFromOld[i]indicates that pointiindatais now pointnewFromOld[i]in the tree’s dataset; that is,node.Dataset().col(newFromOld[i])is the pointdata.col(i).
 
- Construct an 
 
node = Octree<DistanceType, StatisticType, MatType>(data, maxLeafSize=20)node = Octree<DistanceType, StatisticType, MatType>(data, oldFromNew, maxLeafSize=20)node = Octree<DistanceType, StatisticType, MatType>(data, oldFromNew, newFromOld, maxLeafSize=20)- Construct an 
Octreeon the givendata, using custom template parameters to control the behavior of the tree, usingmaxLeafSizeas the maximum number of points held in a leaf. - By default, 
datais copied. Avoid a copy by usingstd::move()(e.g.std::move(data)); when doing this,datawill be set to an empty matrix. - Optionally, construct mappings from old points to new points.  
oldFromNewandnewFromOldwill have lengthdata.n_cols, and:oldFromNew[i]indicates that pointiin the tree’s dataset was originally pointoldFromNew[i]indata; that is,node.Dataset().col(i)is the pointdata.col(oldFromNew[i]).newFromOld[i]indicates that pointiindatais now pointnewFromOld[i]in the tree’s dataset; that is,node.Dataset().col(newFromOld[i])is the pointdata.col(i).
 
- Construct an 
 
node = Octree()- Construct an empty octree with no children and no points.
 
Notes:
- 
    
The name
nodeis used here forOctreeobjects instead oftree, because eachOctreeobject is a single node in the tree. The constructor returns the node that is the root of the tree. - 
    
Inserting individual points or removing individual points from an
Octreeis not supported, because this generally results in a octree with very loose bounding boxes. It is better to simply build a newOctreeon the modified dataset. For trees that support individual insertion and deletions, see theRectangleTreeclass and all its variants (e.g.RTree,RStarTree, etc.). - 
    
See also the developer documentation on tree constructors.
 
🔗 Constructor parameters:
| name | type | description | default | 
|---|---|---|---|
data | 
      arma::mat | 
      Column-major matrix to build the tree on.  Pass with std::move(data) to avoid copying the matrix. | 
      (N/A) | 
maxLeafSize | 
      size_t | 
      Maximum number of points to store in each leaf. | 20 | 
    
oldFromNew | 
      std::vector<size_t> | 
      Mappings from points in node.Dataset() to points in data. | 
      (N/A) | 
newFromOld | 
      std::vector<size_t> | 
      Mappings from points in data to points in node.Dataset(). | 
      (N/A) | 
🔗 Basic tree properties
Once an Octree object is constructed, various properties of the tree can be
accessed or inspected.  Many of these functions are required by the TreeType
API.
🔗 Navigating the tree
- 
    
node.NumChildren()returns the number of children innode. This is0ifnodeis a leaf, and anywhere between1andpow(2, node.Dataset().n_rows)otherwise. - 
    
node.IsLeaf()returns aboolindicating whether or notnodeis a leaf. node.Child(i)returns anOctree&that is theith child.imust be less thannode.NumChildren().- This function should only be called if 
node.NumChildren()is not0(e.g. ifnodeis not a leaf). Note that this returns a validOctree&that can itself be used just like the root node of the tree! 
node.Parent()will return anOctree*that points to the parent ofnode, orNULLifnodeis the root of theOctree.
🔗 Accessing members of a tree
node.Bound()will return anHRectBound&object that represents the hyperrectangle bounding box ofnode. This is the smallest hyperrectangle that encloses all the descendant points ofnode.- Note: mlpack’s 
Octreeimplementation does not require that the bounding structure is a hypercube—each dimension of the bound is allowed to have a different width. 
- Note: mlpack’s 
 - 
    
node.Stat()will return anEmptyStatistic&(or aStatisticType&if a customStatisticTypewas specified as a template parameter) holding the statistics of the node that were computed during tree construction. node.Distance()will return aEuclideanDistance&(or aDistanceType&if a customDistanceTypewas specified as a template parameter).- This function is required by the
TreeType API, but given
that 
Octreerequires anLMetricto be used, andLMetriconly hasstaticfunctions and holds no state, this function is not likely to be useful. 
- This function is required by the
TreeType API, but given
that 
 
See also the developer documentation for basic tree functionality in mlpack.
🔗 Accessing data held in a tree
node.Dataset()will return aconst arma::mat&that is the dataset the tree was built on. Note that this is a permuted version of thedatamatrix passed to the constructor.- If a custom 
MatTypeis being used, the return type will beconst MatType&instead ofconst arma::mat&. 
- If a custom 
 node.NumPoints()returns asize_tindicating the number of points held directly innode.- If 
nodeis not a leaf, this will return0, asOctreeonly holds points directly in its leaves. - If 
nodeis a leaf, then the number of points will be less than or equal to themaxLeafSizethat was specified when the tree was constructed. 
- If 
 node.Point(i)returns asize_tindicating the index of thei‘th point innode.Dataset().imust be in the range[0, node.NumPoints() - 1](inclusive).nodemust be a leaf (as non-leaves do not hold any points).- The 
i‘th point innodecan then be accessed asnode.Dataset().col(node.Point(i)). - In an 
Octree, because of the permutation of points done during construction, point indices are contiguous:node.Point(i + j)is the same asnode.Point(i) + jfor validiandj. - Accessing the actual 
i‘th point itself can be done with, e.g.,node.Dataset().col(node.Point(i)). 
node.NumDescendants()returns asize_tindicating the number of points held in all descendant leaves ofnode.- If 
nodeis the root of the tree, thennode.NumDescendants()will be equal tonode.Dataset().n_cols. 
- If 
 node.Descendant(i)returns asize_tindicating the index of thei‘th descendant point innode.Dataset().imust be in the range[0, node.NumDescendants() - 1](inclusive).nodedoes not need to be a leaf.- The 
i‘th descendant point innodecan then be accessed asnode.Dataset().col(node.Descendant(i)). - In an 
Octree, because of the permutation of points done during construction, point indices are contiguous:node.Descendant(i + j)is the same asnode.Descendant(i) + jfor validiandj. - Accessing the actual 
i‘th descendant itself can be done with, e.g.,node.Dataset().col(node.Descendant(i)). 
node.Begin()returns asize_tindicating the index of the first descendant point ofnode.- This is equivalent to 
node.Descendant(0). 
- This is equivalent to 
 node.Count()returns asize_tindicating the number of descendant points ofnode.- This is equivalent to 
node.NumDescendants(). 
- This is equivalent to 
 
🔗 Accessing computed bound quantities of a tree
The following quantities are cached for each node in an Octree, and so
accessing them does not require any computation.
node.FurthestPointDistance()returns adoublerepresenting the distance between the center of the bounding hyperrectangle ofnodeand the furthest point held bynode.- If 
nodeis not a leaf, this returns 0 (becausenodedoes not hold any points). 
- If 
 - 
    
node.FurthestDescendantDistance()returns adoublerepresenting the distance between the center of the bounding hyperrectangle ofnodeand the furthest descendant point held bynode. node.MinimumBoundDistance()returns adoublerepresenting the minimum possible distance from the center of the node to any edge of the hyperrectangle bound.- This quantity is half the width of the smallest dimension of
node.Bound(). 
- This quantity is half the width of the smallest dimension of
 node.ParentDistance()returns adoublerepresenting the distance between the center of the bounding hyperrectangle ofnodeand the center of the bounding hyperrectangle of its parent.- If 
nodeis the root of the tree,0is returned. 
- If 
 
Notes:
- 
    
If a custom
MatTypewas specified when constructing theOctree, then the return type of each method is the element type of the givenMatTypeinstead ofdouble. (e.g., ifMatTypeisarma::fmat, then the return type isfloat.) - 
    
For more details on each bound quantity, see the developer documentation on bound quantities for trees.
 
🔗 Other functionality
node.Center(center)computes the center of the bounding hyperrectangle ofnodeand stores it incenter.centershould be of typearma::vec&. (If a customMatTypewas specified when constructing theOctree, the type is instead the column vector type for the givenMatType; e.g.,arma::fvec&whenMatTypeisarma::fmat.)centerwill be set to have size equivalent to the dimensionality of the dataset held bynode.- This is equivalent to calling 
node.Bound().Center(center). 
- An 
Octreecan be serialized withdata::Save()anddata::Load(). 
🔗 Bounding distances with the tree
The primary use of trees in mlpack is bounding distances to points or other tree nodes. The following functions can be used for these tasks.
node.GetNearestChild(point)node.GetFurthestChild(point)- Return a 
size_tindicating the index of the child (0for left,1for right) that is closest to (or furthest from)point, with respect to theMinDistance()(orMaxDistance()) function. - If there is a tie, the child with the lowest index is returned.
 - If 
nodeis a leaf,0is returned. pointshould be of typearma::vec. (If a customMatTypewas specified when constructing theOctree, the type is instead the column vector type for the givenMatType; e.g.,arma::fvecwhenMatTypeisarma::fmat.)
- Return a 
 node.GetNearestChild(other)node.GetFurthestChild(other)- Return a 
size_tindicating the index of the child (0for left,1for right) that is closest to (or furthest from) theOctreenodeother, with respect to theMinDistance()(orMaxDistance()) function. - If there is a tie, the child with the lowest index is returned.
 - If 
nodeis a leaf,0is returned. 
- Return a 
 
node.MinDistance(point)node.MinDistance(other)- Return a 
doubleindicating the minimum possible distance betweennodeandpoint, or theOctreenodeother. - This is equivalent to the minimum possible distance between any point
contained in the bounding hyperrectangle of 
nodeandpoint, or between any point contained in the bounding hyperrectangle ofnodeand any point contained in the bounding hyperrectangle ofother. pointshould be of typearma::vec. (If a customMatTypewas specified when constructing theOctree, the type is instead the column vector type for the givenMatType, and the return type is the element type ofMatType; e.g.,pointshould bearma::fvecwhenMatTypeisarma::fmat, and the returned distance isfloat).
- Return a 
 node.MaxDistance(point)node.MaxDistance(other)- Return a 
doubleindicating the maximum possible distance betweennodeandpoint, or theOctreenodeother. - This is equivalent to the maximum possible distance between any point
contained in the bounding hyperrectangle of 
nodeandpoint, or between any point contained in the bounding hyperrectangle ofnodeand any point contained in the bounding hyperrectangle ofother. pointshould be of typearma::vec. (If a customMatTypewas specified when constructing theOctree, the type is instead the column vector type for the givenMatType, and the return type is the element type ofMatType; e.g.,pointshould bearma::fvecwhenMatTypeisarma::fmat, and the returned distance isfloat).
- Return a 
 node.RangeDistance(point)node.RangeDistance(other)- Return a 
Rangewhose lower bound isnode.MinDistance(point)ornode.MinDistance(other), and whose upper bound isnode.MaxDistance(point)ornode.MaxDistance(other). pointshould be of typearma::vec. (If a customMatTypewas specified when constructing theOctree, the type is instead the column vector type for the givenMatType, and the return type is aRangeTypewith element type the same asMatType; e.g.,pointshould bearma::fvecwhenMatTypeisarma::fmat, and the returned type isRangeType<float>).
- Return a 
 
🔗 Tree traversals
Like every mlpack tree, the Octree class provides a single-tree and dual-tree
traversal that can be paired with a
RuleType class to implement a single-tree
or dual-tree algorithm.
Octree::SingleTreeTraverser- Implements a depth-first single-tree traverser.
 
Octree::DualTreeTraverser- Implements a dual-depth-first dual-tree traverser.
 
🔗 Example usage
Build an Octree on the iris dataset and print basic statistics about the
tree.
// See https://datasets.mlpack.org/iris.csv.
arma::mat dataset;
mlpack::data::Load("iris.csv", dataset, true);
// Build the octree with a leaf size of 10.  (This means that nodes are split
// until they contain 10 or fewer points.)
//
// The std::move() means that `dataset` will be empty after this call, and no
// data will be copied during tree building.
//
// Note that the '<>' isn't necessary if C++20 is being used (e.g.
// `mlpack::Octree tree(...)` will work fine in C++20 or newer).
mlpack::Octree<> tree(std::move(dataset));
// Print the bounding box of the root node.
std::cout << "Bounding box of root node:" << std::endl;
for (size_t i = 0; i < tree.Bound().Dim(); ++i)
{
  std::cout << " - Dimension " << i << ": [" << tree.Bound()[i].Lo() << ", "
      << tree.Bound()[i].Hi() << "]." << std::endl;
}
std::cout << std::endl;
// Print the number of descendant points of the root, and of each of its
// children.
std::cout << "Descendant points of root:        "
    << tree.NumDescendants() << "." << std::endl;
for (size_t i = 0; i < tree.NumChildren(); ++i)
{
  std::cout << "Descendant points of child " << i << ": "
      << tree.Child(i).NumDescendants() << "." << std::endl;
}
std::cout << std::endl;
// Compute the center of the octree.
arma::vec center;
tree.Center(center);
std::cout << "Center of octree: " << center.t();
Build two Octrees on subsets of the tiny LCDM dataset (3-dimensional
astronomical data) and compute minimum and maximum distances between different
nodes in the tree.
// See https://datasets.mlpack.org/lcdm_tiny.csv.
arma::mat dataset;
mlpack::data::Load("lcdm_tiny.csv", dataset, true);
// Build octrees on the first half and the second half of points.
mlpack::Octree<> tree1(dataset.cols(0, dataset.n_cols / 2));
mlpack::Octree<> tree2(dataset.cols(dataset.n_cols / 2 + 1,
    dataset.n_cols - 1));
// Compute the maximum distance between the trees.
std::cout << "Maximum distance between tree root nodes: "
    << tree1.MaxDistance(tree2) << "." << std::endl;
// Get the leftmost grandchild of the first tree's root---if it exists.
if (!tree1.IsLeaf() && !tree1.Child(0).IsLeaf())
{
  mlpack::Octree<>& node1 = tree1.Child(0).Child(0);
  // Get the rightmost grandchild of the second tree's root---if it exists.
  if (!tree2.IsLeaf() && !tree2.Child(0).IsLeaf())
  {
    mlpack::Octree<>& node2 = tree2.Child(0).Child(0);
    // Print the minimum and maximum distance between the nodes.
    mlpack::Range dists = node1.RangeDistance(node2);
    std::cout << "Possible distances between two grandchild nodes: ["
        << dists.Lo() << ", " << dists.Hi() << "]." << std::endl;
    // Print the minimum distance between the first node and the first
    // descendant point of the second node.
    const size_t descendantIndex = node2.Descendant(0);
    const double descendantMinDist =
        node1.MinDistance(node2.Dataset().col(descendantIndex));
    std::cout << "Minimum distance between grandchild node and descendant "
        << "point: " << descendantMinDist << "." << std::endl;
    // Which child of node2 is closer to node1?
    const size_t closerIndex = node2.GetNearestChild(node1);
    if (closerIndex == 0)
      std::cout << "The left child of node2 is closer to node1." << std::endl;
    else if (closerIndex == 1)
      std::cout << "The right child of node2 is closer to node1." << std::endl;
    else // closerIndex == 2 in this case.
      std::cout << "Both children of node2 are equally close to node1."
          << std::endl;
    // And which child of node1 is further from node2?
    const size_t furtherIndex = node1.GetFurthestChild(node2);
    if (furtherIndex == 0)
      std::cout << "The left child of node1 is further from node2."
          << std::endl;
    else if (furtherIndex == 1)
      std::cout << "The right child of node1 is further from node2."
          << std::endl;
    else // furtherIndex == 2 in this case.
      std::cout << "Both children of node1 are equally far from node2."
          << std::endl;
  }
}
Build an Octree on 32-bit floating point data and save it to disk.
// See https://datasets.mlpack.org/iris.csv.
arma::fmat dataset;
mlpack::data::Load("iris.csv", dataset);
// Build the Octree using 32-bit floating point data as the matrix type.
// We will still use the default EmptyStatistic and EuclideanDistance
// parameters.  A leaf size of 100 is used here.
mlpack::Octree<mlpack::EuclideanDistance,
               mlpack::EmptyStatistic,
               arma::fmat> tree(std::move(dataset), 100);
// Save the Octree to disk with the name 'tree'.
mlpack::data::Save("tree.bin", "tree", tree);
std::cout << "Saved tree with " << tree.Dataset().n_cols << " points to "
    << "'tree.bin'." << std::endl;
Load a 32-bit floating point Octree from disk, then traverse it manually and
find the number of leaf nodes with fewer than 10 points.
// This assumes the tree has already been saved to 'tree.bin' (as in the example
// above).
// This convenient typedef saves us a long type name!
using TreeType = mlpack::Octree<mlpack::EuclideanDistance,
                                mlpack::EmptyStatistic,
                                arma::fmat>;
TreeType tree;
mlpack::data::Load("tree.bin", "tree", tree);
std::cout << "Tree loaded with " << tree.NumDescendants() << " points."
    << std::endl;
// Recurse in a depth-first manner.  Count both the total number of leaves, and
// the number of leaves with fewer than 10 points.
size_t leafCount = 0;
size_t totalLeafCount = 0;
std::stack<TreeType*> stack;
stack.push(&tree);
while (!stack.empty())
{
  TreeType* node = stack.top();
  stack.pop();
  if (node->NumPoints() < 10)
    ++leafCount;
  ++totalLeafCount;
  if (!node->IsLeaf())
  {
    for (size_t i = 0; i < node->NumChildren(); ++i)
      stack.push(&node->Child(i));
  }
}
// Note that it would be possible to use TreeType::SingleTreeTraverser to
// perform the recursion above, but that is more well-suited for more complex
// tasks that require pruning and other non-trivial behavior; so using a simple
// stack is the better option here.
// Print the results.
std::cout << leafCount << " out of " << totalLeafCount << " leaves have fewer "
  << "than 10 points." << std::endl;
Build an Octree and map between original points and new points.
// See https://datasets.mlpack.org/lcdm_tiny.csv.
arma::mat dataset;
mlpack::data::Load("lcdm_tiny.csv", dataset, true);
// Build the tree.
std::vector<size_t> oldFromNew, newFromOld;
mlpack::Octree<> tree(dataset, oldFromNew, newFromOld);
// oldFromNew and newFromOld will be set to the same size as the dataset.
std::cout << "Number of points in dataset: " << dataset.n_cols << "."
    << std::endl;
std::cout << "Size of oldFromNew: " << oldFromNew.size() << "." << std::endl;
std::cout << "Size of newFromOld: " << newFromOld.size() << "." << std::endl;
std::cout << std::endl;
// See where point 42 in the tree's dataset came from.
std::cout << "Point 42 in the permuted tree's dataset:" << std::endl;
std::cout << "  " << tree.Dataset().col(42).t();
std::cout << "Was originally point " << oldFromNew[42] << ":" << std::endl;
std::cout << "  " << dataset.col(oldFromNew[42]).t();
std::cout << std::endl;
// See where point 7 in the original dataset was mapped.
std::cout << "Point 7 in original dataset:" << std::endl;
std::cout << "  " << dataset.col(7).t();
std::cout << "Mapped to point " << newFromOld[7] << ":" << std::endl;
std::cout << "  " << tree.Dataset().col(newFromOld[7]).t();