Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Typedefs | Functions
rtree.h File Reference
#include <deal.II/base/config.h>
#include <deal.II/base/point.h>
#include <deal.II/base/std_cxx20/iota_view.h>
#include <deal.II/boost_adaptors/bounding_box.h>
#include <deal.II/boost_adaptors/point.h>
#include <deal.II/boost_adaptors/segment.h>
#include <boost/geometry/algorithms/distance.hpp>
#include <boost/geometry/index/rtree.hpp>
#include <boost/geometry/strategies/strategies.hpp>
#include <memory>

Go to the source code of this file.

Classes

class  IndexableGetterFromIndices< Container >
 
struct  ExtractLevelVisitor< Value, Options, Translator, Box, Allocators >
 

Typedefs

template<typename LeafType , typename IndexType = boost::geometry::index::linear<16>, typename IndexableGetter = boost::geometry::index::indexable<LeafType>>
using RTree = boost::geometry::index::rtree< LeafType, IndexType, IndexableGetter >
 

Functions

template<typename IndexType = boost::geometry::index::linear<16>, typename LeafTypeIterator , typename IndexableGetter = boost::geometry::index::indexable< typename LeafTypeIterator::value_type>>
RTree< typename LeafTypeIterator::value_type, IndexType, IndexableGetter > pack_rtree (const LeafTypeIterator &begin, const LeafTypeIterator &end)
 
template<typename IndexType = boost::geometry::index::linear<16>, typename ContainerType , typename IndexableGetter = boost::geometry::index::indexable< typename ContainerType::value_type>>
RTree< typename ContainerType::value_type, IndexType, IndexableGetter > pack_rtree (const ContainerType &container)
 
template<typename IndexType = boost::geometry::index::linear<16>, typename ContainerType >
RTree< typename ContainerType::size_type, IndexType, IndexableGetterFromIndices< ContainerType > > pack_rtree_of_indices (const ContainerType &container)
 
template<typename Rtree >
std::vector< BoundingBox< boost::geometry::dimension< typename Rtree::indexable_type >::value > > extract_rtree_level (const Rtree &tree, const unsigned int level)
 

Typedef Documentation

◆ RTree

template<typename LeafType , typename IndexType = boost::geometry::index::linear<16>, typename IndexableGetter = boost::geometry::index::indexable<LeafType>>
using RTree = boost::geometry::index::rtree<LeafType, IndexType, IndexableGetter>

A wrapper for the boost::geometry::index::rtree class, implementing a self-balancing spatial index (the R-tree) capable of storing various types of values, using different balancing algorithms.

From Wikipedia:

R-trees are tree data structures used for spatial access methods, i.e., for indexing multi-dimensional information such as geographical coordinates, rectangles or polygons. The R-tree was proposed by Antonin Guttman in 1984 and has found significant use in both theoretical and applied contexts. A common real-world usage for an R-tree might be to store spatial objects such as restaurant locations or the polygons that typical maps are made of: streets, buildings, outlines of lakes, coastlines, etc. and then find answers quickly to queries such as "Find all museums within 2 km of my current location", "retrieve all road segments within 2 km of my location" (to display them in a navigation system) or "find the nearest gas station" (although not taking roads into account). The R-tree can also accelerate nearest neighbor search for various distance metrics, including great-circle distance.

The key idea of the data structure is to group nearby objects and represent them with their minimum bounding rectangle in the next higher level of the tree; the "R" in R-tree is for rectangle. Since all objects lie within this bounding rectangle, a query that does not intersect the bounding rectangle also cannot intersect any of the contained objects. At the leaf level, each rectangle describes a single object; at higher levels the aggregation of an increasing number of objects. This can also be seen as an increasingly coarse approximation of the data set.

