libcuspatial  23.12.00
Loading...
Searching...
No Matches
Files | Functions
Spatial Relationship

APIs related to spatial relationship. More...

Files

file  bounding_boxes.hpp
 
file  bounding_boxes.cuh
 
file  point_in_polygon.hpp
 
file  point_in_polygon.cuh
 

Functions

template<typename IdInputIt , typename PointInputIt , typename BoundingBoxOutputIt , typename T = iterator_vec_base_type<PointInputIt>>
BoundingBoxOutputIt cuspatial::point_bounding_boxes (IdInputIt ids_first, IdInputIt ids_last, PointInputIt points_first, BoundingBoxOutputIt bounding_boxes_first, T expansion_radius=T{0}, rmm::cuda_stream_view stream=rmm::cuda_stream_default)
 Compute the spatial bounding boxes of sequences of points.
 
template<class LinestringOffsetIterator , class VertexIterator , class BoundingBoxIterator , class T = iterator_vec_base_type<VertexIterator>, class IndexT = iterator_value_type<LinestringOffsetIterator>>
BoundingBoxIterator cuspatial::linestring_bounding_boxes (LinestringOffsetIterator linestring_offsets_first, LinestringOffsetIterator linestring_offsets_last, VertexIterator linestring_vertices_first, VertexIterator linestring_vertices_last, BoundingBoxIterator bounding_boxes_first, T expansion_radius=T{0}, rmm::cuda_stream_view stream=rmm::cuda_stream_default)
 Compute minimum bounding box for each linestring.
 
template<class PolygonOffsetIterator , class RingOffsetIterator , class VertexIterator , class BoundingBoxIterator , class T = iterator_vec_base_type<VertexIterator>, class IndexT = iterator_value_type<PolygonOffsetIterator>>
BoundingBoxIterator cuspatial::polygon_bounding_boxes (PolygonOffsetIterator polygon_offsets_first, PolygonOffsetIterator polygon_offsets_last, RingOffsetIterator polygon_ring_offsets_first, RingOffsetIterator polygon_ring_offsets_last, VertexIterator polygon_vertices_first, VertexIterator polygon_vertices_last, BoundingBoxIterator bounding_boxes_first, T expansion_radius=T{0}, rmm::cuda_stream_view stream=rmm::cuda_stream_default)
 Compute minimum bounding box for each polygon.
 
std::unique_ptr< cudf::table > cuspatial::linestring_bounding_boxes (cudf::column_view const &linestring_offsets, cudf::column_view const &x, cudf::column_view const &y, double expansion_radius, rmm::mr::device_memory_resource *mr=rmm::mr::get_current_device_resource())
 Compute minimum bounding boxes of a set of linestrings and an expansion radius.
 
std::unique_ptr< cudf::table > cuspatial::polygon_bounding_boxes (cudf::column_view const &poly_offsets, cudf::column_view const &ring_offsets, cudf::column_view const &x, cudf::column_view const &y, double expansion_radius=0.0, rmm::mr::device_memory_resource *mr=rmm::mr::get_current_device_resource())
 Compute minimum bounding box for each polygon in a list.
 
template<class PointRange , class PolygonRange , class OutputIt >
OutputIt cuspatial::point_in_polygon (PointRange points, PolygonRange polygons, OutputIt output, rmm::cuda_stream_view stream=rmm::cuda_stream_default)
 Tests whether the specified points are inside any of the specified polygons.
 
template<typename PointRange , typename PolygonRange , typename OutputIt >
OutputIt cuspatial::pairwise_point_in_polygon (PointRange points, PolygonRange polygons, OutputIt results, rmm::cuda_stream_view stream=rmm::cuda_stream_default)
 Given (point, polygon) pairs, tests whether the point in the pair is in the polygon in the pair.
 
std::unique_ptr< cudf::column > cuspatial::point_in_polygon (cudf::column_view const &test_points_x, cudf::column_view const &test_points_y, cudf::column_view const &poly_offsets, cudf::column_view const &poly_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())
 Tests whether the specified points are inside any of the specified polygons.
 
std::unique_ptr< cudf::column > cuspatial::pairwise_point_in_polygon (cudf::column_view const &test_points_x, cudf::column_view const &test_points_y, cudf::column_view const &poly_offsets, cudf::column_view const &poly_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())
 Given (point, polygon pairs), tests whether the point of each pair is inside the polygon of the pair.
 

