Other
DelaunayTriangulation.distance_to_polygon
— Functiondistance_to_polygon(q, points, boundary_nodes) -> Number
Given a query point q
and a polygon defined by (points, boundary_nodes)
, returns the signed distance from q
to the polygon. The boundary_nodes
must match the specification in the documentation and in check_args
.
See also dist
.
DelaunayTriangulation.number_type
— Functionnumber_type(x) -> DataType
Given a container x
, returns the number type used for storing coordinates.
Examples
julia> using DelaunayTriangulation
julia> DelaunayTriangulation.number_type([1, 2, 3])
Int64
julia> DelaunayTriangulation.number_type((1, 2, 3))
Int64
julia> DelaunayTriangulation.number_type([1.0 2.0 3.0; 4.0 5.0 6.0])
Float64
julia> DelaunayTriangulation.number_type([[[1, 2, 3, 4, 5, 1]], [[6, 8, 9], [9, 10, 11], [11, 12, 6]]])
Int64
julia> DelaunayTriangulation.number_type((1.0f0, 2.0f0))
Float32
julia> DelaunayTriangulation.number_type(Vector{Float64})
Float64
julia> DelaunayTriangulation.number_type(Vector{Vector{Float64}})
Float64
julia> DelaunayTriangulation.number_type(NTuple{2, Float64})
Float64
DelaunayTriangulation.pole_of_inaccessibility
— Functionpole_of_inaccessibility(points, boundary_nodes; precision = one(number_type(points)))
Finds the pole of inaccessibility for the polygon defined by points
and boundary_nodes
. The boundary_nodes
must match the specification given in the documentation. See also check_args
for this specification.
Arguments
points
: The points defining the polygon.boundary_nodes
: The boundary nodes defining the polygon.
Keyword Arguments
precision=one(number_type(points))
: The precision of the returned pole. The default isone(number_type(points))
.
Outputs
x
: The x-coordinate of the pole.y
: The y-coordinate of the pole.
Extended help
The pole of inaccessibility is the point within a polygon that is furthest from an edge. For DelaunayTriangulation.jl, this is useful as it is a representative point for ghost edges that is guaranteed to be inside the polygon, in contrast to for example a centroid which is not always inside the polygon. Some useful links are this blog post and the the original repo. Our implementation is partially based on on the python implementation and this other Julia implementation.
DelaunayTriangulation.clip_polygon
— Functionclip_polygon(vertices, points, clip_vertices, clip_points; predicates::AbstractPredicateKernel=AdaptiveKernel()) -> Vector
Clip a polygon defined by (vertices, points)
to a convex clip polygon defined by (clip_vertices, clip_points)
with the Sutherland-Hodgman algorithm. The polygons should be defined in counter-clockwise order.
Arguments
vertices
: The vertices of the polygon to be clipped.points
: The underlying point set that the vertices are defined over.clip_vertices
: The vertices of the clipping polygon.clip_points
: The underlying point set that the clipping vertices are defined over.
Keyword Arguments
predicates::AbstractPredicateKernel=AdaptiveKernel()
: Method to use for computing predicates. Can be one ofFastKernel
,ExactKernel
, andAdaptiveKernel
. See the documentation for a further discussion of these methods.
Output
clipped_polygon
: The coordinates of the clipped polygon, given in counter-clockwise order andclipped_polygon[begin] == clipped_polygon[end]
.
DelaunayTriangulation.construct_polygon_hierarchy
— Functionconstruct_polygon_hierarchy(points; IntegerType=Int) -> PolygonHierarchy{IntegerType}
Returns a PolygonHierarchy
defining the polygon hierarchy for a given set of points
. This defines a hierarchy with a single polygon.
construct_polygon_hierarchy(points, boundary_nodes; IntegerType=Int) -> PolygonHierarchy{IntegerType}
Returns a PolygonHierarchy
defining the polygon hierarchy for a given set of boundary_nodes
that define a set of piecewise linear curves.
construct_polygon_hierarchy(points, boundary_nodes, boundary_curves; IntegerType=Int, n=4096) -> PolygonHierarchy{IntegerType}
Returns a PolygonHierarchy
defining the polygon hierarchy for a given set of boundary_nodes
that define a curve-bounded domain from the curves in boundary_curves
. Uses polygonise
to fill in the boundary curves.
Arguments
points
: The point set.boundary_nodes
: The boundary nodes. These should be the output fromconvert_boundary_curves!
.boundary_curves
: The boundary curves. These should be the output fromconvert_boundary_curves!
.
Keyword Arguments
IntegerType=Int
: The integer type to use for indexing the polygons.n=4096
: The number of points to use for filling in the boundary curves inpolygonise
.