Loading...
Searching...
No Matches
Files | Functions
Spatial Join

APIs related to spatial join. More...

Files

file  spatial_join.hpp
 

Functions

template<class InputIt , class OutputIt , class T >
OutputIt cuspatial::sinusoidal_projection (InputIt lon_lat_first, InputIt lon_lat_last, OutputIt xy_first, vec_2d< T > origin, rmm::cuda_stream_view stream=rmm::cuda_stream_default)
 Sinusoidal projection of longitude/latitude relative to origin to Cartesian (x/y) coordinates in km.
 
template<class BoundingBoxIterator , class T = typename cuspatial::iterator_vec_base_type<BoundingBoxIterator>>
std::pair< rmm::device_uvector< uint32_t >, rmm::device_uvector< uint32_t > > cuspatial::join_quadtree_and_bounding_boxes (point_quadtree_ref quadtree, BoundingBoxIterator bounding_boxes_first, BoundingBoxIterator bounding_boxes_last, vec_2d< T > const &v_min, T scale, int8_t max_depth, rmm::cuda_stream_view stream=rmm::cuda_stream_default, rmm::device_async_resource_ref mr=rmm::mr::get_current_device_resource())
 Search a quadtree for polygon or linestring bounding box intersections.
 
template<class PolyIndexIterator , class QuadIndexIterator , class PointIndexIterator , class PointIterator , class MultiPolygonRange , class IndexType = iterator_value_type<PointIndexIterator>>
std::pair< rmm::device_uvector< IndexType >, rmm::device_uvector< IndexType > > cuspatial::quadtree_point_in_polygon (PolyIndexIterator poly_indices_first, PolyIndexIterator poly_indices_last, QuadIndexIterator quad_indices_first, point_quadtree_ref quadtree, PointIndexIterator point_indices_first, PointIndexIterator point_indices_last, PointIterator points_first, MultiPolygonRange polygons, rmm::cuda_stream_view stream=rmm::cuda_stream_default, rmm::device_async_resource_ref mr=rmm::mr::get_current_device_resource())
 Test whether the specified points are inside any of the specified polygons.
 
template<class LinestringIndexIterator , class QuadIndexIterator , class PointIndexIterator , class PointIterator , class MultiLinestringRange , typename IndexType = iterator_value_type<PointIndexIterator>, typename T = iterator_vec_base_type<PointIterator>>
std::tuple< rmm::device_uvector< IndexType >, rmm::device_uvector< IndexType >, rmm::device_uvector< T > > cuspatial::quadtree_point_to_nearest_linestring (LinestringIndexIterator linestring_indices_first, LinestringIndexIterator linestring_indices_last, QuadIndexIterator quad_indices_first, point_quadtree_ref quadtree, PointIndexIterator point_indices_first, PointIndexIterator point_indices_last, PointIterator points_first, MultiLinestringRange linestrings, rmm::cuda_stream_view stream=rmm::cuda_stream_default, rmm::device_async_resource_ref mr=rmm::mr::get_current_device_resource())
 Finds the nearest linestring to each point in a quadrant, and computes the distances between each point and linestring.
 
std::unique_ptr< cudf::table > cuspatial::join_quadtree_and_bounding_boxes (cudf::table_view const &quadtree, cudf::table_view const &bbox, double x_min, double x_max, double y_min, double y_max, double scale, int8_t max_depth, rmm::device_async_resource_ref mr=rmm::mr::get_current_device_resource())
 Search a quadtree for polygon or linestring bounding box intersections.
 
std::unique_ptr< cudf::table > cuspatial::quadtree_point_in_polygon (cudf::table_view const &poly_quad_pairs, cudf::table_view const &quadtree, cudf::column_view const &point_indices, cudf::column_view const &point_x, cudf::column_view const &point_y, cudf::column_view const &poly_offsets, cudf::column_view const &ring_offsets, cudf::column_view const &poly_points_x, cudf::column_view const &poly_points_y, rmm::device_async_resource_ref mr=rmm::mr::get_current_device_resource())
 Test whether the specified points are inside any of the specified polygons.
 