Detailed Description

APIs related to spatial relationship.

Function Documentation

◆ linestring_bounding_boxes() [1/2]

std::unique_ptr< cudf::table > cuspatial::linestring_bounding_boxes ( cudf::column_view const &  linestring_offsets,
cudf::column_view const &  x,
cudf::column_view const &  y,
double  expansion_radius,
rmm::mr::device_memory_resource *  mr = rmm::mr::get_current_device_resource() 
)

Compute minimum bounding boxes of a set of linestrings and an expansion radius.

Parameters
linestring_offsetsBegin indices of the first point in each linestring (i.e. prefix-sum)
xLinestring point x-coordinates
yLinestring point y-coordinates
expansion_radiusRadius of each linestring point
Returns
a cudf table of bounding boxes as four columns of the same type as x and y: x_min - the minimum x-coordinate of each bounding box y_min - the minimum y-coordinate of each bounding box x_max - the maximum x-coordinate of each bounding box y_max - the maximum y-coordinate of each bounding box
Precondition
For compatibility with GeoArrow, the size of linestring_offsets should be one more than the number of linestrings to process. The final offset is not used by this function, but the number of offsets determines the output size.

◆ linestring_bounding_boxes() [2/2]

template<class LinestringOffsetIterator , class VertexIterator , class BoundingBoxIterator , class T = iterator_vec_base_type<VertexIterator>, class IndexT = iterator_value_type<LinestringOffsetIterator>>
BoundingBoxIterator cuspatial::linestring_bounding_boxes ( LinestringOffsetIterator  linestring_offsets_first,
LinestringOffsetIterator  linestring_offsets_last,
VertexIterator  linestring_vertices_first,
VertexIterator  linestring_vertices_last,
BoundingBoxIterator  bounding_boxes_first,
expansion_radius = T{0},
rmm::cuda_stream_view  stream = rmm::cuda_stream_default 
)

Compute minimum bounding box for each linestring.

Template Parameters
LinestringOffsetIteratorIterator type to linestring offsets. Must meet the requirements of LegacyRandomAccessIterator and be device-readable.
VertexIteratorIterator type to linestring vertices. Must meet the requirements of LegacyRandomAccessIterator and be device-readable.
BoundingBoxIteratorIterator type to bounding boxes. Must be writable using data of type cuspatial::box<T>. Must meet the requirements of LegacyRandomAccessIterator and be device-writeable.
TThe coordinate data value type.
IndexTThe offset data value type.
Parameters
linestring_offsets_firstIterator to beginning of the range of input polygon offsets.
linestring_offsets_lastIterator to end of the range of input polygon offsets.
linestring_vertices_firstIterator to beginning of the range of input polygon vertices.
linestring_vertices_lastIterator to end of the range of input polygon vertices.
bounding_boxes_firstIterator to beginning of the range of output bounding boxes.
expansion_radiusOptional radius to expand each vertex of the output bounding boxes.
streamthe CUDA stream on which to perform computations and allocate memory.
Returns
An iterator to the end of the range of output bounding boxes.
Precondition
For compatibility with GeoArrow, the number of linestring offsets std::distance(linestring_offsets_first, linestring_offsets_last) should be one more than the number of linestrings. The final offset is not used by this function, but the number of offsets determines the output size.

◆ pairwise_point_in_polygon() [1/2]

std::unique_ptr< cudf::column > cuspatial::pairwise_point_in_polygon ( cudf::column_view const &  test_points_x,
cudf::column_view const &  test_points_y,
cudf::column_view const &  poly_offsets,
cudf::column_view const &  poly_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() 
)

Given (point, polygon pairs), tests whether the point of each pair is inside the polygon of the pair.

Tests that each point is or is not inside of the polygon in the corresponding index. Polygons are a collection of one or more * rings. Rings are a collection of three or more vertices.