The key difficulty of R-tree is to build an efficient tree that on one hand is balanced (so the leaf nodes are at the same height) on the other hand the rectangles do not cover too much empty space and do not overlap too much (so that during search, fewer subtrees need to be processed). For example, the original idea for inserting elements to obtain an efficient tree is to always insert into the subtree that requires least enlargement of its bounding box. Once that page is full, the data is split into two sets that should cover the minimal area each. Most of the research and improvements for R-trees aims at improving the way the tree is built and can be grouped into two objectives: building an efficient tree from scratch (known as bulk-loading) and performing changes on an existing tree (insertion and deletion).

An RTree may store any type of LeafType as long as it is possible to extract an Indexable that the RTree can handle and compare values. An Indexable is a type adapted to the Point, BoundingBox or Segment concept, for which distance and equality comparison are implemented. The deal.II Point, Segment, and BoundingBox classes satisfy this requirement, but you can mix in any geometry object that boost::geometry accepts as indexable.

In particular, given an Indexable type (for example a Point, a BoundingBox, or a Segment), LeafType can by any of Indexable, std::pair<Indexable, T>, boost::tuple<Indexable, ...> or std::tuple<Indexable, ...>.

The optional argument IndexType is used only when adding elements to the tree one by one. If a range insertion is used, then the tree is built using the packing algorithm.

Linear, quadratic, and rstar algorithms are available if one wants to construct the tree sequentially. However, none of these is very efficient, and users should use the packing algorithm when possible.

The packing algorithm constructs the tree all at once, and may be used when you have all the leaves at your disposal.

This class is usually used in combination with one of the two helper functions pack_rtree(), that takes a container or a range of iterators to construct the RTree using the packing algorithm.

An example usage is the following:

std::vector<Point<2>> points = generate_some_points();
auto tree = pack_rtree(points.begin(), points.end());
// or, equivalently:
// auto tree = pack_rtree(points);
RTree< typename LeafTypeIterator::value_type, IndexType, IndexableGetter > pack_rtree(const LeafTypeIterator &begin, const LeafTypeIterator &end)

The tree is accessed by using boost::geometry::index queries. For example, after constructing the tree with the snippet above, one can ask for the closest points to a segment in the following way:

namespace bgi = boost::geometry::index;
Segment<2> segment(Point<2>(0,0), Point<2>(1,1));
std::vector<Point<2>> intersection;
tree.query(bgi::nearest(segment,3), std::back_inserter(intersection));
// Returns the 3 closest points to the Segment defined above.
Definition point.h:112
boost::geometry::model::segment< Point< dim > > Segment
Definition segment.h:34

In general, a tree does not need to store the actual objects, as long as it knows how to access a const reference to an indexable type. This is achieved by passing the optional template argument IndexableGetter, that extracts a const reference to one of the possible indexable types given an object of type LeafType. As an example, one may store points in a container, and only create a tree of the indices within the container, using the IndexableGetterFromIndices class defined below, and the function pack_rtree_of_indices().

Definition at line 143 of file rtree.h.

Function Documentation

◆ pack_rtree() [1/2]

template<typename IndexType = boost::geometry::index::linear<16>, typename LeafTypeIterator , typename IndexableGetter = boost::geometry::index::indexable< typename LeafTypeIterator::value_type>>
RTree< typename LeafTypeIterator::value_type, IndexType, IndexableGetter > pack_rtree ( const LeafTypeIterator &  begin,
const LeafTypeIterator &  end 
)

Construct the correct RTree object by passing an iterator range.

Notice that the order of the parameters is the opposite with respect to the RTree class, since we can automatically infer the LeafType from the arguments.

◆ pack_rtree() [2/2]

template<typename IndexType = boost::geometry::index::linear<16>, typename ContainerType , typename IndexableGetter = boost::geometry::index::indexable< typename ContainerType::value_type>>
RTree< typename ContainerType::value_type, IndexType, IndexableGetter > pack_rtree ( const ContainerType &  container)

Construct an RTree object by passing an STL container type. This function is used in step-70.

Notice that the order of the template parameters is the opposite with respect to the RTree class, since we can automatically infer the LeafType from the arguments, and we only need to specify the IndexType if the default is not adequate.

