Convex Hull

LazySets.API.convex_hullMethod
convex_hull(X::LazySet, Y::LazySet; [algorithm]=nothing,
            [backend]=nothing, [solver]=nothing)

Compute the convex hull of two polytopic sets.

Input

  • X – polytopic set
  • Y – polytopic set
  • algorithm – (optional, default: nothing) the convex-hull algorithm
  • backend – (optional, default: nothing) backend for polyhedral computations (used for higher-dimensional sets)
  • solver – (optional, default: nothing) the linear-programming solver used in the backend

Output

If the sets are empty, the result is an EmptySet. If the convex hull consists of a single point, the result is a Singleton. If the input sets are one-dimensional, the result is an Interval. If the input sets are two-dimensional, the result is a VPolygon. Otherwise the result is a VPolytope.

Algorithm

We compute the vertices of both X and Y using vertices_list and then compute the convex hull of the union of those vertices.

source
LazySets.API.convex_hullMethod
convex_hull(P1::HPoly, P2::HPoly;
           [backend]=default_polyhedra_backend(P1))

Compute the convex hull of the set union of two polyhedra in constraint representation.

Input

  • P1 – polyhedron
  • P2 – polyhedron
  • backend – (optional, default: default_polyhedra_backend(P1)) the backend for polyhedral computations

Output

The HPolyhedron (resp. HPolytope) obtained by the concrete convex hull of P1 and P2.

Notes

For performance reasons, it is suggested to use the CDDLib.Library() backend for the convex_hull.

For further information on the supported backends see Polyhedra's documentation.

source
LazySets.API.convex_hullMethod
convex_hull(points::Vector{VN}; [algorithm]=nothing, [backend]=nothing,
            [solver]=nothing) where {N, VN<:AbstractVector{N}}

Compute the convex hull of a list of points.

Input

  • points – list of points
  • algorithm – (optional, default: nothing) the convex-hull algorithm; see below for valid options
  • backend – (optional, default: nothing) polyhedral computation backend for higher-dimensional point sets
  • solver – (optional, default: nothing) the linear-programming solver used in the backend

Output

The convex hull as a list of points.

Algorithm

A pre-processing step treats the cases with up to two points for one dimension and up to four points for two dimensions. For more points in one resp. two dimensions, we use more general algorithms.

For the one-dimensional case, we return the minimum and maximum points, in that order.

The two-dimensional case is handled with a planar convex-hull algorithm. The following algorithms are available:

  • "monotone_chain" – compute the convex hull of points in the plane using Andrew's monotone-chain method
  • "monotone_chain_sorted" – the same as "monotone_chain" but assuming that the points are already sorted in counter-clockwise fashion

See the reference docstring of each of those algorithms for details.

The higher-dimensional case is treated using the concrete polyhedra library Polyhedra, which gives access to libraries such as CDDLib and ConvexHull.jl. These libraries can be chosen via the backend argument.

Notes

For the in-place version use convex_hull! instead of convex_hull.

Examples

Compute the convex hull of a random set of points:

julia> points = [randn(2) for i in 1:30]; # 30 random points in 2D

julia> hull = convex_hull(points);

julia> typeof(hull)
Vector{Vector{Float64}} (alias for Array{Array{Float64, 1}, 1})
source
LazySets.monotone_chain!Function
monotone_chain!(points::Vector{VN}; sort::Bool=true
               ) where {N, VN<:AbstractVector{N}}

Compute the convex hull of a list of points in the plane using Andrew's monotone-chain method.

Input

  • points – list of 2D vectors; will be sorted in-place inside this method
  • sort – (optional, default: true) flag for sorting the vertices lexicographically; sortedness is required for correctness

Output

List of vectors containing the 2D coordinates of the corner points of the convex hull.

Notes

For large sets of points, it is convenient to use static vectors to get maximum performance. For information on how to convert usual vectors into static vectors, see the type SVector provided by the StaticArrays package.

Algorithm

This method implements Andrew's monotone-chain convex hull algorithm to construct the convex hull of a set of $n$ points in the plane in $O(n \log n)$ time. For further details see Monotone chain

source