Parameters
[in]test_points_xx-coordinates of points to test
[in]test_points_yy-coordinates of points to test
[in]poly_offsetsbeginning index of the first ring in each polygon
[in]poly_ring_offsetsbeginning index of the first point in each ring
[in]poly_points_xx-coordinates of polygon points
[in]poly_points_yy-coordinates of polygon points
Returns
A column of booleans for each point/polygon pair.
Note
Direction of rings does not matter.
Supports open or closed polygon formats.
This algorithm supports the ESRI shapefile format, but assumes all polygons are "clean" (as defined by the format), and does not verify whether the input adheres to the shapefile format.
Overlapping rings negate each other. This behavior is not limited to a single negation, allowing for "islands" within the same polygon.
poly_ring_offsets must contain only the rings that make up the polygons indexed by poly_offsets. If there are rings in poly_ring_offsets that are not part of the polygons in poly_offsets, results are likely to be incorrect and behavior is undefined.
poly w/two rings poly w/four rings
+-----------+ +------------------------+
:███████████: :████████████████████████:
:███████████: :██+------------------+██:
:██████+----:------+ :██: +----+ +----+ :██:
:██████: :██████: :██: :████: :████: :██:
+------;----+██████: :██: :----: :----: :██:
:███████████: :██+------------------+██:
:███████████: :████████████████████████:
+-----------+ +------------------------+

◆ pairwise_point_in_polygon() [2/2]

template<typename PointRange , typename PolygonRange , typename OutputIt >
OutputIt cuspatial::pairwise_point_in_polygon ( PointRange  points,
PolygonRange  polygons,
OutputIt  results,
rmm::cuda_stream_view  stream = rmm::cuda_stream_default 
)

Given (point, polygon) pairs, tests whether the point in the pair is in the polygon in the pair.

Note that the input must be a single geometry column, that is a (multi*)geometry_range initialized with counting iterator as the geometry offsets iterator.

Each input point will map to one uint8_t element in the output.

Template Parameters
PointRangean instance of template type multipoint_range, where GeometryIterator must be a counting iterator
PolygonRangean instance of template type multipolygon_range, where GeometryIterator must be a counting iterator
OutputItiterator type for output array. Must meet the requirements of LegacyRandomAccessIterator, be device-accessible, mutable and iterate on int32_t type.
Parameters
pointsRange of points, one per computed point-in-polygon pair,
polygonsRange of polygons, one per comptued point-in-polygon pair
outputbegin iterator to the output buffer
streamThe CUDA stream to use for kernel launches.
Returns
iterator to one past the last element in the output buffer
Note
Direction of rings does not matter.
The points of the rings must be explicitly closed.
Overlapping rings negate each other. This behavior is not limited to a single negation, allowing for "islands" within the same polygon.
poly w/two rings poly w/four rings
+-----------+ +------------------------+
:███████████: :████████████████████████:
:███████████: :██+------------------+██:
:██████+----:------+ :██: +----+ +----+ :██:
:██████: :██████: :██: :████: :████: :██:
+------;----+██████: :██: :----: :----: :██:
:███████████: :██+------------------+██:
:███████████: :████████████████████████:
+-----------+ +------------------------+
Precondition
Output iterator must be mutable and iterate on uint8_t type.

◆ point_bounding_boxes()

template<typename IdInputIt , typename PointInputIt , typename BoundingBoxOutputIt , typename T = iterator_vec_base_type<PointInputIt>>
BoundingBoxOutputIt cuspatial::point_bounding_boxes ( IdInputIt  ids_first,
IdInputIt  ids_last,
PointInputIt  points_first,
BoundingBoxOutputIt  bounding_boxes_first,
expansion_radius = T{0},
rmm::cuda_stream_view  stream = rmm::cuda_stream_default 
)

Compute the spatial bounding boxes of sequences of points.

Computes a bounding box around all points within each group (consecutive points with the same ID). This function can be applied to trajectory data, polygon vertices, linestring vertices, or any grouped point data.

Before merging bounding boxes, each point may be expanded into a bounding box using an optional expansion_radius. The point is expanded to a box with coordinates (point.x - expansion_radius, point.y - expansion_radius) and (point.x + expansion_radius, point.y + expansion_radius).

