Skip to content

Full GeometryOps API documentation

Warning

This page is still very much WIP!

Documentation for GeometryOps's full API (only for reference!).

apply and associated functions

# GeometryOps.applyFunction.
julia
apply(f, target::Union{TraitTarget, GI.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. These will always be GeoInterface geometries/feature/feature collections. But "deeper" objects may remain unchanged or be whatever GeoInterface compatible objects f returns.

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:

julia
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

source


# GeometryOps.applyreduceFunction.
julia
applyreduce(f, op, target::Union{TraitTarget, GI.AbstractTrait}, obj; threaded)

Apply function f to all objects with the target trait, and reduce the result with an op like +.

The order and grouping of application of op is not guaranteed.

If threaded==true threads will be used over arrays and iterables, feature collections and nested geometries.

source


# GeometryOps.reprojectFunction.
julia
reproject(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.

Tip

The Proj.jl package must be loaded for this method to work, since it is implemented in a package extension.

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.

source


# GeometryOps.transformFunction.
julia
transform(f, obj)

Apply a function f to all the points in obj.

Points will be passed to f as an SVector to allow using CoordinateTransformations.jl and Rotations.jl without hassle.

SVector is also a valid GeoInterface.jl point, so will work in all GeoInterface.jl methods.

Example

julia
julia> import GeoInterface as GI

julia> import GeometryOps as GO

julia> geom = GI.Polygon([GI.LinearRing([(1, 2), (3, 4), (5, 6), (1, 2)]), GI.LinearRing([(3, 4), (5, 6), (6, 7), (3, 4)])]);

julia> f = CoordinateTransformations.Translation(3.5, 1.5)
Translation(3.5, 1.5)

julia> GO.transform(f, geom)
GeoInterface.Wrappers.Polygon{false, false, Vector{GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}}, Nothing, Nothing}(GeoInterface.Wrappers.Linea
rRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}[GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}(StaticArraysCo
re.SVector{2, Float64}[[4.5, 3.5], [6.5, 5.5], [8.5, 7.5], [4.5, 3.5]], nothing, nothing), GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}(StaticA
rraysCore.SVector{2, Float64}[[6.5, 5.5], [8.5, 7.5], [9.5, 8.5], [6.5, 5.5]], nothing, nothing)], nothing, nothing)

With Rotations.jl you need to actuall multiply the Rotation by the SVector point, which is easy using an anonymous function.

julia
julia> using Rotations

julia> GO.transform(p -> one(RotMatrix{2}) * p, geom)
GeoInterface.Wrappers.Polygon{false, false, Vector{GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}}, Nothing, Nothing}(GeoInterface.Wrappers.LinearR
ing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}[GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}(StaticArraysCore.SVe
ctor{2, Int64}[[2, 1], [4, 3], [6, 5], [2, 1]], nothing, nothing), GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}(StaticArraysCore.SVector{2, Int64
}[[4, 3], [6, 5], [7, 6], [4, 3]], nothing, nothing)], nothing, nothing)

source


General geometry methods

OGC methods