std::unique_ptr< cudf::table > cuspatial::quadtree_point_to_nearest_linestring (cudf::table_view const &linestring_quad_pairs, cudf::table_view const &quadtree, cudf::column_view const &point_indices, cudf::column_view const &point_x, cudf::column_view const &point_y, cudf::column_view const &linestring_offsets, cudf::column_view const &linestring_points_x, cudf::column_view const &linestring_points_y, rmm::device_async_resource_ref mr=rmm::mr::get_current_device_resource())
 Finds the nearest linestring to each point in a quadrant, and computes the distances between each point and linestring.
 

Detailed Description

APIs related to spatial join.

Function Documentation

◆ join_quadtree_and_bounding_boxes() [1/2]

std::unique_ptr< cudf::table > cuspatial::join_quadtree_and_bounding_boxes ( cudf::table_view const & quadtree,
cudf::table_view const & bbox,
double x_min,
double x_max,
double y_min,
double y_max,
double scale,
int8_t max_depth,
rmm::device_async_resource_ref mr = rmm::mr::get_current_device_resource() )

Search a quadtree for polygon or linestring bounding box intersections.

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.
Parameters
quadtreecudf table representing a quadtree (key, level, is_internal_node, length, offset).
bboxcudf table of bounding boxes as four columns (x_min, y_min, x_max, y_max).
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 at which to stop testing for intersections.
mrThe optional resource to use for output device memory allocations.
Exceptions
cuspatial::logic_errorIf the quadtree table is malformed
cuspatial::logic_errorIf the bounding box table is malformed
cuspatial::logic_errorIf scale is less than or equal to 0
cuspatial::logic_errorIf x_min is greater than x_max
cuspatial::logic_errorIf y_min is greater than y_max
cuspatial::logic_errorIf max_depth is less than 1 or greater than 15
Returns
A cudf table with two columns:
  • bbox_offset - INT32 column of indices for each polygon/linestring bbox that intersects with the quadtree.
  • quad_offset - INT32 column of indices for each leaf quadrant intersecting with a polygon/linestring bbox.

◆ join_quadtree_and_bounding_boxes() [2/2]

template<class BoundingBoxIterator , class T = typename cuspatial::iterator_vec_base_type<BoundingBoxIterator>>
std::pair< rmm::device_uvector< uint32_t >, rmm::device_uvector< uint32_t > > cuspatial::join_quadtree_and_bounding_boxes ( point_quadtree_ref quadtree,
BoundingBoxIterator bounding_boxes_first,
BoundingBoxIterator bounding_boxes_last,
vec_2d< T > const & v_min,
T scale,
int8_t max_depth,
rmm::cuda_stream_view stream = rmm::cuda_stream_default,
rmm::device_async_resource_ref mr = rmm::mr::get_current_device_resource() )

Search a quadtree for polygon or linestring bounding box intersections.

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.
Parameters
quadtreeReference to a quadtree created using point_quadtree()
bounding_boxes_firststart bounding boxes iterator
bounding_boxes_lastend of bounding boxes iterator
v_minThe lower-left (x, y) corner 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 at which to stop testing for intersections.
streamThe CUDA stream on which to perform computations
mrThe optional resource to use for output device memory allocations.
Returns
A pair of UINT32 bounding box and leaf quadrant offset device vectors:
  • bbox_offset - indices for each polygon/linestring bbox that intersects with the quadtree.
  • quad_offset - indices for each leaf quadrant intersecting with a polygon/linestring bbox.
Exceptions
cuspatial::logic_errorIf scale is less than or equal to 0
cuspatial::logic_errorIf max_depth is less than 1 or greater than 15

◆ quadtree_point_in_polygon() [1/2]