Note
Assumes Object IDs and points are presorted by ID.
Template Parameters
IdInputItIterator over object IDs. Must meet the requirements of LegacyRandomAccessIterator and be device-readable.
PointInputItIterator over points. Must meet the requirements of LegacyRandomAccessIterator and be device-readable.
BoundingBoxOutputItIterator over output bounding boxes. Each element is a tuple of two points representing corners of the axis-aligned bounding box. The type of the points is the same as the value_type of PointInputIt. Must meet the requirements of LegacyRandomAccessIterator and be device-writeable.
Parameters
ids_firstbeginning of the range of input object ids
ids_lastend of the range of input object ids
points_firstbeginning of the range of input point (x,y) coordinates
bounding_boxes_firstbeginning of the range of output bounding boxes, one per trajectory
expansion_radiusradius to add to each point when computing its bounding box.
streamthe CUDA stream on which to perform computations.
Returns
An iterator to the end of the range of output bounding boxes.

◆ point_in_polygon() [1/2]

std::unique_ptr< cudf::column > cuspatial::point_in_polygon ( cudf::column_view const &  test_points_x,
cudf::column_view const &  test_points_y,
cudf::column_view const &  poly_offsets,
cudf::column_view const &  poly_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() 
)

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

Tests whether points are inside at most 31 polygons. Polygons are a collection of one or more rings. Rings are a collection of three or more vertices.

Parameters
[in]test_points_xx-coordinates of points to test
[in]test_points_yy-coordinates of points to test
[in]poly_offsetsbeginning index of the first ring in each polygon
[in]poly_ring_offsetsbeginning index of the first point in each ring
[in]poly_points_xx-coordinates of polygon points
[in]poly_points_yy-coordinates of polygon points
Returns
A column of INT32 containing one element per input point. Each bit (except the sign bit) represents a hit or miss for each of the input polygons in least-significant-bit order. i.e. output[3] & 0b0010 indicates a hit or miss for the 3rd point against the 2nd polygon.
Note
Limit 31 polygons per call. Polygons may contain multiple rings.
Direction of rings does not matter.
This algorithm supports the ESRI shapefile format, but assumes all polygons are "clean" (as defined by the format), and does not verify whether the input adheres to the shapefile format.
Overlapping rings negate each other. This behavior is not limited to a single negation, allowing for "islands" within the same polygon.
poly_ring_offsets must contain only the rings that make up the polygons indexed by poly_offsets. If there are rings in poly_ring_offsets that are not part of the polygons in poly_offsets, results are likely to be incorrect and behavior is undefined.
poly w/two rings poly w/four rings
+-----------+ +------------------------+
:███████████: :████████████████████████:
:███████████: :██+------------------+██:
:██████+----:------+ :██: +----+ +----+ :██:
:██████: :██████: :██: :████: :████: :██:
+------;----+██████: :██: :----: :----: :██:
:███████████: :██+------------------+██:
:███████████: :████████████████████████:
+-----------+ +------------------------+

◆ point_in_polygon() [2/2]

template<class PointRange , class PolygonRange , class OutputIt >
OutputIt cuspatial::point_in_polygon ( PointRange  points,
PolygonRange  polygons,
OutputIt  output,
rmm::cuda_stream_view  stream = rmm::cuda_stream_default 
)

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

Tests whether points are inside at most 31 polygons. Polygons are a collection of one or more rings. Rings are a collection of three or more vertices.

Each input point will map to one int32_t element in the output. Each bit (except the sign bit) represents a hit or miss for each of the input polygons in least-significant-bit order. i.e. output[3] & 0b0010 indicates a hit or miss for the 3rd point against the 2nd polygon.

Note that the input must be a single geometry column, that is a (multi*)geometry_range initialized with counting iterator as the geometry offsets iterator.

Template Parameters
PointRangean instance of template type multipoint_range, where GeometryIterator must be a counting iterator
PolygonRangean instance of template type multipolygon_range, where GeometryIterator must be a counting iterator
OutputItiterator type for output array. Must meet the requirements of LegacyRandomAccessIterator, be device-accessible, mutable and iterate on int32_t type.
Parameters
pointsRange of points, one per computed point-in-polygon pair,
polygonsRange of polygons, one per comptued point-in-polygon pair
outputbegin iterator to the output buffer
streamThe CUDA stream to use for kernel launches.
Returns
iterator to one past the last element in the output buffer
Note
Direction of rings does not matter.
The points of the rings must be explicitly closed.
Overlapping rings negate each other. This behavior is not limited to a single negation, allowing for "islands" within the same polygon.
poly w/two rings poly w/four rings
+-----------+ +------------------------+
:███████████: :████████████████████████:
:███████████: :██+------------------+██:
:██████+----:------+ :██: +----+ +----+ :██:
:██████: :██████: :██: :████: :████: :██:
+------;----+██████: :██: :----: :----: :██:
:███████████: :██+------------------+██:
:███████████: :████████████████████████:
+-----------+ +------------------------+
Precondition
Output iterator must be mutable and iterate on int32_t type.

