GeometryOps
Documentation for GeometryOps.
GeometryOps.AbstractBarycentricCoordinateMethod
GeometryOps.DouglasPeucker
GeometryOps.MeanValue
GeometryOps.RadialDistance
GeometryOps.SimplifyAlg
GeometryOps.VisvalingamWhyatt
GeometryOps._det
GeometryOps._equals_curves
GeometryOps._intersection_point
GeometryOps._line_intersects
GeometryOps._line_intersects
GeometryOps._overlaps
GeometryOps.apply
GeometryOps.area
GeometryOps.centroid
GeometryOps.centroid
GeometryOps.centroid
GeometryOps.centroid_and_area
GeometryOps.centroid_and_area
GeometryOps.centroid_and_area
GeometryOps.centroid_and_area
GeometryOps.centroid_and_length
GeometryOps.centroid_and_length
GeometryOps.contains
GeometryOps.crosses
GeometryOps.disjoint
GeometryOps.distance
GeometryOps.embed_extent
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.euclid_distance
GeometryOps.flatten
GeometryOps.flip
GeometryOps.get_contours
GeometryOps.intersection
GeometryOps.intersection
GeometryOps.intersection
GeometryOps.intersection
GeometryOps.intersection_points
GeometryOps.intersection_points
GeometryOps.intersects
GeometryOps.intersects
GeometryOps.intersects
GeometryOps.isclockwise
GeometryOps.isconcave
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.point_in_polygon
GeometryOps.point_on_line
GeometryOps.polygon_to_line
GeometryOps.polygonize
GeometryOps.rebuild
GeometryOps.reconstruct
GeometryOps.reproject
GeometryOps.signed_area
GeometryOps.signed_distance
GeometryOps.simplify
GeometryOps.t_value
GeometryOps.to_edges
GeometryOps.tuples
GeometryOps.unwrap
GeometryOps.weighted_mean
GeometryOps.AbstractBarycentricCoordinateMethod
— Typeabstract type AbstractBarycentricCoordinateMethod
Abstract supertype for barycentric coordinate methods. The subtypes may serve as dispatch types, or may cache some information about the target polygon.
API
The following methods must be implemented for all subtypes:
barycentric_coordinates!(λs::Vector{<: Real}, method::AbstractBarycentricCoordinateMethod, exterior::Vector{<: Point{2, T1}}, point::Point{2, T2})
barycentric_interpolate(method::AbstractBarycentricCoordinateMethod, exterior::Vector{<: Point{2, T1}}, values::Vector{V}, point::Point{2, T2})::V
barycentric_interpolate(method::AbstractBarycentricCoordinateMethod, exterior::Vector{<: Point{2, T1}}, interiors::Vector{<: Vector{<: Point{2, T1}}} values::Vector{V}, point::Point{2, T2})::V
The rest of the methods will be implemented in terms of these, and have efficient dispatches for broadcasting.
GeometryOps.DouglasPeucker
— TypeDouglasPeucker <: SimplifyAlg
+Home · GeometryOps.jl GeometryOps
Documentation for GeometryOps.
GeometryOps.AbstractBarycentricCoordinateMethod
GeometryOps.DouglasPeucker
GeometryOps.MeanValue
GeometryOps.RadialDistance
GeometryOps.SimplifyAlg
GeometryOps.VisvalingamWhyatt
GeometryOps._det
GeometryOps._equals_curves
GeometryOps._intersection_point
GeometryOps._line_intersects
GeometryOps._line_intersects
GeometryOps._overlaps
GeometryOps.apply
GeometryOps.area
GeometryOps.centroid
GeometryOps.centroid
GeometryOps.centroid
GeometryOps.centroid_and_area
GeometryOps.centroid_and_area
GeometryOps.centroid_and_area
GeometryOps.centroid_and_area
GeometryOps.centroid_and_length
GeometryOps.centroid_and_length
GeometryOps.contains
GeometryOps.crosses
GeometryOps.disjoint
GeometryOps.distance
GeometryOps.embed_extent
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.equals
GeometryOps.euclid_distance
GeometryOps.flatten
GeometryOps.flip
GeometryOps.get_contours
GeometryOps.intersection
GeometryOps.intersection
GeometryOps.intersection
GeometryOps.intersection
GeometryOps.intersection_points
GeometryOps.intersection_points
GeometryOps.intersects
GeometryOps.intersects
GeometryOps.intersects
GeometryOps.isclockwise
GeometryOps.isconcave
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.overlaps
GeometryOps.point_in_polygon
GeometryOps.point_on_line
GeometryOps.polygon_to_line
GeometryOps.polygonize
GeometryOps.rebuild
GeometryOps.reconstruct
GeometryOps.reproject
GeometryOps.signed_area
GeometryOps.signed_distance
GeometryOps.simplify
GeometryOps.t_value
GeometryOps.to_edges
GeometryOps.tuples
GeometryOps.unwrap
GeometryOps.weighted_mean
GeometryOps.AbstractBarycentricCoordinateMethod
— Typeabstract type AbstractBarycentricCoordinateMethod
Abstract supertype for barycentric coordinate methods. The subtypes may serve as dispatch types, or may cache some information about the target polygon.
API
The following methods must be implemented for all subtypes:
barycentric_coordinates!(λs::Vector{<: Real}, method::AbstractBarycentricCoordinateMethod, exterior::Vector{<: Point{2, T1}}, point::Point{2, T2})
barycentric_interpolate(method::AbstractBarycentricCoordinateMethod, exterior::Vector{<: Point{2, T1}}, values::Vector{V}, point::Point{2, T2})::V
barycentric_interpolate(method::AbstractBarycentricCoordinateMethod, exterior::Vector{<: Point{2, T1}}, interiors::Vector{<: Vector{<: Point{2, T1}}} values::Vector{V}, point::Point{2, T2})::V
The rest of the methods will be implemented in terms of these, and have efficient dispatches for broadcasting.
sourceGeometryOps.DouglasPeucker
— TypeDouglasPeucker <: SimplifyAlg
-DouglasPeucker(; number, ratio, tol)
Simplifies geometries by removing points below tol
distance from the line between its neighboring points.
Keywords
ratio
: the fraction of points that should remain after simplify
. Useful as it will generalise for large collections of objects.
number
: the number of points that should remain after simplify
. Less useful for large collections of mixed size objects.
tol
: the minimum distance a point will be from the line joining its neighboring points.
sourceGeometryOps.MeanValue
— TypeMeanValue() <: AbstractBarycentricCoordinateMethod
This method calculates barycentric coordinates using the mean value method.
References
sourceGeometryOps.RadialDistance
— TypeRadialDistance <: SimplifyAlg
Simplifies geometries by removing points less than tol
distance from the line between its neighboring points.
Keywords
ratio
: the fraction of points that should remain after simplify
. Useful as it will generalise for large collections of objects.
number
: the number of points that should remain after simplify
. Less useful for large collections of mixed size objects.
tol
: the minimum distance between points.
sourceGeometryOps.SimplifyAlg
— Typeabstract type SimplifyAlg
Abstract type for simplification algorithms.
API
For now, the algorithm must hold the number
, ratio
and tol
properties.
Simplification algorithm types can hook into the interface by implementing the _simplify(trait, alg, geom)
methods for whichever traits are necessary.
sourceGeometryOps.VisvalingamWhyatt
— TypeVisvalingamWhyatt <: SimplifyAlg
+DouglasPeucker(; number, ratio, tol)
Simplifies geometries by removing points below tol
distance from the line between its neighboring points.
Keywords
ratio
: the fraction of points that should remain after simplify
. Useful as it will generalise for large collections of objects.
number
: the number of points that should remain after simplify
. Less useful for large collections of mixed size objects.
tol
: the minimum distance a point will be from the line joining its neighboring points.
sourceGeometryOps.MeanValue
— TypeMeanValue() <: AbstractBarycentricCoordinateMethod
This method calculates barycentric coordinates using the mean value method.
References
sourceGeometryOps.RadialDistance
— TypeRadialDistance <: SimplifyAlg
Simplifies geometries by removing points less than tol
distance from the line between its neighboring points.
Keywords
ratio
: the fraction of points that should remain after simplify
. Useful as it will generalise for large collections of objects.
number
: the number of points that should remain after simplify
. Less useful for large collections of mixed size objects.
tol
: the minimum distance between points.
sourceGeometryOps.SimplifyAlg
— Typeabstract type SimplifyAlg
Abstract type for simplification algorithms.
API
For now, the algorithm must hold the number
, ratio
and tol
properties.
Simplification algorithm types can hook into the interface by implementing the _simplify(trait, alg, geom)
methods for whichever traits are necessary.
sourceGeometryOps.VisvalingamWhyatt
— TypeVisvalingamWhyatt <: SimplifyAlg
-VisvalingamWhyatt(; kw...)
Simplifies geometries by removing points below tol
distance from the line between its neighboring points.
Keywords
ratio
: the fraction of points that should remain after simplify
. Useful as it will generalise for large collections of objects.
number
: the number of points that should remain after simplify
. Less useful for large collections of mixed size objects.
tol
: the minimum area of a triangle made with a point and its neighboring points.
sourceGeometryOps._det
— Method_det(s1::Point2{T1}, s2::Point2{T2}) where {T1 <: Real, T2 <: Real}
Returns the determinant of the matrix formed by hcat
'ing two points s1
and s2
.
Specifically, this is:
s1[1] * s2[2] - s1[2] * s2[1]
sourceGeometryOps._equals_curves
— Method_equals_curves(c1, c2, closed_type1, closed_type2)::Bool
Two curves are equal if they share the same set of point, representing the same geometry. Both curves must must be composed of the same set of points, however, they do not have to wind in the same direction, or start on the same point to be equivalent. Inputs: c1 first geometry c2 second geometry closedtype1::Bool true if c1 is closed by definition (polygon, linear ring) closedtype2::Bool true if c2 is closed by definition (polygon, linear ring)
sourceGeometryOps._intersection_point
— Method_intersection_point(
+VisvalingamWhyatt(; kw...)
Simplifies geometries by removing points below tol
distance from the line between its neighboring points.
Keywords
ratio
: the fraction of points that should remain after simplify
. Useful as it will generalise for large collections of objects.
number
: the number of points that should remain after simplify
. Less useful for large collections of mixed size objects.
tol
: the minimum area of a triangle made with a point and its neighboring points.
sourceGeometryOps._det
— Method_det(s1::Point2{T1}, s2::Point2{T2}) where {T1 <: Real, T2 <: Real}
Returns the determinant of the matrix formed by hcat
'ing two points s1
and s2
.
Specifically, this is:
s1[1] * s2[2] - s1[2] * s2[1]
sourceGeometryOps._equals_curves
— Method_equals_curves(c1, c2, closed_type1, closed_type2)::Bool
Two curves are equal if they share the same set of point, representing the same geometry. Both curves must must be composed of the same set of points, however, they do not have to wind in the same direction, or start on the same point to be equivalent. Inputs: c1 first geometry c2 second geometry closedtype1::Bool true if c1 is closed by definition (polygon, linear ring) closedtype2::Bool true if c2 is closed by definition (polygon, linear ring)
sourceGeometryOps._intersection_point
— Method_intersection_point(
(a1, a2)::Tuple,
(b1, b2)::Tuple,
-)
Calculates the intersection point between two lines if it exists, and as if the line extended to infinity, and the fractional component of each line from the initial end point to the intersection point. Inputs: (a1, a2)::Tuple{Tuple{::Real, ::Real}, Tuple{::Real, ::Real}} first line (b1, b2)::Tuple{Tuple{::Real, ::Real}, Tuple{::Real, ::Real}} second line Outputs: (x, y)::Tuple{::Real, ::Real} intersection point (t, u)::Tuple{::Real, ::Real} fractional length of lines to intersection Both are ::Nothing if point doesn't exist!
Calculation derivation can be found here: https://stackoverflow.com/questions/563198/
sourceGeometryOps._line_intersects
— Method_line_intersects(
+)
Calculates the intersection point between two lines if it exists, and as if the line extended to infinity, and the fractional component of each line from the initial end point to the intersection point. Inputs: (a1, a2)::Tuple{Tuple{::Real, ::Real}, Tuple{::Real, ::Real}} first line (b1, b2)::Tuple{Tuple{::Real, ::Real}, Tuple{::Real, ::Real}} second line Outputs: (x, y)::Tuple{::Real, ::Real} intersection point (t, u)::Tuple{::Real, ::Real} fractional length of lines to intersection Both are ::Nothing if point doesn't exist!
Calculation derivation can be found here: https://stackoverflow.com/questions/563198/
sourceGeometryOps._line_intersects
— Method_line_intersects(
edge_a::Edge,
edge_b::Edge,
-)::Bool
Returns true if there is at least one intersection between two edges.
sourceGeometryOps._line_intersects
— Method_line_intersects(
+)::Bool
Returns true if there is at least one intersection between two edges.
sourceGeometryOps._line_intersects
— Method_line_intersects(
edges_a::Vector{Edge},
edges_b::Vector{Edge}
-)::Bool
Returns true if there is at least one intersection between edges within the two lists of edges.
sourceGeometryOps._overlaps
— Method_overlaps(
+)::Bool
Returns true if there is at least one intersection between edges within the two lists of edges.
sourceGeometryOps._overlaps
— Method_overlaps(
(a1, a2)::Edge,
(b1, b2)::Edge
-)::Bool
If the edges overlap, meaning that they are colinear but each have one endpoint outside of the other edge, return true. Else false.
sourceGeometryOps.apply
— Methodapply(f, target::Type{<:AbstractTrait}, obj; kw...)
Reconstruct a geometry, feature, feature collection or nested vectors of either using the function f
on the target
trait.
f(target_geom) => x
where x
also has the target
trait, or a trait that can be substituted. For example, swapping PolgonTrait
to MultiPointTrait
will fail if the outer object has MultiPolygonTrait
, but should work if it has FeatureTrait
.
Objects "shallower" than the target trait are always completely rebuilt, like a Vector
of FeatureCollectionTrait
of FeatureTrait
when the target has PolygonTrait
and is held in the features. But "deeper" opjects may remain unchanged - such as points and linear rings if the tartet is the same PolygonTrait
.
The result is an functionally similar geometry with values depending on f
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
Example
Flipped point the order in any feature or geometry, or iterables of either:
```juia import GeoInterface as GI import GeometryOps as GO geom = GI.Polygon([GI.LinearRing([(1, 2), (3, 4), (5, 6), (1, 2)]), GI.LinearRing([(3, 4), (5, 6), (6, 7), (3, 4)])])
flipped_geom = GO.apply(GI.PointTrait, geom) do p (GI.y(p), GI.x(p)) end
sourceGeometryOps.area
— Methodarea(geom)::Real
Returns the area of the geometry. This is computed slighly differently for different geometries: - The area of a point is always zero. - The area of a curve is always zero. - The area of a polygon is the absolute value of the signed area. - The area multi-polygon is the sum of the areas of all of the sub-polygons.
sourceGeometryOps.centroid
— Methodcentroid(trait, geom)::Tuple{T, T}
Returns the centroid of a polygon or multipolygon, which is calculated by weighting edges by their area component
by convention.
sourceGeometryOps.centroid
— Methodcentroid(geom)::Tuple{T, T}
Returns the centroid of a given line segment, linear ring, polygon, or mutlipolygon.
sourceGeometryOps.centroid
— Methodcentroid(
+)::Bool
If the edges overlap, meaning that they are colinear but each have one endpoint outside of the other edge, return true. Else false.
sourceGeometryOps.apply
— Methodapply(f, target::Type{<:AbstractTrait}, obj; kw...)
Reconstruct a geometry, feature, feature collection, or nested vectors of either using the function f
on the target
trait.
f(target_geom) => x
where x
also has the target
trait, or a trait that can be substituted. For example, swapping PolgonTrait
to MultiPointTrait
will fail if the outer object has MultiPolygonTrait
, but should work if it has FeatureTrait
.
Objects "shallower" than the target trait are always completely rebuilt, like a Vector
of FeatureCollectionTrait
of FeatureTrait
when the target has PolygonTrait
and is held in the features. But "deeper" objects may remain unchanged - such as points and linear rings if the target is the same PolygonTrait
.
The result is a functionally similar geometry with values depending on f
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
Example
Flipped point the order in any feature or geometry, or iterables of either:
```juia import GeoInterface as GI import GeometryOps as GO geom = GI.Polygon([GI.LinearRing([(1, 2), (3, 4), (5, 6), (1, 2)]), GI.LinearRing([(3, 4), (5, 6), (6, 7), (3, 4)])])
flipped_geom = GO.apply(GI.PointTrait, geom) do p (GI.y(p), GI.x(p)) end
sourceGeometryOps.area
— Methodarea(geom)::Real
Returns the area of the geometry. This is computed slighly differently for different geometries: - The area of a point is always zero. - The area of a curve is always zero. - The area of a polygon is the absolute value of the signed area. - The area multi-polygon is the sum of the areas of all of the sub-polygons.
sourceGeometryOps.centroid
— Methodcentroid(trait, geom)::Tuple{T, T}
Returns the centroid of a polygon or multipolygon, which is calculated by weighting edges by their area component
by convention.
sourceGeometryOps.centroid
— Methodcentroid(geom)::Tuple{T, T}
Returns the centroid of a given line segment, linear ring, polygon, or mutlipolygon.
sourceGeometryOps.centroid
— Methodcentroid(
trait::Union{GI.LineStringTrait, GI.LinearRingTrait},
geom,
-)::Tuple{T, T}
Returns the centroid of a line string or linear ring, which is calculated by weighting line segments by their length by convention.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(
+)::Tuple{T, T}
Returns the centroid of a line string or linear ring, which is calculated by weighting line segments by their length by convention.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(
::Union{GI.LineStringTrait, GI.LinearRingTrait},
geom,
-)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given geom.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(::GI.MultiPolygonTrait, geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given multipolygon.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(::GI.PolygonTrait, geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given polygon.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(
+)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given geom.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(::GI.MultiPolygonTrait, geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given multipolygon.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(::GI.PolygonTrait, geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given polygon.
sourceGeometryOps.centroid_and_area
— Methodcentroid_and_area(
::Union{GI.LineStringTrait, GI.LinearRingTrait},
geom,
-)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given a line string or a linear ring. Note that this is only valid if the line segment or linear ring is closed.
sourceGeometryOps.centroid_and_length
— Methodcentroid_and_length(geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and length of a given line/ring. Note this is only valid for line strings and linear rings.
sourceGeometryOps.centroid_and_length
— Methodcentroid_and_length(geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and length of a given line/ring. Note this is only valid for line strings and linear rings.
sourceGeometryOps.contains
— Methodcontains(ft1::AbstractGeometry, ft2::AbstractGeometry)::Bool
Return true if the second geometry is completely contained by the first geometry. The interiors of both geometries must intersect and, the interior and boundary of the secondary (geometry b) must not intersect the exterior of the primary (geometry a). contains
returns the exact opposite result of within
.
Examples
import GeometryOps as GO, GeoInterface as GI
+)::(::Tuple{T, T}, ::Real)
Returns the centroid and area of a given a line string or a linear ring. Note that this is only valid if the line segment or linear ring is closed.
sourceGeometryOps.centroid_and_length
— Methodcentroid_and_length(geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and length of a given line/ring. Note this is only valid for line strings and linear rings.
sourceGeometryOps.centroid_and_length
— Methodcentroid_and_length(geom)::(::Tuple{T, T}, ::Real)
Returns the centroid and length of a given line/ring. Note this is only valid for line strings and linear rings.
sourceGeometryOps.contains
— Methodcontains(ft1::AbstractGeometry, ft2::AbstractGeometry)::Bool
Return true if the second geometry is completely contained by the first geometry. The interiors of both geometries must intersect and, the interior and boundary of the secondary (geometry b) must not intersect the exterior of the primary (geometry a). contains
returns the exact opposite result of within
.
Examples
import GeometryOps as GO, GeoInterface as GI
line = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
point = (1, 2)
GO.contains(line, point)
# output
-true
sourceGeometryOps.crosses
— Method crosses(geom1, geom2)::Bool
Return true
if the intersection results in a geometry whose dimension is one less than the maximum dimension of the two source geometries and the intersection set is interior to both source geometries.
TODO: broken
Examples
import GeoInterface as GI, GeometryOps as GO
+true
sourceGeometryOps.crosses
— Method crosses(geom1, geom2)::Bool
Return true
if the intersection results in a geometry whose dimension is one less than the maximum dimension of the two source geometries and the intersection set is interior to both source geometries.
TODO: broken
Examples
import GeoInterface as GI, GeometryOps as GO
line1 = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
line2 = GI.LineString([(-2, 2), (4, 2)])
GO.crosses(line1, line2)
# output
-true
sourceGeometryOps.disjoint
— Methoddisjoint(geom1, geom2)::Bool
Return true
if the intersection of the two geometries is an empty set.
Examples
import GeometryOps as GO, GeoInterface as GI
+true
sourceGeometryOps.disjoint
— Methoddisjoint(geom1, geom2)::Bool
Return true
if the intersection of the two geometries is an empty set.
Examples
import GeometryOps as GO, GeoInterface as GI
poly = GI.Polygon([[(-1, 2), (3, 2), (3, 3), (-1, 3), (-1, 2)]])
point = (1, 1)
GO.disjoint(poly, point)
# output
-true
sourceGeometryOps.distance
— Methoddistance(point, geom)::Real
Calculates the ditance from the geometry g1
to the point
. The distance will always be positive or zero.
The method will differ based on the type of the geometry provided: - The distance from a point to a point is just the Euclidean distance between the points. - The distance from a point to a multipolygon is the shortest distance from a the given point to any point within the multipoint object. - The distance from a point to a line is the minimum distance from the point to the closest point on the given line. - The distance from a point to a linestring is the minimum distance from the point to the closest segment of the linestring. - The distance from a point to a linear ring is the minimum distance from the point to the closest segment of the linear ring. - The distance from a point to a polygon is zero if the point is within the polygon and otherwise is the minimum distance from the point to an edge of the polygon. This includes edges created by holes. - The distance from a point to a multipolygon is zero if the point is within the multipolygon and otherwise is the minimum distance from the point to the closest edge of any of the polygons within the multipolygon. This includes edges created by holes of the polygons as well.
sourceGeometryOps.embed_extent
— Methodembed_extent(obj)
Recursively wrap the object with a GeoInterface.jl geometry, calculating and adding an Extents.Extent
to all objects.
This can improve performance when extents need to be checked multiple times, such when needing to check if many points are in geometries, and using their extents as a quick filter for obviously exterior points.
Keywords
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.
sourceGeometryOps.equals
— Methodequals(trait_a, geom_a, trait_b, geom_b)
Two geometries which are not of the same type cannot be equal so they always return false.
sourceGeometryOps.equals
— Methodequals(geom1, geom2)::Bool
Compare two Geometries return true if they are the same geometry.
Examples
import GeometryOps as GO, GeoInterface as GI
+true
sourceGeometryOps.distance
— Methoddistance(point, geom)::Real
Calculates the ditance from the geometry g1
to the point
. The distance will always be positive or zero.
The method will differ based on the type of the geometry provided: - The distance from a point to a point is just the Euclidean distance between the points. - The distance from a point to a multipolygon is the shortest distance from a the given point to any point within the multipoint object. - The distance from a point to a line is the minimum distance from the point to the closest point on the given line. - The distance from a point to a linestring is the minimum distance from the point to the closest segment of the linestring. - The distance from a point to a linear ring is the minimum distance from the point to the closest segment of the linear ring. - The distance from a point to a polygon is zero if the point is within the polygon and otherwise is the minimum distance from the point to an edge of the polygon. This includes edges created by holes. - The distance from a point to a multipolygon is zero if the point is within the multipolygon and otherwise is the minimum distance from the point to the closest edge of any of the polygons within the multipolygon. This includes edges created by holes of the polygons as well.
sourceGeometryOps.embed_extent
— Methodembed_extent(obj)
Recursively wrap the object with a GeoInterface.jl geometry, calculating and adding an Extents.Extent
to all objects.
This can improve performance when extents need to be checked multiple times, such when needing to check if many points are in geometries, and using their extents as a quick filter for obviously exterior points.
Keywords
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.
sourceGeometryOps.equals
— Methodequals(trait_a, geom_a, trait_b, geom_b)
Two geometries which are not of the same type cannot be equal so they always return false.
sourceGeometryOps.equals
— Methodequals(geom1, geom2)::Bool
Compare two Geometries return true if they are the same geometry.
Examples
import GeometryOps as GO, GeoInterface as GI
poly1 = GI.Polygon([[(0,0), (0,5), (5,5), (5,0), (0,0)]])
poly2 = GI.Polygon([[(0,0), (0,5), (5,5), (5,0), (0,0)]])
GO.equals(poly1, poly2)
# output
-true
sourceGeometryOps.equals
— Methodequals(
+true
sourceGeometryOps.equals
— Methodequals(
::GI.LinearRingTrait, l1,
::GI.LinearRingTrait, l2,
-)::Bool
Two linear rings are equal if they share the same set of points going along the curve. Note that rings are closed by definition, so they can have, but don't need, a repeated last point to be equal.
sourceGeometryOps.equals
— Methodequals(
+)::Bool
Two linear rings are equal if they share the same set of points going along the curve. Note that rings are closed by definition, so they can have, but don't need, a repeated last point to be equal.
sourceGeometryOps.equals
— Methodequals(
::GI.LinearRingTrait, l1,
::Union{GI.LineTrait, GI.LineStringTrait}, l2,
-)::Bool
A linear ring and a line/linestring are equal if they share the same set of points going along the curve. Note that lines aren't closed by defintion, but rings are, so the line must have a repeated last point to be equal
sourceGeometryOps.equals
— Methodequals(::GI.MultiPointTrait, mp1, ::GI.MultiPointTrait, mp2)::Bool
Two multipoints are equal if they share the same set of points.
sourceGeometryOps.equals
— Methodequals(::GI.MultiPointTrait, mp1, ::GI.PointTrait, p2)::Bool
A point and a multipoint are equal if the multipoint is composed of a single point that is equivalent to the given point.
sourceGeometryOps.equals
— Methodequals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool
Two multipolygons are equal if they share the same set of polygons.
sourceGeometryOps.equals
— Methodequals(::GI.MultiPolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool
A polygon and a multipolygon are equal if the multipolygon is composed of a single polygon that is equivalent to the given polygon.
sourceGeometryOps.equals
— Methodequals(::GI.PointTrait, p1, ::GI.MultiPointTrait, mp2)::Bool
A point and a multipoint are equal if the multipoint is composed of a single point that is equivalent to the given point.
sourceGeometryOps.equals
— Methodequals(::GI.PointTrait, p1, ::GI.PointTrait, p2)::Bool
Two points are the same if they have the same x and y (and z if 3D) coordinates.
sourceGeometryOps.equals
— Methodequals(::GI.PolygonTrait, geom_a, ::GI.MultiPolygonTrait, geom_b)::Bool
A polygon and a multipolygon are equal if the multipolygon is composed of a single polygon that is equivalent to the given polygon.
sourceGeometryOps.equals
— Methodequals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool
Two polygons are equal if they share the same exterior edge and holes.
sourceGeometryOps.equals
— Methodequals(
+)::Bool
A linear ring and a line/linestring are equal if they share the same set of points going along the curve. Note that lines aren't closed by defintion, but rings are, so the line must have a repeated last point to be equal
sourceGeometryOps.equals
— Methodequals(::GI.MultiPointTrait, mp1, ::GI.MultiPointTrait, mp2)::Bool
Two multipoints are equal if they share the same set of points.
sourceGeometryOps.equals
— Methodequals(::GI.MultiPointTrait, mp1, ::GI.PointTrait, p2)::Bool
A point and a multipoint are equal if the multipoint is composed of a single point that is equivalent to the given point.
sourceGeometryOps.equals
— Methodequals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool
Two multipolygons are equal if they share the same set of polygons.
sourceGeometryOps.equals
— Methodequals(::GI.MultiPolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool
A polygon and a multipolygon are equal if the multipolygon is composed of a single polygon that is equivalent to the given polygon.
sourceGeometryOps.equals
— Methodequals(::GI.PointTrait, p1, ::GI.MultiPointTrait, mp2)::Bool
A point and a multipoint are equal if the multipoint is composed of a single point that is equivalent to the given point.
sourceGeometryOps.equals
— Methodequals(::GI.PointTrait, p1, ::GI.PointTrait, p2)::Bool
Two points are the same if they have the same x and y (and z if 3D) coordinates.
sourceGeometryOps.equals
— Methodequals(::GI.PolygonTrait, geom_a, ::GI.MultiPolygonTrait, geom_b)::Bool
A polygon and a multipolygon are equal if the multipolygon is composed of a single polygon that is equivalent to the given polygon.
sourceGeometryOps.equals
— Methodequals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool
Two polygons are equal if they share the same exterior edge and holes.
sourceGeometryOps.equals
— Methodequals(
::Union{GI.LineTrait, GI.LineStringTrait}, l1,
::GI.LinearRingTrait, l2,
-)::Bool
A line/linestring and a linear ring are equal if they share the same set of points going along the curve. Note that lines aren't closed by defintion, but rings are, so the line must have a repeated last point to be equal
sourceGeometryOps.equals
— Methodequals(
+)::Bool
A line/linestring and a linear ring are equal if they share the same set of points going along the curve. Note that lines aren't closed by defintion, but rings are, so the line must have a repeated last point to be equal
sourceGeometryOps.equals
— Methodequals(
::Union{GI.LineTrait, GI.LineStringTrait}, l1,
::Union{GI.LineTrait, GI.LineStringTrait}, l2,
-)::Bool
Two lines/linestrings are equal if they share the same set of points going along the curve. Note that lines/linestrings aren't closed by defintion.
sourceGeometryOps.equals
— Methodequals(::T, geom_a, ::T, geom_b)::Bool
Two geometries of the same type, which don't have a equals function to dispatch off of should throw an error.
sourceGeometryOps.euclid_distance
— Methodeuclid_distance(p1::Point, p2::Point)::Real
Returns the Euclidean distance between two points.
sourceGeometryOps.flatten
— Methodflatten(target::Type{<:GI.AbstractTrait}, obj)
-flatten(f, target::Type{<:GI.AbstractTrait}, obj)
Lazily flatten any AbstractArray
, iterator, FeatureCollectionTrait
, FeatureTrait
or AbstractGeometryTrait
object obj
, so that objects with the target
trait are returned by the iterator.
If f
is passed in it will be applied to the target geometries.
sourceGeometryOps.flip
— Methodflip(obj)
Swap all of the x and y coordinates in obj, otherwise keeping the original structure (but not necessarily the original type).
Keywords
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
sourceGeometryOps.get_contours
— Methodget_contours(A::AbstractMatrix)
Returns contours as vectors of CartesianIndex
.
sourceGeometryOps.intersection
— Methodintersection(geom_a, geom_b)::Union{Tuple{::Real, ::Real}, ::Nothing}
Return an intersection point between two geometries. Return nothing if none are found. Else, the return type depends on the input. It will be a union between: a point, a line, a linear ring, a polygon, or a multipolygon
Example
import GeoInterface as GI, GeometryOps as GO
+)::Bool
Two lines/linestrings are equal if they share the same set of points going along the curve. Note that lines/linestrings aren't closed by defintion.
sourceGeometryOps.equals
— Methodequals(::T, geom_a, ::T, geom_b)::Bool
Two geometries of the same type, which don't have a equals function to dispatch off of should throw an error.
sourceGeometryOps.euclid_distance
— Methodeuclid_distance(p1::Point, p2::Point)::Real
Returns the Euclidean distance between two points.
sourceGeometryOps.flatten
— Methodflatten(target::Type{<:GI.AbstractTrait}, obj)
+flatten(f, target::Type{<:GI.AbstractTrait}, obj)
Lazily flatten any AbstractArray
, iterator, FeatureCollectionTrait
, FeatureTrait
or AbstractGeometryTrait
object obj
, so that objects with the target
trait are returned by the iterator.
If f
is passed in it will be applied to the target geometries.
sourceGeometryOps.flip
— Methodflip(obj)
Swap all of the x and y coordinates in obj, otherwise keeping the original structure (but not necessarily the original type).
Keywords
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
sourceGeometryOps.get_contours
— Methodget_contours(A::AbstractMatrix)
Returns contours as vectors of CartesianIndex
.
sourceGeometryOps.intersection
— Methodintersection(geom_a, geom_b)::Union{Tuple{::Real, ::Real}, ::Nothing}
Return an intersection point between two geometries. Return nothing if none are found. Else, the return type depends on the input. It will be a union between: a point, a line, a linear ring, a polygon, or a multipolygon
Example
import GeoInterface as GI, GeometryOps as GO
line1 = GI.Line([(124.584961,-12.768946), (126.738281,-17.224758)])
line2 = GI.Line([(123.354492,-15.961329), (127.22168,-14.008696)])
GO.intersection(line1, line2)
# output
-(125.58375366067547, -14.83572303404496)
sourceGeometryOps.intersection
— Methodintersection(
+(125.58375366067547, -14.83572303404496)
sourceGeometryOps.intersection
— Methodintersection(
::GI.AbstractTrait, geom_a,
::GI.AbstractTrait, geom_b,
)::Union{
::Vector{Vector{Tuple{::Real, ::Real}}}, # is this a good return type?
::Nothing
-}
Calculates the intersection between two line segments. Return nothing if there isn't one.
sourceGeometryOps.intersection
— Methodintersection(
+}
Calculates the intersection between two line segments. Return nothing if there isn't one.
sourceGeometryOps.intersection
— Methodintersection(
::GI.LineTrait, line_a,
::GI.LineTrait, line_b,
)::Union{
::Tuple{::Real, ::Real},
::Nothing
-}
Calculates the intersection between two line segments. Return nothing if there isn't one.
sourceGeometryOps.intersection
— Methodintersection(
+}
Calculates the intersection between two line segments. Return nothing if there isn't one.
sourceGeometryOps.intersection
— Methodintersection(
::GI.PolygonTrait, poly_a,
::GI.PolygonTrait, poly_b,
)::Union{
::Vector{Vector{Tuple{::Real, ::Real}}}, # is this a good return type?
::Nothing
-}
Calculates the intersection between two line segments. Return nothing if there isn't one.
sourceGeometryOps.intersection_points
— Methodintersection_points(
+}
Calculates the intersection between two line segments. Return nothing if there isn't one.
sourceGeometryOps.intersection_points
— Methodintersection_points(
geom_a,
geom_b,
)::Union{
::Vector{::Tuple{::Real, ::Real}},
::Nothing,
-}
Return a list of intersection points between two geometries. If no intersection point was possible given geometry extents, return nothing. If none are found, return an empty list.
sourceGeometryOps.intersection_points
— Methodintersection_points(
+}
Return a list of intersection points between two geometries. If no intersection point was possible given geometry extents, return nothing. If none are found, return an empty list.
sourceGeometryOps.intersection_points
— Methodintersection_points(
::GI.AbstractTrait, geom_a,
::GI.AbstractTrait, geom_b,
)::Union{
::Vector{::Tuple{::Real, ::Real}},
::Nothing,
-}
Calculates the list of intersection points between two geometries, inlcuding line segments, line strings, linear rings, polygons, and multipolygons. If no intersection points were possible given geometry extents, return nothing. If none are found, return an empty list.
sourceGeometryOps.intersects
— Methodintersects(geom1, geom2)::Bool
Check if two geometries intersect, returning true if so and false otherwise.
Example
import GeoInterface as GI, GeometryOps as GO
+}
Calculates the list of intersection points between two geometries, inlcuding line segments, line strings, linear rings, polygons, and multipolygons. If no intersection points were possible given geometry extents, return nothing. If none are found, return an empty list.
sourceGeometryOps.intersects
— Methodintersects(geom1, geom2)::Bool
Check if two geometries intersect, returning true if so and false otherwise.
Example
import GeoInterface as GI, GeometryOps as GO
line1 = GI.Line([(124.584961,-12.768946), (126.738281,-17.224758)])
line2 = GI.Line([(123.354492,-15.961329), (127.22168,-14.008696)])
GO.intersects(line1, line2)
# output
-true
sourceGeometryOps.intersects
— Methodintersects(::GI.AbstractTrait, a, ::GI.AbstractTrait, b)::Bool
Returns true if two geometries intersect with one another and false otherwise. For all geometries but lines, convert the geometry to a list of edges and cross compare the edges for intersections.
sourceGeometryOps.intersects
— Methodintersects(::GI.LineTrait, a, ::GI.LineTrait, b)::Bool
Returns true if two line segments intersect and false otherwise.
sourceGeometryOps.isclockwise
— Methodisclockwise(line::Union{LineString, Vector{Position}})::Bool
Take a ring and return true or false whether or not the ring is clockwise or counter-clockwise.
Example
import GeoInterface as GI, GeometryOps as GO
+true
sourceGeometryOps.intersects
— Methodintersects(::GI.AbstractTrait, a, ::GI.AbstractTrait, b)::Bool
Returns true if two geometries intersect with one another and false otherwise. For all geometries but lines, convert the geometry to a list of edges and cross compare the edges for intersections.
sourceGeometryOps.intersects
— Methodintersects(::GI.LineTrait, a, ::GI.LineTrait, b)::Bool
Returns true if two line segments intersect and false otherwise.
sourceGeometryOps.isclockwise
— Methodisclockwise(line::Union{LineString, Vector{Position}})::Bool
Take a ring and return true or false whether or not the ring is clockwise or counter-clockwise.
Example
import GeoInterface as GI, GeometryOps as GO
ring = GI.LinearRing([(0, 0), (1, 1), (1, 0), (0, 0)])
GO.isclockwise(ring)
# output
-true
sourceGeometryOps.isconcave
— Methodisconcave(poly::Polygon)::Bool
Take a polygon and return true or false as to whether it is concave or not.
Examples
import GeoInterface as GI, GeometryOps as GO
+true
sourceGeometryOps.isconcave
— Methodisconcave(poly::Polygon)::Bool
Take a polygon and return true or false as to whether it is concave or not.
Examples
import GeoInterface as GI, GeometryOps as GO
poly = GI.Polygon([[(0, 0), (0, 1), (1, 1), (1, 0), (0, 0)]])
GO.isconcave(poly)
# output
-false
sourceGeometryOps.overlaps
— Methodoverlaps(geom1, geom2)::Bool
Compare two Geometries of the same dimension and return true if their intersection set results in a geometry different from both but of the same dimension. This means one geometry cannot be within or contain the other and they cannot be equal
Examples
import GeometryOps as GO, GeoInterface as GI
+false
sourceGeometryOps.overlaps
— Methodoverlaps(geom1, geom2)::Bool
Compare two Geometries of the same dimension and return true if their intersection set results in a geometry different from both but of the same dimension. This means one geometry cannot be within or contain the other and they cannot be equal
Examples
import GeometryOps as GO, GeoInterface as GI
poly1 = GI.Polygon([[(0,0), (0,5), (5,5), (5,0), (0,0)]])
poly2 = GI.Polygon([[(1,1), (1,6), (6,6), (6,1), (1,1)]])
GO.overlaps(poly1, poly2)
# output
-true
sourceGeometryOps.overlaps
— Methodoverlaps(::GI.AbstractTrait, geom1, ::GI.AbstractTrait, geom2)::Bool
For any non-specified pair, all have non-matching dimensions, return false.
sourceGeometryOps.overlaps
— Methodoverlaps(::GI.LineTrait, line1, ::GI.LineTrait, line)::Bool
If the lines overlap, meaning that they are colinear but each have one endpoint outside of the other line, return true. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
+true
sourceGeometryOps.overlaps
— Methodoverlaps(::GI.AbstractTrait, geom1, ::GI.AbstractTrait, geom2)::Bool
For any non-specified pair, all have non-matching dimensions, return false.
sourceGeometryOps.overlaps
— Methodoverlaps(::GI.LineTrait, line1, ::GI.LineTrait, line)::Bool
If the lines overlap, meaning that they are colinear but each have one endpoint outside of the other line, return true. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
::GI.MultiPointTrait, points1,
::GI.MultiPointTrait, points2,
-)::Bool
If the multipoints overlap, meaning some, but not all, of the points within the multipoints are shared, return true.
sourceGeometryOps.overlaps
— Methodoverlaps(
+)::Bool
If the multipoints overlap, meaning some, but not all, of the points within the multipoints are shared, return true.
sourceGeometryOps.overlaps
— Methodoverlaps(
::GI.MultiPolygonTrait, polys1,
::GI.MultiPolygonTrait, polys2,
-)::Bool
Return true if at least one pair of polygons from multipolygons overlap. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
+)::Bool
Return true if at least one pair of polygons from multipolygons overlap. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
::GI.MultiPolygonTrait, polys1,
::GI.PolygonTrait, poly2,
-)::Bool
Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
+)::Bool
Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
::GI.PolygonTrait, poly1,
::GI.MultiPolygonTrait, polys2,
-)::Bool
Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
+)::Bool
Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
trait_a::GI.PolygonTrait, poly_a,
trait_b::GI.PolygonTrait, poly_b,
-)::Bool
If the two polygons intersect with one another, but are not equal, return true. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
+)::Bool
If the two polygons intersect with one another, but are not equal, return true. Else false.
sourceGeometryOps.overlaps
— Methodoverlaps(
::Union{GI.LineStringTrait, GI.LinearRing}, line1,
::Union{GI.LineStringTrait, GI.LinearRing}, line2,
-)::Bool
If the curves overlap, meaning that at least one edge of each curve overlaps, return true. Else false.
sourceGeometryOps.point_in_polygon
— Methodpoint_in_polygon(point::Point, polygon::Union{Polygon, MultiPolygon}, ignore_boundary::Bool=false)::Bool
Take a Point and a Polygon and determine if the point resides inside the polygon. The polygon can be convex or concave. The function accounts for holes.
Examples
import GeoInterface as GI, GeometryOps as GO
+)::Bool
If the curves overlap, meaning that at least one edge of each curve overlaps, return true. Else false.
sourceGeometryOps.point_in_polygon
— Methodpoint_in_polygon(point::Point, polygon::Union{Polygon, MultiPolygon}, ignore_boundary::Bool=false)::Bool
Take a Point and a Polygon and determine if the point resides inside the polygon. The polygon can be convex or concave. The function accounts for holes.
Examples
import GeoInterface as GI, GeometryOps as GO
point = (-77.0, 44.0)
poly = GI.Polygon([[(-81, 41), (-81, 47), (-72, 47), (-72, 41), (-81, 41)]])
GO.point_in_polygon(point, poly)
# output
-true
sourceGeometryOps.point_on_line
— Methodpoint_on_line(point::Point, line::LineString; ignore_end_vertices::Bool=false)::Bool
Return true if a point is on a line. Accept a optional parameter to ignore the start and end vertices of the linestring.
Examples
import GeoInterface as GI, GeometryOps as GO
+true
sourceGeometryOps.point_on_line
— Methodpoint_on_line(point::Point, line::LineString; ignore_end_vertices::Bool=false)::Bool
Return true if a point is on a line. Accept a optional parameter to ignore the start and end vertices of the linestring.
Examples
import GeoInterface as GI, GeometryOps as GO
point = (1, 1)
line = GI.LineString([(0, 0), (3, 3), (4, 4)])
GO.point_on_line(point, line)
# output
-true
sourceGeometryOps.polygon_to_line
— Methodpolygon_to_line(poly::Polygon)
Converts a Polygon to LineString or MultiLineString
Examples
import GeometryOps as GO, GeoInterface as GI
+true
sourceGeometryOps.polygon_to_line
— Methodpolygon_to_line(poly::Polygon)
Converts a Polygon to LineString or MultiLineString
Examples
import GeometryOps as GO, GeoInterface as GI
poly = GI.Polygon([[(-2.275543, 53.464547), (-2.275543, 53.489271), (-2.215118, 53.489271), (-2.215118, 53.464547), (-2.275543, 53.464547)]])
GO.polygon_to_line(poly)
# output
-GeoInterface.Wrappers.LineString{false, false, Vector{Tuple{Float64, Float64}}, Nothing, Nothing}([(-2.275543, 53.464547), (-2.275543, 53.489271), (-2.215118, 53.489271), (-2.215118, 53.464547), (-2.275543, 53.464547)], nothing, nothing)
sourceGeometryOps.polygonize
— Methodpolygonize(A; minpoints=10)
-polygonize(xs, ys, A; minpoints=10)
Convert matrix A
to polygons.
If xs
and ys
are passed in they are used as the pixel center points.
Keywords
minpoints
: ignore polygons with less than minpoints
points.
sourceGeometryOps.rebuild
— Methodrebuild(geom, child_geoms)
Rebuild a geometry from child geometries.
By default geometries will be rebuilt as a GeoInterface.Wrappers
geometry, but rebuild
can have methods added to it to dispatch on geometries from other packages and specify how to rebuild them.
(Maybe it should go into GeoInterface.jl)
sourceGeometryOps.reconstruct
— Methodreconstruct(geom, components)
Reconstruct geom
from an iterable of component objects that match its structure.
All objects in components
must have the same GeoInterface.trait
.
Ususally used in combination with flatten
.
sourceGeometryOps.reproject
— Methodreproject(geometry; source_crs, target_crs, transform, always_xy, time)
+GeoInterface.Wrappers.LineString{false, false, Vector{Tuple{Float64, Float64}}, Nothing, Nothing}([(-2.275543, 53.464547), (-2.275543, 53.489271), (-2.215118, 53.489271), (-2.215118, 53.464547), (-2.275543, 53.464547)], nothing, nothing)
sourceGeometryOps.polygonize
— Methodpolygonize(A; minpoints=10)
+polygonize(xs, ys, A; minpoints=10)
Convert matrix A
to polygons.
If xs
and ys
are passed in they are used as the pixel center points.
Keywords
minpoints
: ignore polygons with less than minpoints
points.
sourceGeometryOps.rebuild
— Methodrebuild(geom, child_geoms)
Rebuild a geometry from child geometries.
By default geometries will be rebuilt as a GeoInterface.Wrappers
geometry, but rebuild
can have methods added to it to dispatch on geometries from other packages and specify how to rebuild them.
(Maybe it should go into GeoInterface.jl)
sourceGeometryOps.reconstruct
— Methodreconstruct(geom, components)
Reconstruct geom
from an iterable of component objects that match its structure.
All objects in components
must have the same GeoInterface.trait
.
Ususally used in combination with flatten
.
sourceGeometryOps.reproject
— Methodreproject(geometry; source_crs, target_crs, transform, always_xy, time)
reproject(geometry, source_crs, target_crs; always_xy, time)
-reproject(geometry, transform; always_xy, time)
Reproject any GeoInterface.jl compatible geometry
from source_crs
to target_crs
.
The returned object will be constructed from GeoInterface.WrapperGeometry
geometries, wrapping views of a Vector{Proj.Point{D}}
, where D
is the dimension.
Arguments
geometry
: Any GeoInterface.jl compatible geometries.source_crs
: the source coordinate referece system, as a GeoFormatTypes.jl object or a string.target_crs
: the target coordinate referece system, as a GeoFormatTypes.jl object or a string.
If these a passed as keywords, transform
will take priority. Without it target_crs
is always needed, and source_crs
is needed if it is not retreivable from the geometry with GeoInterface.crs(geometry)
.
Keywords
always_xy
: force x, y coordinate order, true
by default. false
will expect and return points in the crs coordinate order.time
: the time for the coordinates. Inf
by default.threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
sourceGeometryOps.signed_area
— Methodsigned_area(geom)::Real
Returns the signed area of the geometry, based on winding order. This is computed slighly differently for different geometries: - The signed area of a point is always zero. - The signed area of a curve is always zero. - The signed area of a polygon is computed with the shoelace formula and is positive if the polygon coordinates wind clockwise and negative if counterclockwise. - You cannot compute the signed area of a multipolygon as it doesn't have a meaning as each sub-polygon could have a different winding order.
sourceGeometryOps.signed_distance
— Methodsigned_distance(point, geom)::Real
Calculates the signed distance from the geometry geom
to the given point. Points within geom
have a negative signed distance, and points outside of geom
have a positive signed distance. - The signed distance from a point to a point, line, linestring, or linear ring is equal to the distance between the two. - The signed distance from a point to a polygon is negative if the point is within the polygon and is positive otherwise. The value of the distance is the minimum distance from the point to an edge of the polygon. This includes edges created by holes. - The signed distance from a point to a mulitpolygon is negative if the point is within one of the polygons that make up the multipolygon and is positive otherwise. The value of the distance is the minimum distance from the point to an edge of the multipolygon. This includes edges created by holes of the polygons as well.
sourceGeometryOps.simplify
— Methodsimplify(obj; kw...)
+reproject(geometry, transform; always_xy, time)
Reproject any GeoInterface.jl compatible geometry
from source_crs
to target_crs
.
The returned object will be constructed from GeoInterface.WrapperGeometry
geometries, wrapping views of a Vector{Proj.Point{D}}
, where D
is the dimension.
Arguments
geometry
: Any GeoInterface.jl compatible geometries.source_crs
: the source coordinate referece system, as a GeoFormatTypes.jl object or a string.target_crs
: the target coordinate referece system, as a GeoFormatTypes.jl object or a string.
If these a passed as keywords, transform
will take priority. Without it target_crs
is always needed, and source_crs
is needed if it is not retreivable from the geometry with GeoInterface.crs(geometry)
.
Keywords
always_xy
: force x, y coordinate order, true
by default. false
will expect and return points in the crs coordinate order.time
: the time for the coordinates. Inf
by default.threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
sourceGeometryOps.signed_area
— Methodsigned_area(geom)::Real
Returns the signed area of the geometry, based on winding order. This is computed slighly differently for different geometries: - The signed area of a point is always zero. - The signed area of a curve is always zero. - The signed area of a polygon is computed with the shoelace formula and is positive if the polygon coordinates wind clockwise and negative if counterclockwise. - You cannot compute the signed area of a multipolygon as it doesn't have a meaning as each sub-polygon could have a different winding order.
sourceGeometryOps.signed_distance
— Methodsigned_distance(point, geom)::Real
Calculates the signed distance from the geometry geom
to the given point. Points within geom
have a negative signed distance, and points outside of geom
have a positive signed distance. - The signed distance from a point to a point, line, linestring, or linear ring is equal to the distance between the two. - The signed distance from a point to a polygon is negative if the point is within the polygon and is positive otherwise. The value of the distance is the minimum distance from the point to an edge of the polygon. This includes edges created by holes. - The signed distance from a point to a mulitpolygon is negative if the point is within one of the polygons that make up the multipolygon and is positive otherwise. The value of the distance is the minimum distance from the point to an edge of the multipolygon. This includes edges created by holes of the polygons as well.
sourceGeometryOps.simplify
— Methodsimplify(obj; kw...)
simplify(::SimplifyAlg, obj; kw...)
Simplify a geometry, feature, feature collection, or nested vectors or a table of these.
RadialDistance
, DouglasPeucker
, or VisvalingamWhyatt
algorithms are available, listed in order of increasing quality but decreaseing performance.
PoinTrait
and MultiPointTrait
are returned unchanged.
The default behaviour is simplify(DouglasPeucker(; kw...), obj)
. Pass in other SimplifyAlg
to use other algorithms.
Keywords
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
Keywords for DouglasPeucker are allowed when no algorithm is specified:
Keywords
ratio
: the fraction of points that should remain after simplify
. Useful as it will generalise for large collections of objects.number
: the number of points that should remain after simplify
. Less useful for large collections of mixed size objects.
Example
Simplify a polygon to have six points:
import GeoInterface as GI
import GeometryOps as GO
@@ -195,5 +195,5 @@
GI.npoint(simple)
# output
-6
sourceGeometryOps.t_value
— Methodt_value(sᵢ, sᵢ₊₁, rᵢ, rᵢ₊₁)
Returns the "T-value" as described in Hormann's presentation [HormannPresentation] on how to calculate the mean-value coordinate.
Here, sᵢ
is the vector from vertex vᵢ
to the point, and rᵢ
is the norm (length) of sᵢ
. s
must be Point
and r
must be real numbers.
\[tᵢ = \frac{\mathrm{det}\left(sᵢ, sᵢ₊₁\right)}{rᵢ * rᵢ₊₁ + sᵢ ⋅ sᵢ₊₁}\]
```
sourceGeometryOps.to_edges
— Methodto_edges()
Convert any geometry or collection of geometries into a flat vector of Tuple{Tuple{Float64,Float64},Tuple{Float64,Float64}}
edges.
sourceGeometryOps.tuples
— Methodtuples(obj)
Convert all points in obj
to Tuple
s, wherever the are nested.
Returns a similar object or collection of objects using GeoInterface.jl geometries wrapping Tuple
points.
Keywords
threaded
: true
or false
. Whether to use multithreading. Defaults to false
.crs
: The CRS to attach to geometries. Defaults to nothing
.calc_extent
: true
or false
. Whether to calculate the extent. Defaults to false
.
sourceGeometryOps.unwrap
— Functionunwrap(target::Type{<:AbstractTrait}, obj)
-unwrap(f, target::Type{<:AbstractTrait}, obj)
Unwrap the object newst to vectors, down to the target trait.
If f
is passed in it will be applied to the target geometries as they are found.
sourceGeometryOps.weighted_mean
— Methodweighted_mean(weight::Real, x1, x2)
Returns the weighted mean of x1
and x2
, where weight
is the weight of x1
.
Specifically, calculates x1 * weight + x2 * (1 - weight)
.
Note The idea for this method is that you can override this for custom types, like Color types, in extension modules.
source- HormannPresentationK. Hormann and N. Sukumar. Generalized Barycentric Coordinates in Computer Graphics and Computational Mechanics. Taylor & Fancis, CRC Press, 2017.
Settings
This document was generated with Documenter.jl version 0.27.24 on Tuesday 26 December 2023. Using Julia version 1.9.4.
+6
GeometryOps.t_value
— Methodt_value(sᵢ, sᵢ₊₁, rᵢ, rᵢ₊₁)
Returns the "T-value" as described in Hormann's presentation [HormannPresentation] on how to calculate the mean-value coordinate.
Here, sᵢ
is the vector from vertex vᵢ
to the point, and rᵢ
is the norm (length) of sᵢ
. s
must be Point
and r
must be real numbers.
\[tᵢ = \frac{\mathrm{det}\left(sᵢ, sᵢ₊₁\right)}{rᵢ * rᵢ₊₁ + sᵢ ⋅ sᵢ₊₁}\]
```
GeometryOps.to_edges
— Methodto_edges()
Convert any geometry or collection of geometries into a flat vector of Tuple{Tuple{Float64,Float64},Tuple{Float64,Float64}}
edges.
GeometryOps.tuples
— Methodtuples(obj)
Convert all points in obj
to Tuple
s, wherever the are nested.
Returns a similar object or collection of objects using GeoInterface.jl geometries wrapping Tuple
points.
Keywords
threaded
:true
orfalse
. Whether to use multithreading. Defaults tofalse
.crs
: The CRS to attach to geometries. Defaults tonothing
.calc_extent
:true
orfalse
. Whether to calculate the extent. Defaults tofalse
.
GeometryOps.unwrap
— Functionunwrap(target::Type{<:AbstractTrait}, obj)
+unwrap(f, target::Type{<:AbstractTrait}, obj)
Unwrap the object newst to vectors, down to the target trait.
If f
is passed in it will be applied to the target geometries as they are found.
GeometryOps.weighted_mean
— Methodweighted_mean(weight::Real, x1, x2)
Returns the weighted mean of x1
and x2
, where weight
is the weight of x1
.
Specifically, calculates x1 * weight + x2 * (1 - weight)
.
The idea for this method is that you can override this for custom types, like Color types, in extension modules.
- HormannPresentationK. Hormann and N. Sukumar. Generalized Barycentric Coordinates in Computer Graphics and Computational Mechanics. Taylor & Fancis, CRC Press, 2017.