std::unique_ptr< cudf::table > cuspatial::quadtree_point_in_polygon ( cudf::table_view const & poly_quad_pairs,
cudf::table_view const & quadtree,
cudf::column_view const & point_indices,
cudf::column_view const & point_x,
cudf::column_view const & point_y,
cudf::column_view const & poly_offsets,
cudf::column_view const & ring_offsets,
cudf::column_view const & poly_points_x,
cudf::column_view const & poly_points_y,
rmm::device_async_resource_ref mr = rmm::mr::get_current_device_resource() )

Test whether the specified points are inside any of the specified polygons.

Uses the table of (polygon, quadrant) pairs returned by cuspatial::join_quadtree_and_bounding_boxes to ensure only the points in the same quadrant as each polygon are tested for intersection.

This pre-filtering can dramatically reduce number of points tested per polygon, enabling faster intersection-testing at the expense of extra memory allocated to store the quadtree and sorted point_indices.

Parameters
poly_quad_pairscudf table of (polygon, quadrant) index pairs returned by cuspatial::join_quadtree_and_bounding_boxes
quadtreecudf table representing a quadtree (key, level, is_internal_node, length, offset).
point_indicesSorted point indices returned by cuspatial::quadtree_on_points
point_xx-coordinates of points to test
point_yy-coordinates of points to test
poly_offsetsBegin indices of the first ring in each polygon (i.e. prefix-sum).
ring_offsetsBegin indices of the first point in each ring (i.e. prefix-sum).
poly_points_xPolygon point x-coordinates.
poly_points_yPolygon point y-coordinates.
mrThe optional resource to use for output device memory allocations.
Exceptions
cuspatial::logic_errorIf the poly_quad_pairs table is malformed.
cuspatial::logic_errorIf the quadtree table is malformed.
cuspatial::logic_errorIf the number of point indices doesn't match the number of points.
cuspatial::logic_errorIf the number of rings is less than the number of polygons.
cuspatial::logic_errorIf any ring has fewer than three vertices.
cuspatial::logic_errorIf the types of point and polygon vertices are different.
Returns
A cudf table with two columns, where each row represents a point/polygon intersection: polygon_offset - UINT32 column of polygon indices point_offset - UINT32 column of point indices
Note
The returned polygon and point indices are offsets into the poly_quad_pairs inputs and point_indices, respectively.

◆ quadtree_point_in_polygon() [2/2]

template<class PolyIndexIterator , class QuadIndexIterator , class PointIndexIterator , class PointIterator , class MultiPolygonRange , class IndexType = iterator_value_type<PointIndexIterator>>
std::pair< rmm::device_uvector< IndexType >, rmm::device_uvector< IndexType > > cuspatial::quadtree_point_in_polygon ( PolyIndexIterator poly_indices_first,
PolyIndexIterator poly_indices_last,
QuadIndexIterator quad_indices_first,
point_quadtree_ref quadtree,
PointIndexIterator point_indices_first,
PointIndexIterator point_indices_last,
PointIterator points_first,
MultiPolygonRange polygons,
rmm::cuda_stream_view stream = rmm::cuda_stream_default,
rmm::device_async_resource_ref mr = rmm::mr::get_current_device_resource() )

Test whether the specified points are inside any of the specified polygons.

Uses the (polygon, quadrant) pairs returned by cuspatial::join_quadtree_and_bounding_boxes to ensure only the points in the same quadrant as each polygon are tested for intersection.

This pre-filtering can dramatically reduce the number of points tested per polygon, enabling faster intersection testing at the expense of extra memory allocated to store the quadtree and sorted point_indices.