◆ polygon_bounding_boxes() [1/2]

std::unique_ptr< cudf::table > cuspatial::polygon_bounding_boxes ( cudf::column_view const &  poly_offsets,
cudf::column_view const &  ring_offsets,
cudf::column_view const &  x,
cudf::column_view const &  y,
double  expansion_radius = 0.0,
rmm::mr::device_memory_resource *  mr = rmm::mr::get_current_device_resource() 
)

Compute minimum bounding box for each polygon in a list.

Parameters
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)
xPolygon point x-coordinates
yPolygon point y-coordinates
expansion_radiusradius to add to each point when computing its bounding box.
Returns
a cudf table of bounding boxes as four columns of the same type as x and y: x_min - the minimum x-coordinate of each bounding box y_min - the minimum y-coordinate of each bounding box x_max - the maximum x-coordinate of each bounding box y_max - the maximum y-coordinate of each bounding box
Precondition
For compatibility with GeoArrow, the size of poly_offsets should be one more than the number of polygons to process. The size of ring_offsets should be one more than the number of total rings. The final offset in each range is not used by this function, but the number of polygon offsets determines the output size.

◆ polygon_bounding_boxes() [2/2]

template<class PolygonOffsetIterator , class RingOffsetIterator , class VertexIterator , class BoundingBoxIterator , class T = iterator_vec_base_type<VertexIterator>, class IndexT = iterator_value_type<PolygonOffsetIterator>>
BoundingBoxIterator cuspatial::polygon_bounding_boxes ( PolygonOffsetIterator  polygon_offsets_first,
PolygonOffsetIterator  polygon_offsets_last,
RingOffsetIterator  polygon_ring_offsets_first,
RingOffsetIterator  polygon_ring_offsets_last,
VertexIterator  polygon_vertices_first,
VertexIterator  polygon_vertices_last,
BoundingBoxIterator  bounding_boxes_first,
expansion_radius = T{0},
rmm::cuda_stream_view  stream = rmm::cuda_stream_default 
)

Compute minimum bounding box for each polygon.

Template Parameters
PolygonOffsetIteratorIterator type to polygon offsets. Must meet the requirements of LegacyRandomAccessIterator and be device-readable.
RingOffsetIteratorIterator type to polygon ring offsets. Must meet the requirements of LegacyRandomAccessIterator and be device-readable.
VertexIteratorIterator type to polygon vertices. Must meet the requirements of LegacyRandomAccessIterator and be device-readable.
BoundingBoxIteratorIterator type to bounding boxes. Must be writable using data of type cuspatial::box<T>. Must meet the requirements of LegacyRandomAccessIterator and be device-writeable.
TThe coordinate data value type.
IndexTThe offset data value type.
Parameters
polygon_offsets_firstIterator to beginning of the range of input polygon offsets.
polygon_offsets_lastIterator to end of the range of input polygon offsets.
polygon_ring_offsets_firstIterator to beginning of the range of input polygon ring offsets.
polygon_ring_offsets_lastIterator to end of the range of input polygon ring offsets.
polygon_vertices_firstIterator to beginning of the range of input polygon vertices.
polygon_vertices_lastIterator to end of the range of input polygon vertices.
bounding_boxes_firstIterator to beginning of the range of output bounding boxes.
expansion_radiusOptional radius to expand each vertex of the output bounding boxes.
streamthe CUDA stream on which to perform computations and allocate memory.
Returns
An iterator to the end of the range of output bounding boxes.
Precondition
For compatibility with GeoArrow, the number of polygon offsets std::distance(polygon_offsets_first, polygon_offsets_last) should be one more than the number of polygons. The number of ring offsets std::distance(polygon_ring_offsets_first, polygon_ring_offsets_last) should be one more than the number of total rings. The final offset in each range is not used by this function, but the number of polygon offsets determines the output size.