General sets (LazySet)
Every set type in this library is a subtype of the abstract type LazySet
.
LazySets.LazySet
— TypeLazySet{N}
Abstract type for the set types in LazySets.
Notes
LazySet
types should be parameterized with a type N
, typically N<:Real
, for using different numeric types.
Every concrete LazySet
must define the following method:
dim(S::LazySet)
– the ambient dimension ofS
While not strictly required, it is useful to define the following method:
σ(d::AbstractVector, S::LazySet)
– the support vector ofS
in a given directiond
The method
ρ(d::AbstractVector, S::LazySet)
– the support function ofS
in a given directiond
is optional because there is a fallback implementation relying on σ
. However, for potentially unbounded sets (which includes most lazy set types) this fallback cannot be used and an explicit method must be implemented.
The subtypes of LazySet
(including abstract interfaces):
julia> subtypes(LazySet, false)
17-element Vector{Any}:
AbstractAffineMap
AbstractPolynomialZonotope
Bloating
CachedMinkowskiSumArray
CartesianProduct
CartesianProductArray
Complement
ConvexSet
Intersection
IntersectionArray
MinkowskiSum
MinkowskiSumArray
Polygon
QuadraticMap
Rectification
UnionSet
UnionSetArray
If we only consider concrete subtypes, then:
julia> concrete_subtypes = subtypes(LazySet, true);
julia> length(concrete_subtypes)
54
julia> println.(concrete_subtypes);
AffineMap
Ball1
Ball2
BallInf
Ballp
Bloating
CachedMinkowskiSumArray
CartesianProduct
CartesianProductArray
Complement
ConvexHull
ConvexHullArray
DensePolynomialZonotope
Ellipsoid
EmptySet
ExponentialMap
ExponentialProjectionMap
HParallelotope
HPolygon
HPolygonOpt
HPolyhedron
HPolytope
HalfSpace
Hyperplane
Hyperrectangle
Intersection
IntersectionArray
Interval
InverseLinearMap
Line
Line2D
LineSegment
LinearMap
MinkowskiSum
MinkowskiSumArray
Polygon
QuadraticMap
Rectification
ResetMap
RotatedHyperrectangle
SimpleSparsePolynomialZonotope
Singleton
SparsePolynomialZonotope
Star
SymmetricIntervalHull
Tetrahedron
Translation
UnionSet
UnionSetArray
Universe
VPolygon
VPolytope
ZeroSet
Zonotope
Plotting
Plotting via the Plots
package is available for one- or two-dimensional sets. The default algorithm is to plot an outer approximation using the support function (1D) respectively the support vector (2D). This means that (1) plotting will fail if these functionalities are not available (e.g., for lazy Intersection
s) and (2) that plots of non-convex sets can be misleading. The implementation below internally relies on the function plot_recipe
. For some set types (e.g., Intersection
), the default implementation is overridden.
RecipesBase.apply_recipe
— Methodplot_lazyset(X::LazySet{N}, [ε]::Real=N(PLOT_PRECISION); ...) where {N}
Plot a set.
Input
X
– setε
– (optional, default:PLOT_PRECISION
) approximation error bound
Notes
This recipe just defines the default plotting options and then calls the function plot_recipe
, which then implements the set-specific plotting.
The argument ε
is ignored by some set types, e.g., for polyhedra (subtypes of AbstractPolyhedron
).
Examples
julia> B = Ball2(ones(2), 0.1);
julia> plot(B, 1e-3) # default accuracy value (explicitly given for clarity here)
julia> plot(B, 1e-2) # faster but less accurate than the previous call
RecipesBase.apply_recipe
— Methodplot_list(list::AbstractVector{VN}, [ε]::Real=N(PLOT_PRECISION),
[Nφ]::Int=PLOT_POLAR_DIRECTIONS; [same_recipe]=false; ...)
where {N, VN<:LazySet{N}}
Plot a list of sets.
Input
list
– list of sets (1D or 2D)ε
– (optional, default:PLOT_PRECISION
) approximation error boundNφ
– (optional, default:PLOT_POLAR_DIRECTIONS
) number of polar directions (used to plot lazy intersections)same_recipe
– (optional, default:false
) switch for faster plotting but without individual plot recipes (see notes below)
Notes
For each set in the list we apply an individual plot recipe.
The option same_recipe
provides access to a faster plotting scheme where all sets in the list are first converted to polytopes and then plotted in one single run. This, however, is not suitable when plotting flat sets (line segments, singletons) because then the polytope plot recipe does not deliver good results. Hence by default we do not use this option. For plotting a large number of (non-flat) polytopes, we highly advise activating this option.
Examples
julia> B1 = BallInf(zeros(2), 0.4);
julia> B2 = BallInf(ones(2), 0.4);
julia> plot([B1, B2])
Some of the sets in the list may not be plotted precisely but rather overapproximated first. The second argument ε
controls the accuracy of this overapproximation.
julia> Bs = [BallInf(zeros(2), 0.4), Ball2(ones(2), 0.4)];
julia> plot(Bs, 1e-3) # default accuracy value (explicitly given for clarity)
julia> plot(Bs, 1e-2) # faster but less accurate than the previous call
LazySets.plot_vlist
— Methodplot_vlist(X::S, ε::Real) where {S<:LazySet}
Return a list of vertices used for plotting a two-dimensional set.
Input
X
– two-dimensional setε
– precision parameter
Output
A list of vertices of a polygon P
. For convex X
, P
usually satisfies that the Hausdorff distance to X
is less than ε
.
For three-dimensional sets, we support Makie
:
LazySets.plot3d
— Functionplot3d(S::LazySet; [backend]=default_polyhedra_backend(S), [alpha]=1.0,
[color]=:blue, [colormap]=:viridis, [colorrange]=nothing,
[interpolate]=false, [linewidth]=1, [overdraw]=false, [shading]=true,
[transparency]=true, [visible]=true)
Plot a three-dimensional set using Makie
.
Input
S
– setbackend
– (optional, default:default_polyhedra_backend(S)
) backend for polyhedral computationsalpha
– (optional, default:1.0
) float in[0,1]
; the alpha or transparency valuecolor
– (optional, default::blue
)Symbol
orColorant
; the color of the main plot element (markers, lines, etc.), which can be a color symbol/string like:red
colormap
– (optional, default::viridis
) the color map of the main plot; useavailable_gradients()
to see which gradients are available, which can also be used as[:red, :black]
colorrange
– (optional, default:nothing
, which falls back toMakie.Automatic()
) a tuple(min, max)
wheremin
andmax
specify the data range to be used for indexing thecolormap
interpolate
– (optional, default:false
) a boolean for heatmap and images; toggles color interpolation between nearby pixelslinewidth
– (optional, default:1
) a number that specifies the width of the line inline
andlinesegments
plotsoverdraw
– (optional, default:false
)shading
– (optional, default:true
) a boolean that toggles shading (for meshes)transparency
– (optional, default:true
) iftrue
, the set is transparent, otherwise it is displayed as a solid objectvisible
– (optional, default:true
) a boolean that toggles visibility of the plot
For a complete list of attributes and usage see Makie's documentation.
Notes
This plot recipe works by computing the list of constraints of S
and converting to a polytope in H-representation. Then, this polytope is transformed with Polyhedra.Mesh
and plotted using the mesh
function.
If the function constraints_list
is not applicable to your set S
, try overapproximation first; e.g. via
julia> Sapprox = overapproximate(S, SphericalDirections(10))
julia> using Polyhedra, GLMakie
julia> plot3d(Sapprox)
The number 10
above corresponds to the number of directions considered; for better resolution use higher values (but it will take longer).
For efficiency consider using the CDDLib
backend, as in
julia> using CDDLib
julia> plot3d(Sapprox, backend=CDDLib.Library())
Examples
The functionality requires both Polyhedra
and a Makie
backend. After loading LazySets
, do using Polyhedra, GLMakie
(or another Makie backend).
julia> using LazySets, Polyhedra, GLMakie
julia> plot3d(10 * rand(Hyperrectangle, dim=3))
julia> plot3d!(10 * rand(Hyperrectangle, dim=3), color=:red)
LazySets.plot3d!
— Functionplot3d!(S::LazySet; backend=default_polyhedra_backend(S), [alpha]=1.0,
[color]=:blue, [colormap]=:viridis, [colorrange]=nothing,
[interpolate]=false, [linewidth]=1, [overdraw]=false, [shading]=true,
[transparency]=true, [visible]=true)
Plot a three-dimensional set using Makie.
Input
See plot3d
for the description of the inputs. For a complete list of attributes and usage see Makie's documentation.
Notes
See the documentation of plot3d
for examples.
Globally defined set functions
LazySets.○
— Method○(c, a)
Convenience constructor of Ellipsoid
s or Ball2
s depending on the type of a
.
Input
c
– centera
– additional parameter (either a shape matrix forEllipsoid
or a radius forBall2
)
Output
A Ellipsoid
s or Ball2
s depending on the type of a
.
Notes
The function symbol can be typed via \bigcirc[TAB]
.
LazySets.isconvextype
— Methodisconvextype(X::Type{<:LazySet})
Check whether the given LazySet
type is convex.
Input
X
– subtype ofLazySet
Output
true
if the given set type is guaranteed to be convex by using only type information, and false
otherwise.
Notes
Since this operation only acts on types (not on values), it can return false negatives, i.e., there may be instances where the set is convex, even though the answer of this function is false
. The examples below illustrate this point.
Examples
A ball in the infinity norm is always convex, hence we get:
julia> isconvextype(BallInf)
true
For instance, the union (UnionSet
) of two sets may in general be either convex or not. Since convexity cannot be decided by just using type information, isconvextype
returns false
.
julia> isconvextype(UnionSet)
false
However, the type parameters of set operations allow to decide convexity in some cases by falling back to the convexity information of the type of its arguments. Consider for instance the lazy intersection. The intersection of two convex sets is always convex, hence we get:
julia> isconvextype(Intersection{Float64, BallInf{Float64}, BallInf{Float64}})
true
LazySets.low
— Methodlow(X::LazySet)
Return a vector with the lowest coordinates of the set in each canonical direction.
Input
X
– set
Output
A vector with the lower coordinate of the set in each dimension.
Notes
See also low(X::LazySet, i::Int)
.
The result is the lowermost corner of the box approximation, so it is not necessarily contained in X
.
LazySets.high
— Methodhigh(X::LazySet)
Return a vector with the highest coordinate of the set in each canonical direction.
Input
X
– set
Output
A vector with the highest coordinate of the set in each dimension.
Notes
See also high(X::LazySet, i::Int)
.
The result is the uppermost corner of the box approximation, so it is not necessarily contained in X
.
Base.extrema
— Methodextrema(X::LazySet, i::Int)
Return the lower and higher coordinate of a set in a given dimension.
Input
X
– seti
– dimension of interest
Output
The lower and higher coordinate of the set in the given dimension.
Notes
The result is equivalent to (low(X, i), high(X, i))
, but sometimes it can be computed more efficiently.
Algorithm
The bounds are computed with low
and high
.
Base.extrema
— Methodextrema(X::LazySet)
Return two vectors with the lowest and highest coordinate of X
in each canonical direction.
Input
X
– set
Output
Two vectors with the lowest and highest coordinates of X
in each dimension.
Notes
See also extrema(X::LazySet, i::Int)
.
The result is equivalent to (low(X), high(X))
, but sometimes it can be computed more efficiently.
The resulting points are the lowermost and uppermost corners of the box approximation, so they are not necessarily contained in X
.
Algorithm
The bounds are computed with low
and high
by default.
LazySets.convex_hull
— Methodconvex_hull(X::LazySet; kwargs...)
Compute the convex hull of a polytopic set.
Input
X
– polytopic set
Output
The set X
itself if its type indicates that it is convex, or a new set with the list of the vertices describing the convex hull.
Algorithm
For non-convex sets, this method relies on the vertices_list
method.
LazySets.triangulate
— Methodtriangulate(X::LazySet)
Triangulate a three-dimensional polyhedral set.
Input
X
– three-dimensional polyhedral set
Output
A tuple (p, c)
where p
is a matrix, with each column containing a point, and c
is a list of 3-tuples containing the indices of the points in each triangle.
LazySets.basetype
— Functionbasetype(T::Type{<:LazySet})
Return the base type of the given set type (i.e., without type parameters).
Input
T
– set type
Output
The base type of T
.
basetype(S::LazySet)
Return the base type of the given set (i.e., without type parameters).
Input
S
– set
Output
The base type of S
.
Examples
julia> Z = rand(Zonotope);
julia> basetype(Z)
Zonotope
julia> basetype(Z + Z)
MinkowskiSum
julia> basetype(LinearMap(rand(2, 2), Z + Z))
LinearMap
LazySets.isboundedtype
— Methodisboundedtype(T::Type{<:LazySet})
Check whether a set type only represents bounded sets.
Input
T
– set type
Output
true
if the set type only represents bounded sets. Note that some sets may still represent an unbounded set even though their type actually does not (example: HPolytope
, because the construction with non-bounding linear constraints is allowed).
Notes
By default this function returns false
. All set types that can determine boundedness should override this behavior.
LazySets.isbounded
— Methodisbounded(S::LazySet)
Check whether a set is bounded.
Input
S
– setalgorithm
– (optional, default:"support_function"
) algorithm choice, possible options are"support_function"
and"stiemke"
Output
true
iff the set is bounded.
Algorithm
See the documentation of _isbounded_unit_dimensions
or _isbounded_stiemke
for details.
LazySets._isbounded_unit_dimensions
— Method_isbounded_unit_dimensions(S::LazySet)
Check whether a set is bounded in each unit dimension.
Input
S
– set
Output
true
iff the set is bounded in each unit dimension.
Algorithm
This function asks for upper and lower bounds in each ambient dimension.
LazySets.is_polyhedral
— Methodis_polyhedral(S::LazySet)
Trait for polyhedral sets.
Input
S
– set
Output
true
only if the set behaves like an AbstractPolyhedron
.
Notes
The answer is conservative, i.e., may sometimes be false
even if the set is polyhedral.
LazySets.isfeasible
— Functionisfeasible(constraints::AbstractVector{<:HalfSpace}, [witness]::Bool=false;
[solver]=nothing)
Check for feasibility of a list of linear constraints.
Input
constraints
– list of linear constraintswitness
– (optional; default:false
) flag for witness productionsolver
– (optional; default:nothing
) LP solver
Output
true
if the linear constraints are feasible, and false
otherwise.
Algorithm
This implementation solves the corresponding feasibility linear program.
LinearAlgebra.norm
— Functionnorm(S::LazySet, [p]::Real=Inf)
Return the norm of a set. It is the norm of the enclosing ball (of the given $p$-norm) of minimal volume that is centered in the origin.
Input
S
– setp
– (optional, default:Inf
) norm
Output
A real number representing the norm.
IntervalArithmetic.radius
— Functionradius(S::LazySet, [p]::Real=Inf)
Return the radius of a set. It is the radius of the enclosing ball (of the given $p$-norm) of minimal volume with the same center.
Input
S
– setp
– (optional, default:Inf
) norm
Output
A real number representing the radius.
LazySets.diameter
— Functiondiameter(S::LazySet, [p]::Real=Inf)
Return the diameter of a set. It is the maximum distance between any two elements of the set, or, equivalently, the diameter of the enclosing ball (of the given $p$-norm) of minimal volume with the same center.
Input
S
– setp
– (optional, default:Inf
) norm
Output
A real number representing the diameter.
Base.isempty
— Methodisempty(P::LazySet{N}, witness::Bool=false;
[use_polyhedra_interface]::Bool=false, [solver]=nothing,
[backend]=nothing) where {N}
Check whether a polyhedral set is empty.
Input
P
– polyhedral setwitness
– (optional, default:false
) compute a witness if activateduse_polyhedra_interface
– (optional, default:false
) iftrue
, we use thePolyhedra
interface for the emptiness testsolver
– (optional, default:nothing
) LP-solver backend; usesdefault_lp_solver(N)
if not providedbackend
– (optional, default:nothing
) backend for polyhedral computations inPolyhedra
; usesdefault_polyhedra_backend(P)
if not provided
Output
- If
witness
option is deactivated:true
iff $P = ∅$ - If
witness
option is activated:(true, [])
iff $P = ∅$(false, v)
iff $P ≠ ∅$ and $v ∈ P$
Notes
The default value of the backend
is set internally and depends on whether the use_polyhedra_interface
option is set or not. If the option is set, we use default_polyhedra_backend(P)
.
Witness production is not supported if use_polyhedra_interface
is true
.
Algorithm
The algorithm sets up a feasibility LP for the constraints of P
. If use_polyhedra_interface
is true
, we call Polyhedra.isempty
. Otherwise, we set up the LP internally.
LazySets.linear_map
— Methodlinear_map(M::AbstractMatrix, P::LazySet; kwargs...)
Concrete linear map of a polyhedral set.
Input
M
– matrixP
– polyhedral set
Output
A set representing the concrete linear map.
LazySets.linear_map
— Methodlinear_map(a::Number, X::LazySet; kwargs...)
Alias for scale(a, X; kwargs...)
.
LazySets.affine_map
— Methodaffine_map(M, X::LazySet, v::AbstractVector; kwargs...)
Compute the concrete affine map $M·X + v$.
Input
M
– linear mapX
– setv
– translation vector
Output
A set representing the affine map $M·X + v$.
Algorithm
The implementation applies the functions linear_map
and translate
.
LazySets.exponential_map
— Methodexponential_map(M::AbstractMatrix, X::LazySet)
Compute the concrete exponential map of M
and X
, i.e., exp(M) * X
.
Input
M
– matrixX
– set
Output
A set representing the exponential map of M
and X
.
Algorithm
The implementation applies the functions exp
and linear_map
.
LazySets.an_element
— Methodan_element(S::LazySet)
Return some element of a set.
Input
S
– set
Output
An element of a set.
Algorithm
An element of the set is obtained by evaluating its support vector along direction $[1, 0, …, 0]$. This may fail for unbounded sets.
LazySets.tosimplehrep
— Methodtosimplehrep(S::LazySet)
Return the simple constraint representation $Ax ≤ b$ of a polyhedral set from its list of linear constraints.
Input
S
– polyhedral set
Output
The tuple (A, b)
where A
is the matrix of normal directions and b
is the vector of offsets.
Algorithm
This fallback implementation relies on constraints_list(S)
.
LazySets.reflect
— Methodreflect(P::LazySet)
Concrete reflection of a set P
, resulting in the reflected set -P
.
Input
P
– polyhedral set
Output
The set -P
, which is either of type HPolytope
if P
is a polytope (i.e., bounded) or of type HPolyhedron
otherwise.
Algorithm
This function requires that the list of constraints of the set P
is available, i.e., that it can be written as $P = \{z ∈ ℝⁿ: ⋂ sᵢᵀz ≤ rᵢ, i = 1, ..., N\}.$
This function can be used to implement the alternative definition of the Minkowski Difference
\[ A ⊖ B = \{a − b \mid a ∈ A, b ∈ B\} = A ⊕ (-B)\]
by calling minkowski_sum(A, reflect(B))
.
LazySets.is_interior_point
— Methodis_interior_point(d::AbstractVector{N}, X::LazySet{N};
p=N(Inf), ε=_rtol(N)) where {N}
Check whether the point d
is contained in the interior of the set X
.
Input
d
– pointX
– setp
– (optional; default:N(Inf)
) norm of the ball used to apply the error toleranceε
– (optional; default:_rtol(N)
) error tolerance of check
Output
Boolean which indicates if the point d
is contained in X
.
Algorithm
The implementation checks if a Ballp
of norm p
with center d
and radius ε
is contained in the set X
. This is a numerical check for d ∈ interior(X)
with error tolerance ε
.
LazySets.isoperationtype
— Methodisoperationtype(X::Type{<:LazySet})
Check whether the given set type is an operation or not.
Input
X
– set type
Output
true
if the given set type is a set operation and false
otherwise.
Notes
This fallback implementation returns an error that isoperationtype
is not implemented. Subtypes of LazySet
should dispatch on this function as required.
See also isoperation(X<:LazySet)
.
Examples
julia> isoperationtype(BallInf)
false
julia> isoperationtype(LinearMap)
true
LazySets.isoperation
— Methodisoperation(X::LazySet)
Check whether a set is an instance of a set operation or not.
Input
X
– set
Output
true
if X
is an instance of a set-based operation and false
otherwise.
Notes
This fallback implementation checks whether the set type of the input is an operation type using isoperationtype(::Type{<:LazySet})
.
Examples
julia> B = BallInf([0.0, 0.0], 1.0);
julia> isoperation(B)
false
julia> isoperation(B ⊕ B)
true
LazySets.isequivalent
— Methodisequivalent(X::LazySet, Y::LazySet)
Check whether two sets are equal in the mathematical sense, i.e., equivalent.
Input
X
– setY
– set
Output
true
iff X
is equivalent to Y
(up to some precision).
Algorithm
First we check X ≈ Y
, which returns true
if and only if X
and Y
have the same type and approximately the same values (checked with LazySets._isapprox
). If that fails, we check the double inclusion X ⊆ Y && Y ⊆ X
.
Examples
julia> X = BallInf([0.1, 0.2], 0.3);
julia> Y = convert(HPolytope, X);
julia> X == Y
false
julia> isequivalent(X, Y)
true
LazySets.surface
— Methodsurface(X::LazySet)
Compute the surface area of a set.
Input
X
– set
Output
A real number representing the surface area of X
.
LazySets.area
— Methodarea(X::LazySet)
Compute the area of a two-dimensional polytopic set.
Input
X
– two-dimensional polytopic set
Output
A number representing the area of X
.
Notes
This algorithm is applicable to any polytopic set X
whose list of vertices can be computed via vertices_list
.
Algorithm
Let m
be the number of vertices of X
. We consider the following instances:
m = 0, 1, 2
: the output is zero.m = 3
: the triangle case is solved using the Shoelace formula with 3 points.m = 4
: the quadrilateral case is solved by the factored version of the Shoelace formula with 4 points.
Otherwise, the general Shoelace formula is used; for details see the Wikipedia page.
LazySets.concretize
— Methodconcretize(X::LazySet)
Construct a concrete representation of a (possibly lazy) set.
Input
X
– set
Output
A concrete representation of X
(as far as possible).
Notes
Since not every lazy set has a concrete set representation in this library, the result may be partially lazy.
LazySets.complement
— Methodcomplement(X::LazySet)
Return the complement of a polyhedral set.
Input
X
– polyhedral set
Output
A UnionSetArray
of half-spaces, i.e., the output is the union of the linear constraints which are obtained by complementing each constraint of X
.
Algorithm
The principle used in this implementation is that for any pair of sets $(X, Y)$ we have that $(X ∩ Y)^C = X^C ∪ Y^C$. In particular, we can apply this rule for each constraint that defines a polyhedral set. Hence the concrete complement can be represented as the set union of the complement of each constraint.
Polyhedra.polyhedron
— Methodpolyhedron(P::LazySet; [backend]=default_polyhedra_backend(P))
Compute a set representation from Polyhedra.jl
.
Input
P
– polyhedral setbackend
– (optional, default: calldefault_polyhedra_backend(P)
) the polyhedral computations backend
Output
A set representation in the Polyhedra
library.
Notes
For further information on the supported backends see Polyhedra's documentation.
Algorithm
This default implementation uses tosimplehrep
, which computes the constraint representation of P
. Set types preferring the vertex representation should implement their own method.
LazySets.project
— Functionproject(S::LazySet, block::AbstractVector{Int}, [::Nothing=nothing],
[n]::Int=dim(S); [kwargs...])
Project a set to a given block by using a concrete linear map.
Input
S
– setblock
– block structure - a vector with the dimensions of interestnothing
– (default:nothing
)n
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A set representing the projection of the set S
to block block
.
Algorithm
We apply the function linear_map
.
LazySets.project
— Methodproject(S::LazySet, block::AbstractVector{Int}, set_type::Type{TS},
[n]::Int=dim(S); [kwargs...]) where {TS<:LazySet}
Project a set to a given block and set type, possibly involving an overapproximation.
Input
S
– setblock
– block structure - a vector with the dimensions of interestset_type
– target set typen
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A set of type set_type
representing an overapproximation of the projection of S
.
Algorithm
- Project the set
S
withM⋅S
, whereM
is the identity matrix in the block
coordinates and zero otherwise.
- Overapproximate the projected set using
overapproximate
andset_type
.
LazySets.project
— Methodproject(S::LazySet, block::AbstractVector{Int},
set_type_and_precision::Pair{T, N}, [n]::Int=dim(S);
[kwargs...]) where {T<:UnionAll, N<:Real}
Project a set to a given block and set type with a certified error bound.
Input
S
– setblock
– block structure - a vector with the dimensions of interestset_type_and_precision
– pair(T, ε)
of a target set typeT
and an error boundε
for approximationn
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A set representing the epsilon-close approximation of the projection of S
.
Notes
Currently we only support HPolygon
as set type, which implies that the set must be two-dimensional.
Algorithm
- Project the set
S
withM⋅S
, whereM
is the identity matrix in the block
coordinates and zero otherwise.
- Overapproximate the projected set with the given error bound
ε
.
LazySets.project
— Functionproject(S::LazySet, block::AbstractVector{Int}, ε::Real, [n]::Int=dim(S);
[kwargs...])
Project a set to a given block and set type with a certified error bound.
Input
S
– setblock
– block structure - a vector with the dimensions of interestε
– error bound for approximationn
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A set representing the epsilon-close approximation of the projection of S
.
Algorithm
- Project the set
S
withM⋅S
, whereM
is the identity matrix in the block
coordinates and zero otherwise.
- Overapproximate the projected set with the given error bound
ε
.
The target set type is chosen automatically.
ReachabilityBase.Arrays.rectify
— Functionrectify(X::LazySet, [concrete_intersection]::Bool=false)
Concrete rectification of a set.
Input
X
– setconcrete_intersection
– (optional, default:false
) flag to compute concrete intersections for intermediate results
Output
A set corresponding to the rectification of X
, which is in general a union of linear maps of intersections.
Algorithm
For each dimension in which X
is both positive and negative, we split X
into these two parts. Additionally we project the negative part to zero.
SparseArrays.permute
— Functionpermute(X::LazySet, p::AbstractVector{Int})
Permute the dimensions of a set according to a given permutation vector.
Input
X
– setp
– permutation vector
Output
A new set corresponding to X
where the dimensions have been permuted according to p
.
Base.rationalize
— Methodrationalize(::Type{T}, X::LazySet{<:AbstractFloat}, tol::Real)
where {T<:Integer}
Approximate a set of floating-point numbers as a set whose entries are rationals of the given integer type.
Input
T
– (optional, default:Int
) integer type to represent the rationalsX
– set which has floating-point componentstol
– (optional, default:eps(N)
) tolerance of the result; each rationalized component will differ by no more thantol
with respect to the floating-point value
Output
A set of the same base type of X
where each numerical component is of type Rational{T}
.
LazySets.singleton_list
— Methodsingleton_list(P::LazySet)
Return the vertices of a polytopic set as a list of singletons.
Input
P
– polytopic set
Output
A list of the vertices of P
as Singleton
s.
Notes
This function relies on vertices_list
, which raises an error if the set is not polytopic (e.g., unbounded).
LazySets.constraints
— Methodconstraints(X::LazySet)
Construct an iterator over the constraints of a polyhedral set.
Input
X
– polyhedral set
Output
An iterator over the constraints of X
.
LazySets.vertices
— Methodvertices(X::LazySet)
Construct an iterator over the vertices of a polytopic set.
Input
X
– polytopic set
Output
An iterator over the vertices of X
.
MiniQhull.delaunay
— Functiondelaunay(X::LazySet)
Compute the Delaunay triangulation of the given polytopic set.
Input
X
– polytopic setcompute_triangles_3d
– (optional; default:false
) flag to compute the 2D triangulation of a 3D set
Output
A union of polytopes in vertex representation.
Notes
This implementation requires the package MiniQhull.jl, which uses the library Qhull.
The method works in arbitrary dimension and the requirement is that the list of vertices of X
can be obtained.
LazySets.chebyshev_center_radius
— Methodchebyshev_center_radius(P::LazySet{N};
[backend]=default_polyhedra_backend(P),
[solver]=default_lp_solver_polyhedra(N; presolve=true)
) where {N}
Compute a Chebyshev center and the corresponding radius of a polytopic set.
Input
P
– polytopic setbackend
– (optional; default:default_polyhedra_backend(P)
) the backend for polyhedral computationssolver
– (optional; default:default_lp_solver_polyhedra(N; presolve=true)
) the LP solver passed toPolyhedra
Output
The pair (c, r)
where c
is a Chebyshev center of P
and r
is the radius of the largest ball with center c
enclosed by P
.
Notes
The Chebyshev center is the center of a largest Euclidean ball enclosed by P
. In general, the center of such a ball is not unique, but the radius is.
Algorithm
We call Polyhedra.chebyshevcenter
.
LazySets.scale
— Methodscale(α::Real, X::LazySet)
Concrete scaling of a set.
Input
α
– scalarX
– set
Output
A new set of the same type (if possible).
Notes
This fallback method calls scale!
on a copy of X
.
LazySets.plot_recipe
— Methodplot_recipe(X::LazySet, [ε])
Convert a compact convex set to a pair (x, y)
of points for plotting.
Input
X
– compact convex setε
– approximation-error bound
Output
A pair (x, y)
of points that can be plotted.
Notes
We do not support three-dimensional or higher-dimensional sets at the moment.
Algorithm
One-dimensional sets are converted to an Interval
.
For two-dimensional sets, we first compute a polygonal overapproximation. The second argument, ε
, corresponds to the error in Hausdorff distance between the overapproximating set and X
. On the other hand, if you only want to produce a fast box-overapproximation of X
, pass ε=Inf
.
Finally, we use the plot recipe for the constructed set (interval or polygon).
The following methods are also defined for LazySet
but cannot be documented due to a bug in the documentation package.
LazySets.low
— Methodlow(X::ConvexSet{N}, i::Int) where {N}
Return the lower coordinate of a convex set in a given dimension.
Input
X
– convex seti
– dimension of interest
Output
The lower coordinate of the set in the given dimension.
LazySets.high
— Methodhigh(X::ConvexSet{N}, i::Int) where {N}
Return the higher coordinate of a convex set in a given dimension.
Input
X
– convex seti
– dimension of interest
Output
The higher coordinate of the set in the given dimension.
LazySets.an_element
— Methodan_element(P::AbstractPolyhedron{N};
[solver]=default_lp_solver(N)) where {N}
Return some element of a polyhedron.
Input
P
– polyhedronsolver
– (optional, default:default_lp_solver(N)
) LP solver
Output
An element of the polyhedron, or an error if the polyhedron is empty.
Algorithm
An element is obtained by solving a feasibility linear program.
an_element(U::Universe{N}) where {N}
Return some element of a universe.
Input
U
– universe
Output
The origin.
Support function and support vector
Every LazySet
type must define a function σ
to compute the support vector. The support function, ρ
, can optionally be defined; otherwise, a fallback definition based on σ
is used.
LazySets.σ
— Functionσ
Function to compute the support vector σ.
LazySets.support_vector
— Functionsupport_vector
Alias for the support vector σ.
LazySets.ρ
— Methodρ(d::AbstractVector, S::LazySet)
Evaluate the support function of a set in a given direction.
Input
d
– directionS
– set
Output
The evaluation of the support function of the set S
for the direction d
.
LazySets.support_function
— Functionsupport_function
Alias for the support function ρ.
Set functions that override Base functions
Base.:==
— Method==(X::LazySet, Y::LazySet)
Check whether two sets use exactly the same set representation.
Input
X
– setY
– set
Output
true
iffX
is equal toY
.
Notes
The check is purely syntactic and the sets need to have the same base type. For instance, X::VPolytope == Y::HPolytope
returns false
even if X
and Y
represent the same polytope. However X::HPolytope{Int64} == Y::HPolytope{Float64}
is a valid comparison.
Algorithm
We recursively compare the fields of X
and Y
until a mismatch is found.
Examples
julia> HalfSpace([1], 1) == HalfSpace([1], 1)
true
julia> HalfSpace([1], 1) == HalfSpace([1.0], 1.0)
true
julia> Ball1([0.0], 1.0) == Ball2([0.0], 1.0)
false
Base.:≈
— Method≈(X::LazySet, Y::LazySet)
Check whether two sets of the same type are approximately equal.
Input
X
– setY
– set of the same base type asX
Output
true
iffX
is equal toY
.
Notes
The check is purely syntactic and the sets need to have the same base type. For instance, X::VPolytope ≈ Y::HPolytope
returns false
even if X
and Y
represent the same polytope. However X::HPolytope{Int64} ≈ Y::HPolytope{Float64}
is a valid comparison.
Algorithm
We recursively compare the fields of X
and Y
until a mismatch is found.
Examples
julia> HalfSpace([1], 1) ≈ HalfSpace([1], 1)
true
julia> HalfSpace([1], 1) ≈ HalfSpace([1.00000001], 0.99999999)
true
julia> Ball1([0.0], 1.0) ≈ Ball2([0.0], 1.0)
false
Base.copy
— Methodcopy(S::LazySet)
Return a copy of a set by copying its values recursively.
Input
S
– set
Output
A copy of S
.
Notes
This function computes a copy
of each field in S
. See the documentation of ?copy
for further details.
Base.eltype
— Functioneltype(::Type{<:LazySet{N}}) where {N}
Return the numeric type (N
) of the given set type.
Input
T
– set type
Output
The numeric type of T
.
eltype(::LazySet{N}) where {N}
Return the numeric type (N
) of the given set.
Input
X
– set
Output
The numeric type of X
.
Aliases for set types
LazySets.CompactSet
— TypeCompactSet
An alias for compact set types.
Notes
Most lazy operations are not captured by this alias because whether their result is compact or not depends on the argument(s).
LazySets.NonCompactSet
— TypeNonCompactSet
An alias for non-compact set types.
Notes
Most lazy operations are not captured by this alias because whether their result is non-compact or not depends on the argument(s).
Implementations
Concrete set representations:
Lazy set operations:
- Affine map (AffineMap)
- Linear map (LinearMap)
- Exponential map (ExponentialMap)
- Exponential projection map (ExponentialProjectionMap)
- Reset map (ResetMap)
- Translation
- Bloating
- Binary Cartesian product (CartesianProduct)
- $n$-ary Cartesian product (CartesianProductArray)
- Binary convex hull (ConvexHull)
- $n$-ary convex hull (ConvexHullArray)
- Binary intersection
- $n$-ary intersection (IntersectionArray)
- Binary Minkowski sum (MinkowskiSum)
- $n$-ary Minkowski sum (MinkowskiSumArray)
- $n$-ary Minkowski sum with cache (CachedMinkowskiSumArray)
- Binary set union (UnionSet)
- $n$-ary set union (UnionSetArray)
- Complement
- Rectification