◆ pack_rtree_of_indices()

template<typename IndexType = boost::geometry::index::linear<16>, typename ContainerType >
RTree< typename ContainerType::size_type, IndexType, IndexableGetterFromIndices< ContainerType > > pack_rtree_of_indices ( const ContainerType &  container)

Construct a RTree object that stores the indices of an existing container of indexable types. The only requirement on the container is that it supports operator[] for any index between 0 and the size of the container (i.e., a std::vector, or an std::array will do, however an std::map won't).

Differently from the object created by the pack_rtree() function, in this case we don't store the actual geometrical types, but just their indices:

namespace bgi = boost::geometry::index;
std::vector<Point<dim>> some_points = fill();
auto tree = pack_rtree(points); // the tree contains a copy of the points
auto index_tree = pack_rtree_of_indices(points); // the tree contains only
// the indices of the
// points
BoundingBox<dim> box = build_a_box();
for(const auto &p: tree | bgi::adaptors::queried(bgi::intersects(box)))
std::cout << "Point p: " << p << std::endl;
for(const auto &i: index_tree | bgi::adaptors::queried(bgi::intersects(box)))
std::cout << "Point p: " << some_points[i] << std::endl;
STL namespace.
RTree< typename ContainerType::size_type, IndexType, IndexableGetterFromIndices< ContainerType > > pack_rtree_of_indices(const ContainerType &container)

The leaves stored in the tree are the indices of the corresponding entry in the container. A reference to the external container is stored internally, but keep in mind that if you change the container, you should rebuild the tree.

◆ extract_rtree_level()

template<typename Rtree >
std::vector< BoundingBox< boost::geometry::dimension< typename Rtree::indexable_type >::value > > extract_rtree_level ( const Rtree &  tree,
const unsigned int  level 
)
inline

Given a RTree object rtree, and a target level level, return a vector of BoundingBox objects containing all the bounding boxes that make the given level of the rtree. This function is a convenient wrapper around the ExtractLevelVisitor class. It is used in step-70.

Since an RTree object is a balanced tree, you can expect each entry of the resulting vector to contain roughly the same number of children, and ultimately, the same number of leaf objects. If you request for a level that is not present in the RTree, the last level is returned.

A typical usage of this function is in the context of parallel::distributed::Triangulation objects, where one would like to construct a rough representation of the area which is covered by the locally owned cells of the active process, and exchange this information with other processes. The finest level of information is given by the leaves, which in this context would be the collection of all the bounding boxes associated to the locally owned cells of the triangulation. Exchanging this information with all participating processes would defeat the purpuse of parallel computations. If however one constructs an RTree containing these bounding boxes (for example, by calling GridTools::Cache::get_cell_bounding_boxes_rtree()), and then extracts one of the first levels of the RTree, only a handful of BoundingBox objects would be returned, allowing the user to have a very efficient description of the geometry of the domain, and of its distribution among processes.

An example usage is given by the following snippet of code:

std::vector<BoundingBox<2>> all_boxes(tria.n_locally_owned_active_cells());
unsigned int i = 0;
for (const auto &cell : tria.active_cell_iterators())
if (cell->is_locally_owned())
all_boxes[i++] = cell->bounding_box();
const auto tree = pack_rtree(all_boxes);
const auto boxes = extract_rtree_level(tree, 1);
void refine_global(const unsigned int times=1)
void hyper_ball(Triangulation< dim > &tria, const Point< dim > &center=Point< dim >(), const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false)
std::vector< BoundingBox< boost::geometry::dimension< typename Rtree::indexable_type >::value > > extract_rtree_level(const Rtree &tree, const unsigned int level)
const ::Triangulation< dim, spacedim > & tria

When run on three processes, the complete set of the BoundingBox objects surrounding only the locally owned cells and the second level of the rtree constructed with those boxes would look like in the following pictures (one image per process):