libcuspatial
24.02.00
|
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::mr::device_memory_resource *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::mr::device_memory_resource *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::mr::device_memory_resource *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::mr::device_memory_resource *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::mr::device_memory_resource *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::mr::device_memory_resource *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. | |
APIs related to spatial join.
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::mr::device_memory_resource * | mr = rmm::mr::get_current_device_resource() ) |
Search a quadtree for polygon or linestring bounding box intersections.
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.quadtree | cudf table representing a quadtree (key, level, is_internal_node, length, offset). |
bbox | cudf table of bounding boxes as four columns (x_min, y_min, x_max, y_max). |
x_min | The lower-left x-coordinate of the area of interest bounding box. |
x_max | The upper-right x-coordinate of the area of interest bounding box. |
y_min | The lower-left y-coordinate of the area of interest bounding box. |
y_max | The upper-right y-coordinate of the area of interest bounding box. |
scale | Scale to apply to each x and y distance from x_min and y_min. |
max_depth | Maximum quadtree depth at which to stop testing for intersections. |
mr | The optional resource to use for output device memory allocations. |
cuspatial::logic_error | If the quadtree table is malformed |
cuspatial::logic_error | If the bounding box table is malformed |
cuspatial::logic_error | If scale is less than or equal to 0 |
cuspatial::logic_error | If x_min is greater than x_max |
cuspatial::logic_error | If y_min is greater than y_max |
cuspatial::logic_error | If max_depth is less than 1 or greater than 15 |
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::mr::device_memory_resource * | mr = rmm::mr::get_current_device_resource() ) |
Search a quadtree for polygon or linestring bounding box intersections.
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.quadtree | Reference to a quadtree created using point_quadtree() |
bounding_boxes_first | start bounding boxes iterator |
bounding_boxes_last | end of bounding boxes iterator |
v_min | The lower-left (x, y) corner of the area of interest bounding box. |
scale | Scale to apply to each x and y distance from x_min and y_min. |
max_depth | Maximum quadtree depth at which to stop testing for intersections. |
stream | The CUDA stream on which to perform computations |
mr | The optional resource to use for output device memory allocations. |
cuspatial::logic_error | If scale is less than or equal to 0 |
cuspatial::logic_error | If max_depth is less than 1 or greater than 15 |
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::mr::device_memory_resource * | 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.
poly_quad_pairs | cudf table of (polygon, quadrant) index pairs returned by cuspatial::join_quadtree_and_bounding_boxes |
quadtree | cudf table representing a quadtree (key, level, is_internal_node, length, offset). |
point_indices | Sorted point indices returned by cuspatial::quadtree_on_points |
point_x | x-coordinates of points to test |
point_y | y-coordinates of points to test |
poly_offsets | Begin indices of the first ring in each polygon (i.e. prefix-sum). |
ring_offsets | Begin indices of the first point in each ring (i.e. prefix-sum). |
poly_points_x | Polygon point x-coordinates. |
poly_points_y | Polygon point y-coordinates. |
mr | The optional resource to use for output device memory allocations. |
cuspatial::logic_error | If the poly_quad_pairs table is malformed. |
cuspatial::logic_error | If the quadtree table is malformed. |
cuspatial::logic_error | If the number of point indices doesn't match the number of points. |
cuspatial::logic_error | If the number of rings is less than the number of polygons. |
cuspatial::logic_error | If any ring has fewer than three vertices. |
cuspatial::logic_error | If the types of point and polygon vertices are different. |
poly_quad_pairs
inputs and point_indices
, respectively. 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::mr::device_memory_resource * | 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.
poly_indices_first | iterator to beginning of sequence of polygon indices returned by cuspatial::join_quadtree_and_bounding_boxes |
poly_indices_first | iterator to end of sequence of polygon indices returned by cuspatial::join_quadtree_and_bounding_boxes |
quad_indices_first | iterator to beginning of sequence of quadrant indices returned by cuspatial::join_quadtree_and_bounding_boxes |
quadtree | Reference to a quadtree created using point_quadtree() |
point_indices_first | iterator to beginning of sequence of point indices returned by cuspatial::quadtree_on_points |
point_indices_last | iterator to end of sequence of point indices returned by cuspatial::quadtree_on_points |
points_first | iterator to beginning of sequence of (x, y) points to test |
polygons | multipolygon_range of polygons. |
stream | The CUDA stream on which to perform computations |
mr | The optional resource to use for output device memory allocations. |
cuspatial::logic_error | If the number of rings is less than the number of polygons. |
cuspatial::logic_error | If any ring has fewer than four vertices. |
cuspatial::logic_error | if the number of multipolygons does not equal the total number of multipolygons (one polygon per multipolygon) |
poly_quad_pairs
input range and point_indices
range, respectively. 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::mr::device_memory_resource * | 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.
linestring_quad_pairs | cudf table of (linestring, quadrant) index pairs returned by cuspatial::join_quadtree_and_bounding_boxes |
quadtree | cudf table representing a quadtree (key, level, is_internal_node, length, offset). |
point_indices | Sorted point indices returned by cuspatial::quadtree_on_points |
point_x | x-coordinates of points to test |
point_y | y-coordinates of points to test |
linestring_offsets | Begin indices of the first point in each linestring (i.e. prefix-sum) |
linestring_points_x | Linestring point x-coordinates |
linestring_points_y | Linestring point y-coordinates |
mr | The optional resource to use for output device memory allocations. |
cuspatial::logic_error | If the linestring_quad_pairs table is malformed. |
cuspatial::logic_error | If the quadtree table is malformed. |
cuspatial::logic_error | If the number of point indices doesn't match the number of points. |
cuspatial::logic_error | If any linestring has fewer than two vertices. |
cuspatial::logic_error | If the types of point and linestring vertices are different. |
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
point_indices
and linestring_quad_pairs
inputs, respectively. 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::mr::device_memory_resource * | 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.
linestring_quad_pairs | cudf table of (linestring, quadrant) index pairs returned by cuspatial::join_quadtree_and_bounding_boxes |
quadtree | cudf table representing a quadtree (key, level, is_internal_node, length, offset). |
point_indices | Sorted point indices returned by cuspatial::quadtree_on_points |
point_x | x-coordinates of points to test |
point_y | y-coordinates of points to test |
linestring_offsets | Begin indices of the first point in each linestring (i.e. prefix-sum) |
linestring_points_x | Linestring point x-coordinates |
linestring_points_y | Linestring point y-coordinates |
mr | The optional resource to use for output device memory allocations. |
cuspatial::logic_error | If the linestring_quad_pairs table is malformed. |
cuspatial::logic_error | If the quadtree table is malformed. |
cuspatial::logic_error | If the number of point indices doesn't match the number of points. |
cuspatial::logic_error | If any linestring has fewer than two vertices. |
cuspatial::logic_error | If the types of point and linestring vertices are different. |
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
point_indices
and linestring_quad_pairs
inputs, respectively. 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.
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).[in] | lon_lat_first | beginning of range of input longitude/latitude coordinates. |
[in] | lon_lat_last | end of range of input longitude/latitude coordinates. |
[in] | origin | longitude and latitude of origin. |
[out] | xy_first | beginning of range of output x/y coordinates. |
[in] | stream | The CUDA stream on which to perform computations and allocate memory. |
InputIt | Iterator over longitude/latitude locations. Must meet the requirements of LegacyRandomAccessIterator and be device-accessible. |
OutputIt | Iterator over Cartesian output points. Must meet the requirements of LegacyRandomAccessIterator and be device-accessible and mutable. |
T | the floating-point coordinate value type of input longitude/latitude coordinates. |
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))
.