Parameters
poly_indices_firstiterator to beginning of sequence of polygon indices returned by cuspatial::join_quadtree_and_bounding_boxes
poly_indices_firstiterator to end of sequence of polygon indices returned by cuspatial::join_quadtree_and_bounding_boxes
quad_indices_firstiterator to beginning of sequence of quadrant indices returned by cuspatial::join_quadtree_and_bounding_boxes
quadtreeReference to a quadtree created using point_quadtree()
point_indices_firstiterator to beginning of sequence of point indices returned by cuspatial::quadtree_on_points
point_indices_lastiterator to end of sequence of point indices returned by cuspatial::quadtree_on_points
points_firstiterator to beginning of sequence of (x, y) points to test
polygonsmultipolygon_range of polygons.
streamThe CUDA stream on which to perform computations
mrThe optional resource to use for output device memory allocations.
Exceptions
cuspatial::logic_errorIf the number of rings is less than the number of polygons.
cuspatial::logic_errorIf any ring has fewer than four vertices.
cuspatial::logic_errorif the number of multipolygons does not equal the total number of multipolygons (one polygon per multipolygon)
Returns
A pair of rmm::device_uvectors where each row represents a point/polygon intersection: polygon_offset - uint32_t polygon indices point_offset - uint32_t point indices
Note
Currently only supports single-polygon multipolygons.
The returned polygon and point indices are offsets into the poly_quad_pairs input range and point_indices range, respectively.

◆ quadtree_point_to_nearest_linestring() [1/2]

std::unique_ptr< cudf::table > cuspatial::quadtree_point_to_nearest_linestring ( cudf::table_view const & linestring_quad_pairs,
cudf::table_view const & quadtree,
cudf::column_view const & point_indices,
cudf::column_view const & point_x,
cudf::column_view const & point_y,
cudf::column_view const & linestring_offsets,
cudf::column_view const & linestring_points_x,
cudf::column_view const & linestring_points_y,
rmm::device_async_resource_ref mr = rmm::mr::get_current_device_resource() )

Finds the nearest linestring to each point in a quadrant, and computes the distances between each point and linestring.

Uses the table of (linestring, quadrant) pairs returned by cuspatial::join_quadtree_and_bounding_boxes to ensure distances are computed only for the points in the same quadrant as each linestring.

Parameters
linestring_quad_pairscudf table of (linestring, quadrant) index pairs returned by cuspatial::join_quadtree_and_bounding_boxes
quadtreecudf table representing a quadtree (key, level, is_internal_node, length, offset).
point_indicesSorted point indices returned by cuspatial::quadtree_on_points
point_xx-coordinates of points to test
point_yy-coordinates of points to test
linestring_offsetsBegin indices of the first point in each linestring (i.e. prefix-sum)
linestring_points_xLinestring point x-coordinates
linestring_points_yLinestring point y-coordinates
mrThe optional resource to use for output device memory allocations.
Exceptions
cuspatial::logic_errorIf the linestring_quad_pairs table is malformed.
cuspatial::logic_errorIf the quadtree table is malformed.
cuspatial::logic_errorIf the number of point indices doesn't match the number of points.
cuspatial::logic_errorIf any linestring has fewer than two vertices.
cuspatial::logic_errorIf the types of point and linestring vertices are different.
Returns
A cudf table with three columns, where each row represents a point/linestring pair and the distance between the two:

point_offset - UINT32 column of point indices linestring_offset - UINT32 column of linestring indices distance - FLOAT or DOUBLE column (based on input point data type) of distances between each point and linestring

Note
The returned point and linestring indices are offsets into the point_indices and linestring_quad_pairs inputs, respectively.

◆ quadtree_point_to_nearest_linestring() [2/2]