# GeometryOps.containsFunction.
julia
contains(g1::AbstractGeometry, g2::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 (g2) must not intersect the exterior of the first (g1).

contains returns the exact opposite result of within.

Examples

julia
import GeometryOps as GO, GeoInterface as GI
line = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
point = GI.Point((1, 2))

GO.contains(line, point)
# output
true

source


# GeometryOps.coveredbyFunction.
julia
coveredby(g1, g2)::Bool

Return true if the first geometry is completely covered by the second geometry. The interior and boundary of the primary geometry (g1) must not intersect the exterior of the secondary geometry (g2).

Furthermore, coveredby returns the exact opposite result of covers. They are equivalent with the order of the arguments swapped.

Examples

julia
import GeometryOps as GO, GeoInterface as GI
p1 = GI.Point(0.0, 0.0)
p2 = GI.Point(1.0, 1.0)
l1 = GI.Line([p1, p2])

GO.coveredby(p1, l1)
# output
true

source


# GeometryOps.coversFunction.
julia
covers(g1::AbstractGeometry, g2::AbstractGeometry)::Bool

Return true if the first geometry is completely covers the second geometry, The exterior and boundary of the second geometry must not be outside of the interior and boundary of the first geometry. However, the interiors need not intersect.

covers returns the exact opposite result of coveredby.

Examples

julia
import GeometryOps as GO, GeoInterface as GI
l1 = GI.LineString([(1.0, 1.0), (1.0, 2.0), (1.0, 3.0), (1.0, 4.0)])
l2 = GI.LineString([(1.0, 1.0), (1.0, 2.0)])

GO.covers(l1, l2)
# output
true

source


# GeometryOps.crossesFunction.
julia
 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

julia
import GeoInterface as GI, GeometryOps as GO
# TODO: Add working example

source


# GeometryOps.disjointFunction.
julia
disjoint(geom1, geom2)::Bool

Return true if the first geometry is disjoint from the second geometry.

Return true if the first geometry is disjoint from the second geometry. The interiors and boundaries of both geometries must not intersect.

Examples

julia
import GeometryOps as GO, GeoInterface as GI

line = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
point = (2, 2)
GO.disjoint(point, line)

# output
true

source


# GeometryOps.intersectsFunction.
julia
intersects(geom1, geom2)::Bool

Return true if the interiors or boundaries of the two geometries interact.

intersects returns the exact opposite result of disjoint.

Example

julia
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

source


# GeometryOps.overlapsFunction.
julia
overlaps(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

julia
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

source

julia
overlaps(::GI.AbstractTrait, geom1, ::GI.AbstractTrait, geom2)::Bool

For any non-specified pair, all have non-matching dimensions, return false.

source

julia
overlaps(
    ::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.

source

julia
overlaps(::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.

source

julia
overlaps(
    ::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.

source

julia
overlaps(
    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.

source

julia
overlaps(
    ::GI.PolygonTrait, poly1,
    ::GI.MultiPolygonTrait, polys2,
)::Bool

Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.

source

julia
overlaps(
    ::GI.MultiPolygonTrait, polys1,
    ::GI.PolygonTrait, poly2,
)::Bool

Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.

source

julia
overlaps(
    ::GI.MultiPolygonTrait, polys1,
    ::GI.MultiPolygonTrait, polys2,
)::Bool

Return true if at least one pair of polygons from multipolygons overlap. Else false.

source


# GeometryOps.touchesFunction.
julia
touches(geom1, geom2)::Bool

Return true if the first geometry touches the second geometry. In other words, the two interiors cannot interact, but one of the geometries must have a boundary point that interacts with either the other geometies interior or boundary.

Examples

julia
import GeometryOps as GO, GeoInterface as GI

l1 = GI.Line([(0.0, 0.0), (1.0, 0.0)])
l2 = GI.Line([(1.0, 1.0), (1.0, -1.0)])

GO.touches(l1, l2)
# output
true

source


# GeometryOps.withinFunction.
julia
within(geom1, geom2)::Bool

Return true if the first geometry is completely within the second geometry. The interiors of both geometries must intersect and the interior and boundary of the primary geometry (geom1) must not intersect the exterior of the secondary geometry (geom2).

Furthermore, within returns the exact opposite result of contains.

Examples

julia
import GeometryOps as GO, GeoInterface as GI

line = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
point = (1, 2)
GO.within(point, line)

# output
true

source


Other general methods

# GeometryOps.equalsFunction.
julia
equals(geom1, geom2)::Bool

Compare two Geometries return true if they are the same geometry.

Examples

julia
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

source

julia
equals(::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.

source

julia
equals(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.

source

julia
equals(::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.

source

julia
equals(::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.

source

julia
equals(::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.

source

julia
equals(::GI.MultiPointTrait, mp1, ::GI.MultiPointTrait, mp2)::Bool

Two multipoints are equal if they share the same set of points.

source

julia
equals(
    ::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.

source

julia
equals(
    ::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

source

julia
equals(
    ::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

source

julia
equals(
    ::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.

source

julia
equals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool

Two polygons are equal if they share the same exterior edge and holes.

source

julia
equals(::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.

source

julia
equals(::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.

source

julia
equals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool

Two multipolygons are equal if they share the same set of polygons.

source


# GeometryOps.centroidFunction.
julia
centroid(geom, [T=Float64])::Tuple{T, T}

Returns the centroid of a given line segment, linear ring, polygon, or mutlipolygon.

source


# GeometryOps.distanceFunction.
julia
distance(point, geom, ::Type{T} = Float64)::T

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 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 multigeometry or a geometry collection is the minimum distance between the point and any of the sub-geometries.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.signed_distanceFunction.
julia
signed_distance(point, geom, ::Type{T} = Float64)::T

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 multigeometry or a geometry collection is the minimum signed distance between the point and any of the sub-geometries.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.areaFunction.
julia
area(geom, [T = Float64])::T

Returns the area of a geometry or collection of geometries. This is computed slightly differently for different geometries:

- The area of a point/multipoint is always zero.
- The area of a curve/multicurve 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.
- The area of a geometry collection, feature collection of array/iterable 
    is the sum of the areas of all of the sub-geometries.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.signed_areaFunction.
julia
signed_area(geom, [T = Float64])::T

Returns the signed area of a single 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.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.anglesFunction.
julia
angles(geom, ::Type{T} = Float64)

Returns the angles of a geometry or collection of geometries. This is computed differently for different geometries:

- The angles of a point is an empty vector.
- The angles of a single line segment is an empty vector.
- The angles of a linestring or linearring is a vector of angles formed by the curve.
- The angles of a polygin is a vector of vectors of angles formed by each ring.
- The angles of a multi-geometry collection is a vector of the angles of each of the
    sub-geometries as defined above.

Result will be a Vector, or nested set of vectors, of type T where an optional argument with a default value of Float64.

source


# GeometryOps.embed_extentFunction.
julia
embed_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.

source


Barycentric coordinates

# GeometryOps.barycentric_coordinatesFunction.
julia
barycentric_coordinates(method = MeanValue(), polygon, point)

Returns the barycentric coordinates of point in polygon using the barycentric coordinate method method.

source


# GeometryOps.barycentric_coordinates!Function.
julia
barycentric_coordinates!(λs::Vector{<: Real}, method::AbstractBarycentricCoordinateMethod, polygon, point)

Loads the barycentric coordinates of point in polygon into λs using the barycentric coordinate method method.

λs must be of the length of the polygon plus its holes.

Tip

Use this method to avoid excess allocations when you need to calculate barycentric coordinates for many points.

source


# GeometryOps.barycentric_interpolateFunction.
julia
barycentric_interpolate(method = MeanValue(), polygon, values::AbstractVector{V}, point)

Returns the interpolated value at point within polygon using the barycentric coordinate method method. values are the per-point values for the polygon which are to be interpolated.

Returns an object of type V.

Warning

Barycentric interpolation is currently defined only for 2-dimensional polygons. If you pass a 3-D polygon in, the Z coordinate will be used as per-vertex value to be interpolated (the M coordinate in GIS parlance).

source


Other methods

# GeometryOps.AbstractBarycentricCoordinateMethodType.
julia
abstract 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.

source


# GeometryOps.ClosedRingType.
julia
ClosedRing() <: GeometryCorrection

This correction ensures that a polygon's exterior and interior rings are closed.

It can be called on any geometry correction as usual.

See also GeometryCorrection.

source


# GeometryOps.DiffIntersectingPolygonsType.
julia
DiffIntersectingPolygons() <: GeometryCorrection

This correction ensures that the polygons included in a multipolygon aren't intersecting. If any polygon's are intersecting, they will be made nonintersecting through the difference operation to create a unique set of disjoint (other than potentially connections by a single point) polygons covering the same area. See also GeometryCorrection, UnionIntersectingPolygons.

source


# GeometryOps.DouglasPeuckerType.
julia
DouglasPeucker <: 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.

Note: user input tol is squared to avoid uneccesary computation in algorithm.

source


# GeometryOps.GEOSType.
julia
GEOS(; params...)

A struct which instructs the method it's passed to as an algorithm to use the appropriate GEOS function via LibGEOS.jl for the operation.

Dispatch is generally carried out using the names of the keyword arguments. For example, segmentize will only accept a GEOS struct with only a max_distance keyword, and no other.

It's generally a lot slower than the native Julia implementations, since it must convert to the LibGEOS implementation and back - so be warned!

source


# GeometryOps.GeodesicSegmentsType.
julia
GeodesicSegments(; max_distance::Real, equatorial_radius::Real=6378137, flattening::Real=1/298.257223563)

A method for segmentizing geometries by adding extra vertices to the geometry so that no segment is longer than a given distance. This method calculates the distance between points on the geodesic, and assumes input in lat/long coordinates.

Warning

Any input geometries must be in lon/lat coordinates! If not, the method may fail or error.

Arguments

  • max_distance::Real: The maximum distance, in meters, between vertices in the geometry.

  • equatorial_radius::Real=6378137: The equatorial radius of the Earth, in meters. Passed to Proj.geod_geodesic.

  • flattening::Real=1/298.257223563: The flattening of the Earth, which is the ratio of the difference between the equatorial and polar radii to the equatorial radius. Passed to Proj.geod_geodesic.

One can also omit the equatorial_radius and flattening keyword arguments, and pass a geodesic object directly to the eponymous keyword.

This method uses the Proj/GeographicLib API for geodesic calculations.

source


# GeometryOps.GeometryCorrectionType.
julia
abstract type GeometryCorrection

This abstract type represents a geometry correction.

Interface

Any GeometryCorrection must implement two functions: * application_level(::GeometryCorrection)::AbstractGeometryTrait: This function should return the GeoInterface trait that the correction is intended to be applied to, like PointTrait or LineStringTrait or PolygonTrait. * (::GeometryCorrection)(::AbstractGeometryTrait, geometry)::(some_geometry): This function should apply the correction to the given geometry, and return a new geometry.

source


# GeometryOps.LineOrientationType.
julia
Enum LineOrientation

Enum for the orientation of a line with respect to a curve. A line can be line_cross (crossing over the curve), line_hinge (crossing the endpoint of the curve), line_over (colinear with the curve), or line_out (not interacting with the curve).

source


# GeometryOps.LinearSegmentsType.
julia
LinearSegments(; max_distance::Real)

A method for segmentizing geometries by adding extra vertices to the geometry so that no segment is longer than a given distance.

Here, max_distance is a purely nondimensional quantity and will apply in the input space. This is to say, that if the polygon is provided in lat/lon coordinates then the max_distance will be in degrees of arc. If the polygon is provided in meters, then the max_distance will be in meters.

source


# GeometryOps.MeanValueType.
julia
MeanValue() <: AbstractBarycentricCoordinateMethod

This method calculates barycentric coordinates using the mean value method.

References

source


# GeometryOps.PointOrientationType.
julia
Enum PointOrientation

Enum for the orientation of a point with respect to a curve. A point can be point_in the curve, point_on the curve, or point_out of the curve.

source


# GeometryOps.RadialDistanceType.
julia
RadialDistance <: 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.

Note: user input tol is squared to avoid uneccesary computation in algorithm.

source


# GeometryOps.SimplifyAlgType.
julia
abstract 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.

source


# GeometryOps.TraitTargetType.
julia
TraitTarget{T}

This struct holds a trait parameter or a union of trait parameters.

It is primarily used for dispatch into methods which select trait levels, like apply, or as a parameter to target.

Constructors

julia
TraitTarget(GI.PointTrait())
TraitTarget(GI.LineStringTrait(), GI.LinearRingTrait()) # and other traits as you may like
TraitTarget(TraitTarget(...))
# There are also type based constructors available, but that's not advised.
TraitTarget(GI.PointTrait)
TraitTarget(Union{GI.LineStringTrait, GI.LinearRingTrait})
# etc.

source


# GeometryOps.UnionIntersectingPolygonsType.
julia
UnionIntersectingPolygons() <: GeometryCorrection

This correction ensures that the polygon's included in a multipolygon aren't intersecting. If any polygon's are intersecting, they will be combined through the union operation to create a unique set of disjoint (other than potentially connections by a single point) polygons covering the same area.

See also GeometryCorrection.

source


# GeometryOps.VisvalingamWhyattType.
julia
VisvalingamWhyatt <: 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.

Note: user input tol is doubled to avoid uneccesary computation in algorithm.

source


# GeometryOps._detMethod.
julia
_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:

julia
s1[1] * s2[2] - s1[2] * s2[1]

source


# GeometryOps._equals_curvesMethod.
julia
_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 closed_type1::Bool true if c1 is closed by definition (polygon, linear ring) closed_type2::Bool true if c2 is closed by definition (polygon, linear ring)

source


# GeometryOps.anglesMethod.
julia
angles(geom, ::Type{T} = Float64)

Returns the angles of a geometry or collection of geometries. This is computed differently for different geometries:

- The angles of a point is an empty vector.
- The angles of a single line segment is an empty vector.
- The angles of a linestring or linearring is a vector of angles formed by the curve.
- The angles of a polygin is a vector of vectors of angles formed by each ring.
- The angles of a multi-geometry collection is a vector of the angles of each of the
    sub-geometries as defined above.

Result will be a Vector, or nested set of vectors, of type T where an optional argument with a default value of Float64.

source


# GeometryOps.applyMethod.
julia
apply(f, target::Union{TraitTarget, GI.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. These will always be GeoInterface geometries/feature/feature collections. But "deeper" objects may remain unchanged or be whatever GeoInterface compatible objects f returns.

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:

julia
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

source


# GeometryOps.applyreduceMethod.
julia
applyreduce(f, op, target::Union{TraitTarget, GI.AbstractTrait}, obj; threaded)

Apply function f to all objects with the target trait, and reduce the result with an op like +.

The order and grouping of application of op is not guaranteed.

If threaded==true threads will be used over arrays and iterables, feature collections and nested geometries.

source


# GeometryOps.areaMethod.
julia
area(geom, [T = Float64])::T

Returns the area of a geometry or collection of geometries. This is computed slightly differently for different geometries:

- The area of a point/multipoint is always zero.
- The area of a curve/multicurve 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.
- The area of a geometry collection, feature collection of array/iterable 
    is the sum of the areas of all of the sub-geometries.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.barycentric_coordinates!Method.
julia
barycentric_coordinates!(λs::Vector{<: Real}, method::AbstractBarycentricCoordinateMethod, polygon, point)

Loads the barycentric coordinates of point in polygon into λs using the barycentric coordinate method method.

λs must be of the length of the polygon plus its holes.

Tip

Use this method to avoid excess allocations when you need to calculate barycentric coordinates for many points.

source


# GeometryOps.barycentric_coordinatesMethod.
julia
barycentric_coordinates(method = MeanValue(), polygon, point)

Returns the barycentric coordinates of point in polygon using the barycentric coordinate method method.

source


# GeometryOps.barycentric_interpolateMethod.
julia
barycentric_interpolate(method = MeanValue(), polygon, values::AbstractVector{V}, point)

Returns the interpolated value at point within polygon using the barycentric coordinate method method. values are the per-point values for the polygon which are to be interpolated.

Returns an object of type V.

Warning

Barycentric interpolation is currently defined only for 2-dimensional polygons. If you pass a 3-D polygon in, the Z coordinate will be used as per-vertex value to be interpolated (the M coordinate in GIS parlance).

source


# GeometryOps.centroidMethod.
julia
centroid(geom, [T=Float64])::Tuple{T, T}

Returns the centroid of a given line segment, linear ring, polygon, or mutlipolygon.

source


# GeometryOps.centroid_and_areaMethod.
julia
centroid_and_area(geom, [T=Float64])::(::Tuple{T, T}, ::Real)

Returns the centroid and area of a given geometry.

source


# GeometryOps.centroid_and_lengthMethod.
julia
centroid_and_length(geom, [T=Float64])::(::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.

source


# GeometryOps.containsMethod.
julia
contains(g1::AbstractGeometry, g2::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 (g2) must not intersect the exterior of the first (g1).

contains returns the exact opposite result of within.

Examples

julia
import GeometryOps as GO, GeoInterface as GI
line = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
point = GI.Point((1, 2))

GO.contains(line, point)
# output
true

source


# GeometryOps.coverageMethod.
julia
coverage(geom, xmin, xmax, ymin, ymax, [T = Float64])::T

Returns the area of intersection between given geometry and grid cell defined by its minimum and maximum x and y-values. This is computed 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 calculated by tracing along its edges and switching to the cell edges if needed.

  • The coverage of a geometry collection, multi-geometry, feature collection of array/iterable is the sum of the coverages of all of the sub-geometries.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.coveredbyMethod.
julia
coveredby(g1, g2)::Bool

Return true if the first geometry is completely covered by the second geometry. The interior and boundary of the primary geometry (g1) must not intersect the exterior of the secondary geometry (g2).

Furthermore, coveredby returns the exact opposite result of covers. They are equivalent with the order of the arguments swapped.

Examples

julia
import GeometryOps as GO, GeoInterface as GI
p1 = GI.Point(0.0, 0.0)
p2 = GI.Point(1.0, 1.0)
l1 = GI.Line([p1, p2])

GO.coveredby(p1, l1)
# output
true

source


# GeometryOps.coversMethod.
julia
covers(g1::AbstractGeometry, g2::AbstractGeometry)::Bool

Return true if the first geometry is completely covers the second geometry, The exterior and boundary of the second geometry must not be outside of the interior and boundary of the first geometry. However, the interiors need not intersect.

covers returns the exact opposite result of coveredby.

Examples

julia
import GeometryOps as GO, GeoInterface as GI
l1 = GI.LineString([(1.0, 1.0), (1.0, 2.0), (1.0, 3.0), (1.0, 4.0)])
l2 = GI.LineString([(1.0, 1.0), (1.0, 2.0)])

GO.covers(l1, l2)
# output
true

source


# GeometryOps.crossesMethod.
julia
 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

julia
import GeoInterface as GI, GeometryOps as GO
# TODO: Add working example

source


# GeometryOps.cutMethod.
julia
cut(geom, line, [T::Type])

Return given geom cut by given line as a list of geometries of the same type as the input geom. Return the original geometry as only list element if none are found. Line must cut fully through given geometry or the original geometry will be returned.

Note: This currently doesn't work for degenerate cases there line crosses through vertices.

Example

julia
import GeoInterface as GI, GeometryOps as GO

poly = GI.Polygon([[(0.0, 0.0), (10.0, 0.0), (10.0, 10.0), (0.0, 10.0), (0.0, 0.0)]])
line = GI.Line([(5.0, -5.0), (5.0, 15.0)])
cut_polys = GO.cut(poly, line)
GI.coordinates.(cut_polys)

# output
2-element Vector{Vector{Vector{Vector{Float64}}}}:
 [[[0.0, 0.0], [5.0, 0.0], [5.0, 10.0], [0.0, 10.0], [0.0, 0.0]]]
 [[[5.0, 0.0], [10.0, 0.0], [10.0, 10.0], [5.0, 10.0], [5.0, 0.0]]]

source


# GeometryOps.differenceMethod.
julia
difference(geom_a, geom_b, [T::Type]; target::Type, fix_multipoly = UnionIntersectingPolygons())

Return the difference between two geometries as a list of geometries. Return an empty list if none are found. The type of the list will be constrained as much as possible given the input geometries. Furthermore, the user can provide a taget type as a keyword argument and a list of target geometries found in the difference will be returned. The user can also provide a float type that they would like the points of returned geometries to be. If the user is taking a intersection involving one or more multipolygons, and the multipolygon might be comprised of polygons that intersect, if fix_multipoly is set to an IntersectingPolygons correction (the default is UnionIntersectingPolygons()), then the needed multipolygons will be fixed to be valid before performing the intersection to ensure a correct answer. Only set fix_multipoly to false if you know that the multipolygons are valid, as it will avoid unneeded computation.

Example

julia
import GeoInterface as GI, GeometryOps as GO

poly1 = GI.Polygon([[[0.0, 0.0], [5.0, 5.0], [10.0, 0.0], [5.0, -5.0], [0.0, 0.0]]])
poly2 = GI.Polygon([[[3.0, 0.0], [8.0, 5.0], [13.0, 0.0], [8.0, -5.0], [3.0, 0.0]]])
diff_poly = GO.difference(poly1, poly2; target = GI.PolygonTrait())
GI.coordinates.(diff_poly)

# output
1-element Vector{Vector{Vector{Vector{Float64}}}}:
 [[[6.5, 3.5], [5.0, 5.0], [0.0, 0.0], [5.0, -5.0], [6.5, -3.5], [3.0, 0.0], [6.5, 3.5]]]

source


# GeometryOps.disjointMethod.
julia
disjoint(geom1, geom2)::Bool

Return true if the first geometry is disjoint from the second geometry.

Return true if the first geometry is disjoint from the second geometry. The interiors and boundaries of both geometries must not intersect.

Examples

julia
import GeometryOps as GO, GeoInterface as GI

line = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
point = (2, 2)
GO.disjoint(point, line)

# output
true

source


# GeometryOps.distanceMethod.
julia
distance(point, geom, ::Type{T} = Float64)::T

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 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 multigeometry or a geometry collection is the minimum distance between the point and any of the sub-geometries.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.embed_extentMethod.
julia
embed_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.

source


# GeometryOps.enforceMethod.
julia
enforce(alg::GO.GEOS, kw::Symbol, f)

Enforce the presence of a keyword argument in a GEOS algorithm, and return alg.params[kw].

Throws an error if the key is not present, and mentions f in the error message (since there isn't a good way to get the name of the function that called this method).

source


# GeometryOps.equalsMethod.
julia
equals(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.

source


# GeometryOps.equalsMethod.
julia
equals(geom1, geom2)::Bool

Compare two Geometries return true if they are the same geometry.

Examples

julia
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

source


# GeometryOps.equalsMethod.
julia
equals(
    ::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.

source


# GeometryOps.equalsMethod.
julia
equals(
    ::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

source


# GeometryOps.equalsMethod.
julia
equals(::GI.MultiPointTrait, mp1, ::GI.MultiPointTrait, mp2)::Bool

Two multipoints are equal if they share the same set of points.

source


# GeometryOps.equalsMethod.
julia
equals(::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.

source


# GeometryOps.equalsMethod.
julia
equals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool

Two multipolygons are equal if they share the same set of polygons.

source


# GeometryOps.equalsMethod.
julia
equals(::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.

source


# GeometryOps.equalsMethod.
julia
equals(::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.

source


# GeometryOps.equalsMethod.
julia
equals(::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.

source


# GeometryOps.equalsMethod.
julia
equals(::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.

source


# GeometryOps.equalsMethod.
julia
equals(::GI.PolygonTrait, geom_a, ::GI.PolygonTrait, geom_b)::Bool

Two polygons are equal if they share the same exterior edge and holes.

source


# GeometryOps.equalsMethod.
julia
equals(
    ::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

source


# GeometryOps.equalsMethod.
julia
equals(
    ::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.

source


# GeometryOps.equalsMethod.
julia
equals(::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.

source


# GeometryOps.flattenMethod.
julia
flatten(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.

source


# GeometryOps.flipMethod.
julia
flip(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.

source


# GeometryOps.intersectionMethod.
julia
intersection(geom_a, geom_b, [T::Type]; target::Type, fix_multipoly = UnionIntersectingPolygons())

Return the intersection between two geometries as a list of geometries. Return an empty list if none are found. The type of the list will be constrained as much as possible given the input geometries. Furthermore, the user can provide a target type as a keyword argument and a list of target geometries found in the intersection will be returned. The user can also provide a float type that they would like the points of returned geometries to be. If the user is taking a intersection involving one or more multipolygons, and the multipolygon might be comprised of polygons that intersect, if fix_multipoly is set to an IntersectingPolygons correction (the default is UnionIntersectingPolygons()), then the needed multipolygons will be fixed to be valid before performing the intersection to ensure a correct answer. Only set fix_multipoly to nothing if you know that the multipolygons are valid, as it will avoid unneeded computation.

Example

julia
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)])
inter_points = GO.intersection(line1, line2; target = GI.PointTrait())
GI.coordinates.(inter_points)

# output
1-element Vector{Vector{Float64}}:
 [125.58375366067548, -14.83572303404496]

source


# GeometryOps.intersection_pointsMethod.
julia
intersection_points(
    geom_a,
    geom_b,
)::Union{
    ::Vector{::Tuple{::Real, ::Real}},
    ::Nothing,
}

Return a list of intersection points between two geometries of type GI.Point. If no intersection point was possible given geometry extents, returns an empty list.

source


# GeometryOps.intersectsMethod.
julia
intersects(geom1, geom2)::Bool

Return true if the interiors or boundaries of the two geometries interact.

intersects returns the exact opposite result of disjoint.

Example

julia
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

source


# GeometryOps.isclockwiseMethod.
julia
isclockwise(line::Union{LineString, Vector{Position}})::Bool

Take a ring and return true if the line goes clockwise, or false if the line goes counter-clockwise. "Going clockwise" means, mathematically,

(i=2n(xixi1)(yi+yi1))>0

Example

julia
julia> import GeoInterface as GI, GeometryOps as GO
julia> ring = GI.LinearRing([(0, 0), (1, 1), (1, 0), (0, 0)]);
julia> GO.isclockwise(ring)
# output
true

source


# GeometryOps.isconcaveMethod.
julia
isconcave(poly::Polygon)::Bool

Take a polygon and return true or false as to whether it is concave or not.

Examples

julia
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

source


# GeometryOps.overlapsMethod.
julia
overlaps(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

julia
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

source


# GeometryOps.overlapsMethod.
julia
overlaps(::GI.AbstractTrait, geom1, ::GI.AbstractTrait, geom2)::Bool

For any non-specified pair, all have non-matching dimensions, return false.

source


# GeometryOps.overlapsMethod.
julia
overlaps(::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.

source


# GeometryOps.overlapsMethod.
julia
overlaps(
    ::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.

source


# GeometryOps.overlapsMethod.
julia
overlaps(
    ::GI.MultiPolygonTrait, polys1,
    ::GI.MultiPolygonTrait, polys2,
)::Bool

Return true if at least one pair of polygons from multipolygons overlap. Else false.

source


# GeometryOps.overlapsMethod.
julia
overlaps(
    ::GI.MultiPolygonTrait, polys1,
    ::GI.PolygonTrait, poly2,
)::Bool

Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.

source


# GeometryOps.overlapsMethod.
julia
overlaps(
    ::GI.PolygonTrait, poly1,
    ::GI.MultiPolygonTrait, polys2,
)::Bool

Return true if polygon overlaps with at least one of the polygons within the multipolygon. Else false.

source


# GeometryOps.overlapsMethod.
julia
overlaps(
    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.

source


# GeometryOps.overlapsMethod.
julia
overlaps(
    ::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.

source


# GeometryOps.polygon_to_lineMethod.
julia
polygon_to_line(poly::Polygon)

Converts a Polygon to LineString or MultiLineString

Examples

julia
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)

source


# GeometryOps.polygonizeMethod.
julia
polygonize(A::AbstractMatrix{Bool}; kw...)
polygonize(f, A::AbstractMatrix; kw...)
polygonize(xs, ys, A::AbstractMatrix{Bool}; kw...)
polygonize(f, xs, ys, A::AbstractMatrix; kw...)

Polygonize an AbstractMatrix of values, currently to a single class of polygons.

Returns a MultiPolygon for Bool values and f return values, and a FeatureCollection of Features holding MultiPolygon for all other values.

Function f should return either true or false or a transformation of values into simpler groups, especially useful for floating point arrays.

If xs and ys are ranges, they are used as the pixel/cell center points. If they are Vector of Tuple they are used as the lower and upper bounds of each pixel/cell.

Keywords

  • minpoints: ignore polygons with less than minpoints points.

  • values: the values to turn into polygons. By default these are union(A), If function f is passed these refer to the return values of f, by default union(map(f, A). If values Bool, false is ignored and a single MultiPolygon is returned rather than a FeatureCollection.

Example

julia
using GeometryOps
A = rand(100, 100)
multipolygon = polygonize(>(0.5), A);

source


# GeometryOps.rebuildMethod.
julia
rebuild(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)

source


# GeometryOps.reconstructMethod.
julia
reconstruct(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.

source


# GeometryOps.segmentizeMethod.
julia
segmentize([method = LinearSegments()], geom; max_distance::Real, threaded)

Segmentize a geometry by adding extra vertices to the geometry so that no segment is longer than a given distance. This is useful for plotting geometries with a limited number of vertices, or for ensuring that a geometry is not too "coarse" for a given application.

Arguments

  • method::SegmentizeMethod = LinearSegments(): The method to use for segmentizing the geometry. At the moment, only LinearSegments and GeodesicSegments are available.

  • geom: The geometry to segmentize. Must be a LineString, LinearRing, or greater in complexity.

  • max_distance::Real: The maximum distance, in the input space, between vertices in the geometry. Only used if you don't explicitly pass a method.

Returns a geometry of similar type to the input geometry, but resampled.

source


# GeometryOps.signed_areaMethod.
julia
signed_area(geom, [T = Float64])::T

Returns the signed area of a single 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.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.signed_distanceMethod.
julia
signed_distance(point, geom, ::Type{T} = Float64)::T

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 multigeometry or a geometry collection is the minimum signed distance between the point and any of the sub-geometries.

Result will be of type T, where T is an optional argument with a default value of Float64.

source


# GeometryOps.simplifyMethod.
julia
simplify(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

  • prefilter_alg: SimplifyAlg algorithm used to pre-filter object before using primary filtering algorithm.

  • 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.

  • tol: the minimum distance a point will be from the line joining its neighboring points.

Example

Simplify a polygon to have six points:

julia
import GeoInterface as GI
import GeometryOps as GO

poly = GI.Polygon([[
    [-70.603637, -33.399918],
    [-70.614624, -33.395332],
    [-70.639343, -33.392466],
    [-70.659942, -33.394759],
    [-70.683975, -33.404504],
    [-70.697021, -33.419406],
    [-70.701141, -33.434306],
    [-70.700454, -33.446339],
    [-70.694274, -33.458369],
    [-70.682601, -33.465816],
    [-70.668869, -33.472117],
    [-70.646209, -33.473835],
    [-70.624923, -33.472117],
    [-70.609817, -33.468107],
    [-70.595397, -33.458369],
    [-70.587158, -33.442901],
    [-70.587158, -33.426283],
    [-70.590591, -33.414248],
    [-70.594711, -33.406224],
    [-70.603637, -33.399918]]])

simple = GO.simplify(poly; number=6)
GI.npoint(simple)

# output
6

source


# GeometryOps.t_valueMethod.
julia
t_value(sᵢ, sᵢ₊₁, rᵢ, rᵢ₊₁)

Returns the "T-value" as described in Hormann's presentation [1] 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=det(s,s)rr+ss


[source](https://github.com/JuliaGeo/GeometryOps.jl/blob/v0.1.9/src/methods/barycentric.jl#L289-L305)

</div>
<br>
<div style='border-width:1px; border-style:solid; border-color:black; padding: 1em; border-radius: 25px;'>
<a id='GeometryOps.to_edges-Union{Tuple{Any}, Tuple{T}, Tuple{Any, Type{T}}} where T' href='#GeometryOps.to_edges-Union{Tuple{Any}, Tuple{T}, Tuple{Any, Type{T}}} where T'>#</a>&nbsp;<b><u>GeometryOps.to_edges</u></b> &mdash; <i>Method</i>.




```julia
to_edges()

Convert any geometry or collection of geometries into a flat vector of Tuple{Tuple{Float64,Float64},Tuple{Float64,Float64}} edges.

source


# GeometryOps.touchesMethod.
julia
touches(geom1, geom2)::Bool

Return true if the first geometry touches the second geometry. In other words, the two interiors cannot interact, but one of the geometries must have a boundary point that interacts with either the other geometies interior or boundary.

Examples

julia
import GeometryOps as GO, GeoInterface as GI

l1 = GI.Line([(0.0, 0.0), (1.0, 0.0)])
l2 = GI.Line([(1.0, 1.0), (1.0, -1.0)])

GO.touches(l1, l2)
# output
true

source


# GeometryOps.transformMethod.
julia
transform(f, obj)

Apply a function f to all the points in obj.

Points will be passed to f as an SVector to allow using CoordinateTransformations.jl and Rotations.jl without hassle.

SVector is also a valid GeoInterface.jl point, so will work in all GeoInterface.jl methods.

Example

julia
julia> import GeoInterface as GI

julia> import GeometryOps as GO

julia> geom = GI.Polygon([GI.LinearRing([(1, 2), (3, 4), (5, 6), (1, 2)]), GI.LinearRing([(3, 4), (5, 6), (6, 7), (3, 4)])]);

julia> f = CoordinateTransformations.Translation(3.5, 1.5)
Translation(3.5, 1.5)

julia> GO.transform(f, geom)
GeoInterface.Wrappers.Polygon{false, false, Vector{GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}}, Nothing, Nothing}(GeoInterface.Wrappers.Linea
rRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}[GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}(StaticArraysCo
re.SVector{2, Float64}[[4.5, 3.5], [6.5, 5.5], [8.5, 7.5], [4.5, 3.5]], nothing, nothing), GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Float64}}, Nothing, Nothing}(StaticA
rraysCore.SVector{2, Float64}[[6.5, 5.5], [8.5, 7.5], [9.5, 8.5], [6.5, 5.5]], nothing, nothing)], nothing, nothing)

With Rotations.jl you need to actuall multiply the Rotation by the SVector point, which is easy using an anonymous function.

julia
julia> using Rotations

julia> GO.transform(p -> one(RotMatrix{2}) * p, geom)
GeoInterface.Wrappers.Polygon{false, false, Vector{GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}}, Nothing, Nothing}(GeoInterface.Wrappers.LinearR
ing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}[GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}(StaticArraysCore.SVe
ctor{2, Int64}[[2, 1], [4, 3], [6, 5], [2, 1]], nothing, nothing), GeoInterface.Wrappers.LinearRing{false, false, Vector{StaticArraysCore.SVector{2, Int64}}, Nothing, Nothing}(StaticArraysCore.SVector{2, Int64
}[[4, 3], [6, 5], [7, 6], [4, 3]], nothing, nothing)], nothing, nothing)

source


# GeometryOps.tuplesMethod.
julia
tuples(obj)

Convert all points in obj to Tuples, 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.

source


# GeometryOps.unionMethod.
julia
union(geom_a, geom_b, [::Type{T}]; target::Type, fix_multipoly = UnionIntersectingPolygons())

Return the union between two geometries as a list of geometries. Return an empty list if none are found. The type of the list will be constrained as much as possible given the input geometries. Furthermore, the user can provide a taget type as a keyword argument and a list of target geometries found in the difference will be returned. The user can also provide a float type 'T' that they would like the points of returned geometries to be. If the user is taking a intersection involving one or more multipolygons, and the multipolygon might be comprised of polygons that intersect, if fix_multipoly is set to an IntersectingPolygons correction (the default is UnionIntersectingPolygons()), then the needed multipolygons will be fixed to be valid before performing the intersection to ensure a correct answer. Only set fix_multipoly to false if you know that the multipolygons are valid, as it will avoid unneeded computation.

Calculates the union between two polygons.

Example

julia
import GeoInterface as GI, GeometryOps as GO

p1 = GI.Polygon([[(0.0, 0.0), (5.0, 5.0), (10.0, 0.0), (5.0, -5.0), (0.0, 0.0)]])
p2 = GI.Polygon([[(3.0, 0.0), (8.0, 5.0), (13.0, 0.0), (8.0, -5.0), (3.0, 0.0)]])
union_poly = GO.union(p1, p2; target = GI.PolygonTrait())
GI.coordinates.(union_poly)

# output
1-element Vector{Vector{Vector{Vector{Float64}}}}:
 [[[6.5, 3.5], [5.0, 5.0], [0.0, 0.0], [5.0, -5.0], [6.5, -3.5], [8.0, -5.0], [13.0, 0.0], [8.0, 5.0], [6.5, 3.5]]]

source


# GeometryOps.unwrapFunction.
julia
unwrap(target::Type{<:AbstractTrait}, obj)
unwrap(f, target::Type{<:AbstractTrait}, obj)

Unwrap the object to vectors, down to the target trait.

If f is passed in it will be applied to the target geometries as they are found.

source


# GeometryOps.weighted_meanMethod.
julia
weighted_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


# GeometryOps.withinMethod.
julia
within(geom1, geom2)::Bool

Return true if the first geometry is completely within the second geometry. The interiors of both geometries must intersect and the interior and boundary of the primary geometry (geom1) must not intersect the exterior of the secondary geometry (geom2).

Furthermore, within returns the exact opposite result of contains.

Examples

julia
import GeometryOps as GO, GeoInterface as GI

line = GI.LineString([(1, 1), (1, 2), (1, 3), (1, 4)])
point = (1, 2)
GO.within(point, line)

# output
true

source



  1. K. Hormann and N. Sukumar. Generalized Barycentric Coordinates in Computer Graphics and Computational Mechanics. Taylor & Fancis, CRC Press, 2017. ↩︎