libcuspatial  23.12.00
Loading...
Searching...
No Matches
Files | Classes | Functions
Spatial Indexing

APIs related to spatial indexing. More...

Files

file  point_quadtree.hpp
 

Classes

struct  cuspatial::point_quadtree
 
struct  cuspatial::point_quadtree_ref
 

Functions

template<class PointIterator , class T = typename cuspatial::iterator_value_type<PointIterator>>
std::pair< rmm::device_uvector< uint32_t >, point_quadtreecuspatial::quadtree_on_points (PointIterator points_first, PointIterator points_last, vec_2d< T > vertex_1, vec_2d< T > vertex_2, T scale, int8_t max_depth, int32_t max_size, rmm::cuda_stream_view stream=rmm::cuda_stream_default, rmm::mr::device_memory_resource *mr=rmm::mr::get_current_device_resource())
 Construct a quadtree structure from points.
 
std::pair< std::unique_ptr< cudf::column >, std::unique_ptr< cudf::table > > cuspatial::quadtree_on_points (cudf::column_view const &x, cudf::column_view const &y, double x_min, double x_max, double y_min, double y_max, double scale, int8_t max_depth, cudf::size_type max_size, rmm::mr::device_memory_resource *mr=rmm::mr::get_current_device_resource())
 Construct a quadtree structure from points.
 

Detailed Description

APIs related to spatial indexing.

Function Documentation

◆ quadtree_on_points() [1/2]

std::pair< std::unique_ptr< cudf::column >, std::unique_ptr< cudf::table > > cuspatial::quadtree_on_points ( cudf::column_view const &  x,
cudf::column_view const &  y,
double  x_min,
double  x_max,
double  y_min,
double  y_max,
double  scale,
int8_t  max_depth,
cudf::size_type  max_size,
rmm::mr::device_memory_resource *  mr = rmm::mr::get_current_device_resource() 
)

Construct a quadtree structure from points.

See also
http://www.adms-conf.org/2019-camera-ready/zhang_adms19.pdf for details.
Note
scale is applied to (x - x_min) and (y - y_min) to convert coordinates into a Morton code in 2D space.
max_depth should be less than 16, since Morton codes are represented as uint32_t. The eventual number of levels may be less than max_depth if the number of points is small or max_size is large.
All intermediate quadtree nodes will have fewer than max_size number of points. Leaf nodes are permitted (but not guaranteed) to have >= max_size number of points.
Parameters
xColumn of x-coordinates for each point.
yColumn of y-coordinates for each point.
x_minThe lower-left x-coordinate of the area of interest bounding box.
x_maxThe upper-right x-coordinate of the area of interest bounding box.
y_minThe lower-left y-coordinate of the area of interest bounding box.
y_maxThe upper-right y-coordinate of the area of interest bounding box.
scaleScale to apply to each x and y distance from x_min and y_min.
max_depthMaximum quadtree depth.
max_sizeMaximum number of points allowed in a node before it's split into 4 leaf nodes.
mrThe optional resource to use for output device memory allocations.
Exceptions
cuspatial::logic_errorIf the x and y column sizes are different
Returns
Pair of INT32 column of sorted keys to point indices, and cudf table with five columns for a complete quadtree: key - UINT32 column of quad node keys level - UINT8 column of quadtree levels is_internal_node - BOOL8 column indicating whether the node is a quad (true) or leaf (false) length - UINT32 column for the number of child nodes (if is_internal_node), or number of points offset - UINT32 column for the first child position (if is_internal_node), or first point position

◆ quadtree_on_points() [2/2]

template<class PointIterator , class T = typename cuspatial::iterator_value_type<PointIterator>>
std::pair< rmm::device_uvector< uint32_t >, point_quadtree > cuspatial::quadtree_on_points ( PointIterator  points_first,
PointIterator  points_last,
vec_2d< T >  vertex_1,
vec_2d< T >  vertex_2,
scale,
int8_t  max_depth,
int32_t  max_size,
rmm::cuda_stream_view  stream = rmm::cuda_stream_default,
rmm::mr::device_memory_resource *  mr = rmm::mr::get_current_device_resource() 
)

Construct a quadtree structure from points.

See also
http://www.adms-conf.org/2019-camera-ready/zhang_adms19.pdf for details.
Note
2D coordinates are converted into a 1D Morton code by dividing each x/y by the scale: ((x - min_x) / scale and (y - min_y) / scale).
max_depth should be less than 16, since Morton codes are represented as uint32_t. The eventual number of levels may be less than max_depth if the number of points is small or max_size is large.
All intermediate quadtree nodes will have fewer than max_size number of points. Leaf nodes are permitted (but not guaranteed) to have >= max_size number of points.
Template Parameters
PointIteratorIterator over x/y points. Must meet the requirements of LegacyRandomAccessIterator and be device-accessible.
Tthe floating-point coordinate value type of the input x/y points.
Parameters
points_firstIterator to the beginning of the range of (x, y) points.
points_lastIterator to the end of the range of (x, y) points.
vertex_1Vertex of the area of interest bounding box
vertex_2Vertex of the area of interest bounding box opposite vertex_1
scaleScale to apply to each x and y distance from min.x and min.y.
max_depthMaximum quadtree depth.
max_sizeMaximum number of points allowed in a node before it's split into 4 leaf nodes.
streamThe CUDA stream on which to perform computations
mrThe optional resource to use for output device memory allocations.
Returns
Pair of UINT32 column of sorted keys to point indices and a point_quadtree
Precondition
Point iterators must have the same vec_2d value type, with the same underlying floating-point coordinate type (e.g. cuspatial::vec_2d<float>).