template<class LinestringIndexIterator , class QuadIndexIterator , class PointIndexIterator , class PointIterator , class MultiLinestringRange , typename IndexType = iterator_value_type<PointIndexIterator>, typename T = iterator_vec_base_type<PointIterator>>
std::tuple< rmm::device_uvector< IndexType >, rmm::device_uvector< IndexType >, rmm::device_uvector< T > > cuspatial::quadtree_point_to_nearest_linestring ( LinestringIndexIterator linestring_indices_first,
LinestringIndexIterator linestring_indices_last,
QuadIndexIterator quad_indices_first,
point_quadtree_ref quadtree,
PointIndexIterator point_indices_first,
PointIndexIterator point_indices_last,
PointIterator points_first,
MultiLinestringRange linestrings,
rmm::cuda_stream_view stream = rmm::cuda_stream_default,
rmm::device_async_resource_ref mr = rmm::mr::get_current_device_resource() )

Finds the nearest linestring to each point in a quadrant, and computes the distances between each point and linestring.

Uses the (linestring, quadrant) pairs returned by cuspatial::join_quadtree_and_bounding_boxes to ensure distances are computed only for the points in the same quadrant as each linestring.

Parameters
linestring_quad_pairscudf table of (linestring, quadrant) index pairs returned by cuspatial::join_quadtree_and_bounding_boxes
quadtreecudf table representing a quadtree (key, level, is_internal_node, length, offset).
point_indicesSorted point indices returned by cuspatial::quadtree_on_points
point_xx-coordinates of points to test
point_yy-coordinates of points to test
linestring_offsetsBegin indices of the first point in each linestring (i.e. prefix-sum)
linestring_points_xLinestring point x-coordinates
linestring_points_yLinestring point y-coordinates
mrThe optional resource to use for output device memory allocations.
Exceptions
cuspatial::logic_errorIf the linestring_quad_pairs table is malformed.
cuspatial::logic_errorIf the quadtree table is malformed.
cuspatial::logic_errorIf the number of point indices doesn't match the number of points.
cuspatial::logic_errorIf any linestring has fewer than two vertices.
cuspatial::logic_errorIf the types of point and linestring vertices are different.
Returns
A cudf table with three columns, where each row represents a point/linestring pair and the distance between the two:

point_offset - UINT32 column of point indices linestring_offset - UINT32 column of linestring indices distance - FLOAT or DOUBLE column (based on input point data type) of distances between each point and linestring

Note
The returned point and linestring indices are offsets into the point_indices and linestring_quad_pairs inputs, respectively.

◆ sinusoidal_projection()

template<class InputIt , class OutputIt , class T >
OutputIt cuspatial::sinusoidal_projection ( InputIt lon_lat_first,
InputIt lon_lat_last,
OutputIt xy_first,
vec_2d< T > origin,
rmm::cuda_stream_view stream = rmm::cuda_stream_default )

Sinusoidal projection of longitude/latitude relative to origin to Cartesian (x/y) coordinates in km.

Can be used to approximately convert longitude/latitude coordinates to Cartesian coordinates given that all points are near the origin. Error increases with distance from the origin. See Sinusoidal Projection for more detail.

Note
All input iterators must have a value_type of cuspatial::vec_2d<T> (Lat/Lon coordinates), and the output iterator must be able to accept for storage values of type cuspatial::vec_2d<T> (Cartesian coordinates).
Parameters
[in]lon_lat_firstbeginning of range of input longitude/latitude coordinates.
[in]lon_lat_lastend of range of input longitude/latitude coordinates.
[in]originlongitude and latitude of origin.
[out]xy_firstbeginning of range of output x/y coordinates.
[in]streamThe CUDA stream on which to perform computations and allocate memory.
Template Parameters
InputItIterator over longitude/latitude locations. Must meet the requirements of LegacyRandomAccessIterator and be device-accessible.
OutputItIterator over Cartesian output points. Must meet the requirements of LegacyRandomAccessIterator and be device-accessible and mutable.
Tthe floating-point coordinate value type of input longitude/latitude coordinates.
Precondition
lonlat_first may equal xy_first, but the range [lonlat_first, lonlat_last) shall not otherwise overlap the range [xy_first, xy_first + std::distance(lonlat_first, lonlat_last)).
Returns
Output iterator to the element past the last x/y coordinate computed.