Common Set Operations

Common Set Operations

This section of the manual describes the basic symbolic types describing operations between sets.

Cartesian Product

Binary Cartesian Product

CartesianProduct{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents a Cartesian product of two convex sets.

Fields

  • X – first convex set
  • Y – second convex set

Notes

The Cartesian product of three elements is obtained recursively. See also CartesianProductArray for an implementation of a Cartesian product of many sets without recursion, instead using an array.

The EmptySet is the absorbing element for CartesianProduct.

Constructors:

  • CartesianProduct{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}}(X1::S1, X2::S2) – default constructor

  • CartesianProduct(Xarr::Vector{S}) where {S<:LazySet} – constructor from an array of convex sets

source
LinearAlgebra.:×Method.
×

Alias for the binary Cartesian product.

source
Base.:*Method.
    *(X::LazySet, Y::LazySet)

Alias for the binary Cartesian product.

source
LazySets.dimMethod.
dim(cp::CartesianProduct)::Int

Return the dimension of a Cartesian product.

Input

  • cp – Cartesian product

Output

The ambient dimension of the Cartesian product.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, cp::CartesianProduct{N}) where {N<:Real}

Return the support function of a Cartesian product.

Input

  • d – direction
  • cp – Cartesian product

Output

The support function in the given direction. If the direction has norm zero, the result depends on the wrapped sets.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cp::CartesianProduct{N}) where {N<:Real}

Return the support vector of a Cartesian product.

Input

  • d – direction
  • cp – Cartesian product

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the wrapped sets.

source
LazySets.isboundedMethod.
isbounded(cp::CartesianProduct)::Bool

Determine whether a Cartesian product is bounded.

Input

  • cp – Cartesian product

Output

true iff both wrapped sets are bounded.

source
Base.:∈Method.
∈(x::AbstractVector{N}, cp::CartesianProduct{N})::Bool where {N<:Real}

Check whether a given point is contained in a Cartesian product.

Input

  • x – point/vector
  • cp – Cartesian product

Output

true iff $x ∈ cp$.

source
Base.isemptyMethod.
isempty(cp::CartesianProduct)::Bool

Return if a Cartesian product is empty or not.

Input

  • cp – Cartesian product

Output

true iff any of the sub-blocks is empty.

source
constraints_list(cp::CartesianProduct{N}
                )::Vector{LinearConstraint{N}} where {N<:Real}

Return the list of constraints of a (polytopic) Cartesian product.

Input

  • cp – Cartesian product

Output

A list of constraints.

source
vertices_list(cp::CartesianProduct{N})::Vector{Vector{N}} where {N<:Real}

Return the list of vertices of a (polytopic) Cartesian product.

Input

  • cp – Cartesian product

Output

A list of vertices.

Algorithm

We assume that the underlying sets are polytopic. Then the high-dimensional set of vertices is just the Cartesian product of the low-dimensional sets of vertices.

source

Inherited from LazySet:

$n$-ary Cartesian Product

CartesianProductArray{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the Cartesian product of a finite number of convex sets.

Fields

  • array – array of sets

Notes

The EmptySet is the absorbing element for CartesianProductArray.

Constructors:

  • CartesianProductArray(array::Vector{<:LazySet}) – default constructor

  • CartesianProductArray([n]::Int=0, [N]::Type=Float64) – constructor for an empty product with optional size hint and numeric type

source
LazySets.dimMethod.
dim(cpa::CartesianProductArray)::Int

Return the dimension of a Cartesian product of a finite number of convex sets.

Input

  • cpa – Cartesian product array

Output

The ambient dimension of the Cartesian product of a finite number of convex sets.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, cp::CartesianProductArray{N}) where {N<:Real}

Return the support function of a Cartesian product array.

Input

  • d – direction
  • cpa – Cartesian product array

Output

The support function in the given direction. If the direction has norm zero, the result depends on the wrapped sets.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cpa::CartesianProductArray{N}) where {N<:Real}

Support vector of a Cartesian product array.

Input

  • d – direction
  • cpa – Cartesian product array

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the product sets.

source
LazySets.isboundedMethod.
isbounded(cpa::CartesianProductArray)::Bool

Determine whether a Cartesian product of a finite number of convex sets is bounded.

Input

  • cpa – Cartesian product of a finite number of convex sets

Output

true iff all wrapped sets are bounded.

source
Base.:∈Method.
∈(x::AbstractVector{N}, cpa::CartesianProductArray{N}
 )::Bool where {N<:Real}

Check whether a given point is contained in a Cartesian product of a finite number of sets.

Input

  • x – point/vector
  • cpa – Cartesian product array

Output

true iff $x ∈ \text{cpa}$.

source
Base.isemptyMethod.
isempty(cpa::CartesianProductArray)::Bool

Return if a Cartesian product is empty or not.

Input

  • cp – Cartesian product

Output

true iff any of the sub-blocks is empty.

source
constraints_list(cpa::CartesianProductArray{N}
                )::Vector{LinearConstraint{N}} where {N<:Real}

Return the list of constraints of a (polytopic) Cartesian product of a finite number of sets.

Input

  • cpa – Cartesian product array

Output

A list of constraints.

source
vertices_list(cpa::CartesianProductArray{N}
             )::Vector{Vector{N}} where {N<:Real}

Return the list of vertices of a (polytopic) Cartesian product of a finite number of sets.

Input

  • cpa – Cartesian product array

Output

A list of vertices.

Algorithm

We assume that the underlying sets are polytopic. Then the high-dimensional set of vertices is just the Cartesian product of the low-dimensional sets of vertices.

source
LazySets.arrayMethod.
array(cpa::CartesianProductArray{N, S}
     )::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Cartesian product of a finite number of convex sets.

Input

  • cpa – Cartesian product array

Output

The array of a Cartesian product of a finite number of convex sets.

source
array(cha::ConvexHullArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a convex hull of a finite number of convex sets.

Input

  • cha – convex hull array

Output

The array of a convex hull of a finite number of convex sets.

source
array(ia::IntersectionArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of an intersection of a finite number of convex sets.

Input

  • ia – intersection of a finite number of convex sets

Output

The array of an intersection of a finite number of convex sets.

source
array(msa::MinkowskiSumArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Minkowski sum of a finite number of convex sets.

Input

  • msa – Minkowski sum array

Output

The array of a Minkowski sum of a finite number of convex sets.

source
array(cms::CacheMinkowskiSum{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a caching Minkowski sum.

Input

  • cms – caching Minkowski sum

Output

The array of a caching Minkowski sum.

source
array(cup::UnionSetArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

The array that holds the union of a finite number of convex sets.

source

Inherited from LazySet:

Convex Hull

Binary Convex Hull

ConvexHull{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents the convex hull of the union of two convex sets.

Fields

  • X – convex set
  • Y – convex set

Notes

The EmptySet is the neutral element for ConvexHull.

Examples

Convex hull of two 100-dimensional Euclidean balls:

julia> b1, b2 = Ball2(zeros(100), 0.1), Ball2(4*ones(100), 0.2);

julia> c = ConvexHull(b1, b2);

julia> typeof(c)
ConvexHull{Float64,Ball2{Float64},Ball2{Float64}}
source
LazySets.CHType.
CH

Alias for ConvexHull.

source
LazySets.dimMethod.
dim(ch::ConvexHull)::Int

Return the dimension of a convex hull of two convex sets.

Input

  • ch – convex hull of two convex sets

Output

The ambient dimension of the convex hull of two convex sets.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, ch::ConvexHull{N}) where {N<:Real}

Return the support function of a convex hull of two convex sets in a given direction.

Input

  • d – direction
  • ch – convex hull of two convex sets

Output

The support function of the convex hull in the given direction.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, ch::ConvexHull{N}) where {N<:Real}

Return the support vector of a convex hull of two convex sets in a given direction.

Input

  • d – direction
  • ch – convex hull of two convex sets

Output

The support vector of the convex hull in the given direction.

source
LazySets.isboundedMethod.
isbounded(ch::ConvexHull)::Bool

Determine whether a convex hull of two convex sets is bounded.

Input

  • ch – convex hull of two convex sets

Output

true iff both wrapped sets are bounded.

source
Base.isemptyMethod.
isempty(ch::ConvexHull)::Bool

Return if a convex hull of two convex sets is empty or not.

Input

  • ch – convex hull

Output

true iff both wrapped sets are empty.

source

Inherited from LazySet:

$n$-ary Convex Hull

ConvexHullArray{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the symbolic convex hull of a finite number of convex sets.

Fields

  • array – array of sets

Notes

The EmptySet is the neutral element for ConvexHullArray.

Constructors:

  • ConvexHullArray(array::Vector{<:LazySet}) – default constructor

  • ConvexHullArray([n]::Int=0, [N]::Type=Float64) – constructor for an empty hull with optional size hint and numeric type

Examples

Convex hull of 100 two-dimensional balls whose centers follows a sinusoidal:

julia> b = [Ball2([2*pi*i/100, sin(2*pi*i/100)], 0.05) for i in 1:100];

julia> c = ConvexHullArray(b);
source
CHArray

Alias for ConvexHullArray.

source
LazySets.dimMethod.
dim(cha::ConvexHullArray)::Int

Return the dimension of the convex hull of a finite number of convex sets.

Input

  • cha – convex hull array

Output

The ambient dimension of the convex hull of a finite number of convex sets.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, cha::ConvexHullArray{N}) where {N<:Real}

Return the support function of a convex hull array in a given direction.

Input

  • d – direction
  • cha – convex hull array

Output

The support function of the convex hull array in the given direction.

Algorithm

This algorihm calculates the maximum over all $ρ(d, X_i)$ where the $X_1, …, X_k$ are the sets in the array cha.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cha::ConvexHullArray{N}) where {N<:Real}

Return the support vector of a convex hull array in a given direction.

Input

  • d – direction
  • cha – convex hull array
source
LazySets.isboundedMethod.
isbounded(cha::ConvexHullArray)::Bool

Determine whether a convex hull of a finite number of convex sets is bounded.

Input

  • cha – convex hull of a finite number of convex sets

Output

true iff all wrapped sets are bounded.

source
LazySets.arrayMethod.
array(cpa::CartesianProductArray{N, S}
     )::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Cartesian product of a finite number of convex sets.

Input

  • cpa – Cartesian product array

Output

The array of a Cartesian product of a finite number of convex sets.

source
array(cha::ConvexHullArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a convex hull of a finite number of convex sets.

Input

  • cha – convex hull array

Output

The array of a convex hull of a finite number of convex sets.

source
array(ia::IntersectionArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of an intersection of a finite number of convex sets.

Input

  • ia – intersection of a finite number of convex sets

Output

The array of an intersection of a finite number of convex sets.

source
array(msa::MinkowskiSumArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Minkowski sum of a finite number of convex sets.

Input

  • msa – Minkowski sum array

Output

The array of a Minkowski sum of a finite number of convex sets.

source
array(cms::CacheMinkowskiSum{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a caching Minkowski sum.

Input

  • cms – caching Minkowski sum

Output

The array of a caching Minkowski sum.

source
array(cup::UnionSetArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

The array that holds the union of a finite number of convex sets.

source
Base.isemptyMethod.
isempty(cha::ConvexHullArray)::Bool

Return if a convex hull array is empty or not.

Input

  • cha – convex hull array

Output

true iff all wrapped sets are empty.

source

Inherited from LazySet:

Convex Hull Algorithms

LazySets.convex_hullFunction.
convex_hull(P1::HPoly{N}, P2::HPoly{N};
           [backend]=default_polyhedra_backend(P1, N)) where {N}

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

Input

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

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
convex_hull(P::VPolygon{N}, Q::VPolygon{N};
            [algorithm]::String="monotone_chain")::VPolygon{N} where {N<:Real}

Return the convex hull of two polygons in vertex representation.

Input

  • P – polygon in vertex representation
  • Q – another polygon in vertex representation
  • algorithm – (optional, default: "monotone_chain") the algorithm used to compute the convex hull

Output

A new polygon such that its vertices are the convex hull of the given two polygons.

Algorithm

A convex hull algorithm is used to compute the convex hull of the vertices of the given input polygons P and Q; see ?convex_hull for details on the available algorithms. The vertices of the output polygon are sorted in counter-clockwise fashion.

source
convex_hull(P1::VPolytope{N}, P2::VPolytope{N};
            [backend]=default_polyhedra_backend(P1, N)) where {N}

Compute the convex hull of the set union of two polytopes in V-representation.

Input

  • P1 – polytope
  • P2 – another polytope
  • backend – (optional, default: default_polyhedra_backend(P1, N)) the polyhedral computations backend, see Polyhedra's documentation for further information

Output

The VPolytope 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.

source
convex_hull(points::Vector{VN};
            [algorithm]=default_convex_hull_algorithm(points),
            [backend]=nothing
            )::Vector{VN} where {N<:Real, VN<:AbstractVector{N}}

Compute the convex hull of the given points.

Input

  • points – list of vectors
  • algorithm – (optional, default: depends on the dimension) the convex hull algorithm, see valid options below
  • backend – (optional, default: "nothing") polyhedral computation backend for higher-dimensional point sets

Output

The convex hull as a list of vectors with the coordinates of the points.

Algorithm

A pre-processing step treats the cases with 0, 1 and 2 points in any dimension. For more than 3 points, the algorithm used depends on the dimension.

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, that gives access to libraries such as CDDLib and ConvexHull.jl. These libraries can be chosen from 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)
Array{Array{Float64,1},1}

Plot both the random points and the computed convex hull polygon:

julia> using Plots;

julia> plot([Tuple(pi) for pi in points], seriestype=:scatter);

julia> plot!(VPolygon(hull), alpha=0.2);
source
LazySets.right_turnFunction.
right_turn(O::AbstractVector{N}, A::AbstractVector{N}, B::AbstractVector{N}
          )::N where {N<:Real}

Determine if the acute angle defined by the three points O, A, B in the plane is a right turn (counter-clockwise) with respect to the center O.

Input

  • O – 2D center point
  • A – 2D one point
  • B – 2D another point

Output

Scalar representing the rotation.

Algorithm

The cross product is used to determine the sense of rotation. If the result is 0, the points are collinear; if it is positive, the three points constitute a positive angle of rotation around O from A to B; otherwise they constitute a negative angle.

source
monotone_chain!(points::Vector{VN}; sort::Bool=true
               )::Vector{VN} where {N<:Real, VN<:AbstractVector{N}}

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

Input

  • points – list of 2D vectors; is sorted in-place inside this function
  • 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 function 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

Intersection

Binary Intersection

Intersection{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents the intersection of two convex sets.

Fields

  • X – convex set
  • Y – convex set
  • cache – internal cache for avoiding recomputation; see IntersectionCache

Examples

Create an expression, $Z$, which lazily represents the intersection of two squares $X$ and $Y$:

julia> X, Y = BallInf([0,0.], 0.5), BallInf([1,0.], 0.65);

julia> Z = X ∩ Y;

julia> typeof(Z)
Intersection{Float64,BallInf{Float64},BallInf{Float64}}

julia> dim(Z)
2

We can check if the intersection is empty with isempty:

julia> isempty(Z)
false

Do not confuse Intersection with the concrete operation, which is computed with the lowercase intersection function:

julia> W = intersection(X, Y)
Hyperrectangle{Float64}([0.425, 0.0], [0.075, 0.5])
source
Base.:∩Method.

Alias for Intersection.

source
LazySets.dimMethod.
dim(cap::Intersection)::Int

Return the dimension of an intersection of two convex sets.

Input

  • cap – intersection of two convex sets

Output

The ambient dimension of the intersection of two convex sets.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, cap::Intersection{N}) where {N<:Real}

Return an upper bound on the support function of the intersection of two convex sets in a given direction.

Input

  • d – direction
  • cap – intersection of two convex sets

Output

An uper bound on the support function in the given direction.

Algorithm

The support function of an intersection of $X$ and $Y$ is upper bounded by the minimum of the support functions of $X$ and $Y$.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N},
  cap::Intersection{N, S1, S2};
  [algorithm]::String="line_search",
  [kwargs...]) where {N<:Real,
                      S1<:LazySet{N},
                      S2<:Union{HalfSpace{N}, Hyperplane{N}, Line{N}}}

Return the support function of the intersection of a compact set and a half-space/hyperplane/line in a given direction.

Input

  • d – direction

  • cap – lazy intersection of a compact set and a half-space/hyperplane/ line

  • algorithm – (optional, default: "line_search"): the algorithm to calculate the support function; valid options are:

    • "line_search" – solve the associated univariate optimization problem using a line search method (either Brent or the Golden Section method)
    • "projection" – only valid for intersection with a hyperplane; evaluates the support function by reducing the problem to the 2D intersection of a rank 2 linear transformation of the given compact set in the plane generated by the given direction d and the hyperplane's normal vector n
    • "simple" – take the $\min$ of the support function evaluation of each operand

Output

The scalar value of the support function of the set cap in the given direction.

Notes

It is assumed that the set cap.X is compact.

Any additional number of arguments to the algorithm backend can be passed as keyword arguments.

Algorithm

The algorithms are based on solving the associated optimization problem

\[\min_\{ λ ∈ D_h \} ρ(ℓ - λa, X) + λb.\]

where $D_h = \{ λ : λ ≥ 0 \}$ if $H$ is a half-space or $D_h = \{ λ : λ ∈ \mathbb{R} \}$ if $H$ is a hyperplane.

For additional information we refer to:

source
LazySets.ρMethod.
ρ(d::AbstractVector{N},
  cap::Intersection{N, S1, S2};
  kwargs...) where {N<:Real, S1<:LazySet{N}, S2<:AbstractPolytope{N}}

Return an upper bound of the intersection between a compact set and a polytope along a given direction.

Input

  • d – direction
  • cap – intersection of a compact set and a polytope
  • kwargs – additional arguments that are passed to the support function algorithm

Output

An upper bound of the support function of the given intersection.

Algorithm

The idea is to solve the univariate optimization problem ρ(di, X ∩ Hi) for each half-space in the set P and then take the minimum. This gives an overapproximation of the exact support function.

This algorithm is inspired from G. Frehse, R. Ray. Flowpipe-Guard Intersection for Reachability Computations with Support Functions.

Notes

This method relies on having available the constraints_list of the polytope P.

This method of overapproximation can return a non-empty set even if the original intersection is empty.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cap::Intersection{N}) where {N<:Real}

Return the support vector of an intersection of two convex sets in a given direction.

Input

  • d – direction
  • cap – intersection of two convex sets

Output

The support vector in the given direction.

source
LazySets.isboundedMethod.
isbounded(cap::Intersection)::Bool

Determine whether an intersection of two convex sets is bounded.

Input

  • cap – intersection of two convex sets

Output

true iff the intersection is bounded.

Algorithm

We first check if any of the wrapped sets is bounded. Otherwise, we check boundedness via isbounded_unit_dimensions.

source
Base.isemptyMethod.
isempty(cap::Intersection)::Bool

Return if the intersection is empty or not.

Input

  • cap – intersection of two convex sets

Output

true iff the intersection is empty.

Notes

The result will be cached, so a second query will be fast.

source
Base.:∈Method.
∈(x::AbstractVector{N}, cap::Intersection{N})::Bool where {N<:Real}

Check whether a given point is contained in an intersection of two convex sets.

Input

  • x – point/vector
  • cap – intersection of two convex sets

Output

true iff $x ∈ cap$.

source
constraints_list(cap::Intersection{N}) where {N<:Real}

Return the list of constraints of an intersection of two (polyhedral) sets.

Input

  • cap – intersection of two (polyhedral) sets

Output

The list of constraints of the intersection.

Notes

We assume that the underlying sets are polyhedral, i.e., offer a method constraints_list.

Algorithm

We create the polyhedron by taking the intersection of the constraints_lists of the sets and remove redundant constraints.

This function ignores the boolean output from the in-place remove_redundant_constraints!, which may inform the user that the constraints are infeasible. In that case, the list of constraints at the moment when the infeasibility was detected is returned.

source
isempty_known(cap::Intersection)

Ask whether the status of emptiness is known.

Input

  • cap – intersection of two convex sets

Output

true iff the emptiness status is known. In this case, isempty(cap) can be used to obtain the status.

source
set_isempty!(cap::Intersection, isempty::Bool)

Set the status of emptiness in the cache.

Input

  • cap – intersection of two convex sets
  • isempty – new status of emptiness
source
LazySets.swapMethod.
swap(cap::Intersection{N, S1, S2})::Intersection{N} where {N<:Real, S1, S2}

Return a new Intersection object with the arguments swapped.

Input

  • cap – intersection of two convex sets

Output

A new Intersection object with the arguments swapped. The old cache is shared between the old and new objects.

Notes

The advantage of using this function instead of manually swapping the arguments is that the cache is shared.

source
use_precise_ρ(cap::Intersection{N})::Bool where {N<:Real}

Determine whether a precise algorithm for computing $ρ$ shall be applied.

Input

  • cap – intersection of two convex sets

Output

true if a precise algorithm shall be applied.

Notes

The default implementation always returns true.

If the result is false, a coarse approximation of the support function is returned.

This function can be overwritten by the user to control the policy.

source
LazySets._line_searchFunction.
_line_search(ℓ, X, H; [kwargs...])

Given a compact and convex set $X$ and a halfspace $H = \{x: a^T x ≤ b \}$ or a hyperplane $H = \{x: a^T x = b \}$, calculate:

\[\min_\{ λ ∈ D_h \} ρ(ℓ - λa, X) + λb.\]

where $D_h = \{ λ : λ ≥ 0 \}$ if $H$ is a half-space or $D_h = \{ λ : λ ∈ \mathbb{R} \}$ if $H$ is a hyperplane.

Input

  • – direction
  • X – set
  • H – halfspace or hyperplane

Output

The tuple (fmin, λmin), where fmin is the minimum value of the function $f(λ) = ρ(ℓ - λa) + λb$ over the feasible set $λ ≥ 0$, and $λmin$ is the minimizer.

Notes

This function requires the Optim package, and relies on the univariate optimization interface Optim.optimize(...).

Additional arguments to the optimize backend can be passed as keyword arguments. The default method is Optim.Brent().

Examples

julia> X = Ball1(zeros(2), 1.0);

julia> H = HalfSpace([-1.0, 0.0], -1.0); # x >= 0 

julia> using Optim

julia> using LazySets: _line_search

julia> _line_search([1.0, 0.0], X, H) # uses Brent's method by default
(1.0, 999999.9849478417)

We can specify the upper bound in Brent's method:

julia> _line_search([1.0, 0.0], X, H, upper=1e3)
(1.0, 999.9999849478418)

Instead of using Brent, we use the Golden Section method:

julia> _line_search([1.0, 0.0], X, H, upper=1e3, method=GoldenSection())
(1.0, 381.9660112501051)
source
LazySets._projectionFunction.
_projection(ℓ, X, H::Union{Hyperplane{N}, Line{N}};
            [lazy_linear_map]=false,
            [lazy_2d_intersection]=true,
            [algorithm_2d_intersection]=nothing,
            [kwargs...]) where {N}

Given a compact and convex set $X$ and a hyperplane $H = \{x: n ⋅ x = γ \}$, calculate the support function of the intersection between the rank-2 projection $Π_{nℓ} X$ and the line $Lγ = \{(x, y): x = γ \}$.

Input

  • – direction
  • X – set
  • H – hyperplane
  • lazy_linear_map – (optional, default: false) to perform the projection lazily or concretely
  • lazy_2d_intersection – (optional, default: true) to perform the 2D intersection between the projected set and the line lazily or concretely
  • algorithm_2d_intersection – (optional, default: nothing) if given, fixes the support function algorithm used for the intersection in 2D; otherwise the default is implied

Output

The support function of $X ∩ H$ along direction $ℓ$.

Algorithm

This projection method is based on Prop. 8.2, page 103, C. Le Guernic. Reachability Analysis of Hybrid Systems with Linear Continuous Dynamics, PhD thesis.

In the original algorithm, Section 8.2 of Le Guernic's thesis, the linear map is performed concretely and the intersection is performed lazily (these are the default options in this algorithm, but here the four combinations are available). If the set $X$ is a zonotope, its concrete projection is again a zonotope (sometimes called "zonogon"). The intersection between this zonogon and the line can be taken efficiently in a lazy way (see Section 8.2.2 of Le Guernic's thesis), if one uses dispatch on ρ(y_dir, Sℓ⋂Lγ; kwargs...) given that Sℓ is itself a zonotope.

Notes

This function depends itself on the calculation of the support function of another set in two dimensions. Obviously one doesn't want to use again algorithm="projection" for this second calculation. The option algorithm_2d_intersection is such that, if it is not given, the default support function algorithm is used (e.g. "line_search"). You can still pass additional arguments to the "line_search" backend through the kwargs.

source
linear_map(M::AbstractMatrix{N}, cap::Intersection{N}) where {N}

Return the concrete linear map of a lazy intersection.

Input

  • M – matrix
  • cap – lazy intersection

Output

The set obtained by applying the given linear map to the lazy intersection.

Notes

This function relies on computing cap concretely (i.e. as a set representation), and then applying the linear map.

source

Inherited from LazySet:

Intersection cache

IntersectionCache

Container for information cached by a lazy Intersection object.

Fields

  • isempty – is the intersection empty? There are three possible states, encoded as Int8 values -1, 0, 1:

    • $-1$ - it is currently unknown whether the intersection is empty or not
    • $0$ - intersection is not empty
    • $1$ - intersection is empty
source

$n$-ary Intersection

IntersectionArray{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the intersection of a finite number of convex sets.

Fields

  • array – array of convex sets

Notes

This type assumes that the dimensions of all elements match.

The EmptySet is the absorbing element for IntersectionArray.

Constructors:

  • IntersectionArray(array::Vector{<:LazySet}) – default constructor

  • IntersectionArray([n]::Int=0, [N]::Type=Float64) – constructor for an empty sum with optional size hint and numeric type

source
LazySets.dimMethod.
dim(ia::IntersectionArray)::Int

Return the dimension of an intersection of a finite number of sets.

Input

  • ia – intersection of a finite number of convex sets

Output

The ambient dimension of the intersection of a finite number of sets.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, ia::IntersectionArray{N})::Vector{N} where {N<:Real}

Return the support vector of an intersection of a finite number of sets in a given direction.

Input

  • d – direction
  • ia – intersection of a finite number of convex sets

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the individual sets.

source
LazySets.isboundedMethod.
isbounded(ia::IntersectionArray)::Bool

Determine whether an intersection of a finite number of convex sets is bounded.

Input

  • ia – intersection of a finite number of convex sets

Output

true iff the intersection is bounded.

Algorithm

We first check if any of the wrapped sets is bounded. Otherwise, we check boundedness via isbounded_unit_dimensions.

source
Base.:∈Method.
∈(x::AbstractVector{N}, ia::IntersectionArray{N})::Bool where {N<:Real}

Check whether a given point is contained in an intersection of a finite number of convex sets.

Input

  • x – point/vector
  • ia – intersection of a finite number of convex sets

Output

true iff $x ∈ ia$.

source
LazySets.arrayMethod.
array(cpa::CartesianProductArray{N, S}
     )::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Cartesian product of a finite number of convex sets.

Input

  • cpa – Cartesian product array

Output

The array of a Cartesian product of a finite number of convex sets.

source
array(cha::ConvexHullArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a convex hull of a finite number of convex sets.

Input

  • cha – convex hull array

Output

The array of a convex hull of a finite number of convex sets.

source
array(ia::IntersectionArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of an intersection of a finite number of convex sets.

Input

  • ia – intersection of a finite number of convex sets

Output

The array of an intersection of a finite number of convex sets.

source
array(msa::MinkowskiSumArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Minkowski sum of a finite number of convex sets.

Input

  • msa – Minkowski sum array

Output

The array of a Minkowski sum of a finite number of convex sets.

source
array(cms::CacheMinkowskiSum{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a caching Minkowski sum.

Input

  • cms – caching Minkowski sum

Output

The array of a caching Minkowski sum.

source
array(cup::UnionSetArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

The array that holds the union of a finite number of convex sets.

source
constraints_list(ia::IntersectionArray{N}) where {N<:Real}

Return the list of constraints of an intersection of a finite number of (polyhedral) sets.

Input

  • ia – intersection of a finite number of (polyhedral) sets

Output

The list of constraints of the intersection.

Notes

We assume that the underlying sets are polyhedral, i.e., offer a method constraints_list.

Algorithm

We create the polyhedron from the constraints_lists of the sets and remove redundant constraints.

source

Inherited from LazySet:

Minkowski Sum

Binary Minkowski Sum

MinkowskiSum{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents the Minkowski sum of two convex sets.

Fields

  • X – first convex set
  • Y – second convex set

Notes

The ZeroSet is the neutral element and the EmptySet is the absorbing element for MinkowskiSum.

source
LazySets.:⊕Method.
⊕(X::LazySet, Y::LazySet)

Unicode alias constructor ⊕ (oplus) for the lazy Minkowski sum operator.

source
Base.:+Method.
X + Y

Convenience constructor for Minkowski sum.

Input

  • X – a convex set
  • Y – another convex set

Output

The symbolic Minkowski sum of $X$ and $Y$.

source
LazySets.dimMethod.
dim(ms::MinkowskiSum)::Int

Return the dimension of a Minkowski sum.

Input

  • ms – Minkowski sum

Output

The ambient dimension of the Minkowski sum.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, ms::MinkowskiSum{N}) where {N<:Real}

Return the support function of a Minkowski sum.

Input

  • d – direction
  • ms – Minkowski sum

Output

The support function in the given direction.

Algorithm

The support function in direction $d$ of the Minkowski sum of two sets $X$ and $Y$ is the sum of the support functions of $X$ and $Y$ in direction $d$.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, ms::MinkowskiSum{N}) where {N<:Real}

Return the support vector of a Minkowski sum.

Input

  • d – direction
  • ms – Minkowski sum

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the summand sets.

Algorithm

The support vector in direction $d$ of the Minkowski sum of two sets $X$ and $Y$ is the sum of the support vectors of $X$ and $Y$ in direction $d$.

source
LazySets.isboundedMethod.
isbounded(ms::MinkowskiSum)::Bool

Determine whether a Minkowski sum is bounded.

Input

  • ms – Minkowski sum

Output

true iff both wrapped sets are bounded.

source
Base.isemptyMethod.
isempty(ms::MinkowskiSum)::Bool

Return if a Minkowski sum is empty or not.

Input

  • ms – Minkowski sum

Output

true iff any of the wrapped sets are empty.

source

Inherited from LazySet:

$n$-ary Minkowski Sum

MinkowskiSumArray{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the Minkowski sum of a finite number of convex sets.

Fields

  • array – array of convex sets

Notes

This type assumes that the dimensions of all elements match.

The ZeroSet is the neutral element and the EmptySet is the absorbing element for MinkowskiSumArray.

Constructors:

  • MinkowskiSumArray(array::Vector{<:LazySet}) – default constructor

  • MinkowskiSumArray([n]::Int=0, [N]::Type=Float64) – constructor for an empty sum with optional size hint and numeric type

source
LazySets.dimMethod.
dim(msa::MinkowskiSumArray)::Int

Return the dimension of a Minkowski sum of a finite number of sets.

Input

  • msa – Minkowski sum array

Output

The ambient dimension of the Minkowski sum of a finite number of sets.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, msa::MinkowskiSumArray{N}) where {N<:Real}

Return the support function of a Minkowski sum array of a finite number of sets in a given direction.

Input

  • d – direction
  • msa – Minkowski sum array

Output

The support function in the given direction.

Algorithm

The support function of the Minkowski sum of sets is the sum of the support functions of each set.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, msa::MinkowskiSumArray{N}) where {N<:Real}

Return the support vector of a Minkowski sum of a finite number of sets in a given direction.

Input

  • d – direction
  • msa – Minkowski sum array

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the summand sets.

source
LazySets.isboundedMethod.
isbounded(msa::MinkowskiSumArray)::Bool

Determine whether a Minkowski sum of a finite number of convex sets is bounded.

Input

  • msa – Minkowski sum of a finite number of convex sets

Output

true iff all wrapped sets are bounded.

source
Base.isemptyMethod.
isempty(msa::MinkowskiSumArray)::Bool

Return if a Minkowski sum array is empty or not.

Input

  • msa – Minkowski sum array

Output

true iff any of the wrapped sets are empty.

source
LazySets.arrayMethod.
array(cpa::CartesianProductArray{N, S}
     )::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Cartesian product of a finite number of convex sets.

Input

  • cpa – Cartesian product array

Output

The array of a Cartesian product of a finite number of convex sets.

source
array(cha::ConvexHullArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a convex hull of a finite number of convex sets.

Input

  • cha – convex hull array

Output

The array of a convex hull of a finite number of convex sets.

source
array(ia::IntersectionArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of an intersection of a finite number of convex sets.

Input

  • ia – intersection of a finite number of convex sets

Output

The array of an intersection of a finite number of convex sets.

source
array(msa::MinkowskiSumArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Minkowski sum of a finite number of convex sets.

Input

  • msa – Minkowski sum array

Output

The array of a Minkowski sum of a finite number of convex sets.

source
array(cms::CacheMinkowskiSum{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a caching Minkowski sum.

Input

  • cms – caching Minkowski sum

Output

The array of a caching Minkowski sum.

source
array(cup::UnionSetArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

The array that holds the union of a finite number of convex sets.

source

Inherited from LazySet:

$n$-ary Minkowski Sum with cache

CacheMinkowskiSum{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the Minkowski sum of a finite number of convex sets. Support vector queries are cached.

Fields

  • array – array of convex sets
  • cache – cache of support vector query results

Notes

This type assumes that the dimensions of all elements match.

The ZeroSet is the neutral element and the EmptySet is the absorbing element for CacheMinkowskiSum.

The cache (field cache) is implemented as dictionary whose keys are directions and whose values are pairs (k, s) where k is the number of elements in the array array when the support vector was evaluated last time, and s is the support vector that was obtained. Thus this type assumes that array is not modified except by adding new sets at the end.

Constructors:

  • CacheMinkowskiSum(array::Vector{<:LazySet}) – default constructor

  • CacheMinkowskiSum([n]::Int=0, [N]::Type=Float64) – constructor for an empty sum with optional size hint and numeric type

source
LazySets.dimMethod.
dim(cms::CacheMinkowskiSum)::Int

Return the dimension of a caching Minkowski sum.

Input

  • cms – caching Minkowski sum

Output

The ambient dimension of the caching Minkowski sum.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cms::CacheMinkowskiSum{N}) where {N<:Real}

Return the support vector of a caching Minkowski sum in a given direction.

Input

  • d – direction
  • cms – caching Minkowski sum

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the summand sets.

Notes

The result is cached, i.e., any further query with the same direction runs in constant time. When sets are added to the caching Minkowski sum, the query is only performed for the new sets.

source
LazySets.isboundedMethod.
isbounded(cms::CacheMinkowskiSum)::Bool

Determine whether a caching Minkowski sum is bounded.

Input

  • cms – caching Minkowski sum

Output

true iff all wrapped sets are bounded.

source
Base.isemptyMethod.
isempty(cms::CacheMinkowskiSum)::Bool

Return if a caching Minkowski sum array is empty or not.

Input

  • cms – caching Minkowski sum

Output

true iff any of the wrapped sets are empty.

Notes

Forgotten sets cannot be checked anymore. Usually they have been empty because otherwise the support vector query should have crashed before. In that case, the caching Minkowski sum should not be used further.

source
LazySets.arrayMethod.
array(cpa::CartesianProductArray{N, S}
     )::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Cartesian product of a finite number of convex sets.

Input

  • cpa – Cartesian product array

Output

The array of a Cartesian product of a finite number of convex sets.

source
array(cha::ConvexHullArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a convex hull of a finite number of convex sets.

Input

  • cha – convex hull array

Output

The array of a convex hull of a finite number of convex sets.

source
array(ia::IntersectionArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of an intersection of a finite number of convex sets.

Input

  • ia – intersection of a finite number of convex sets

Output

The array of an intersection of a finite number of convex sets.

source
array(msa::MinkowskiSumArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a Minkowski sum of a finite number of convex sets.

Input

  • msa – Minkowski sum array

Output

The array of a Minkowski sum of a finite number of convex sets.

source
array(cms::CacheMinkowskiSum{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a caching Minkowski sum.

Input

  • cms – caching Minkowski sum

Output

The array of a caching Minkowski sum.

source
array(cup::UnionSetArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

The array that holds the union of a finite number of convex sets.

source
forget_sets!(cms::CacheMinkowskiSum)::Int

Tell a caching Minkowski sum to forget the stored sets (but not the support vectors). Only those sets are forgotten such that for each cached direction the support vector has been computed before.

Input

  • cms – caching Minkowski sum

Output

The number of sets that have been forgotten.

Notes

This function should only be used under the assertion that no new directions are queried in the future; otherwise such support vector results will be incorrect.

This implementation is optimistic and first tries to remove all sets. However, it also checks that for all cached directions the support vector has been computed before. If it finds that this is not the case, the implementation identifies the biggest index $k$ such that the above holds for the $k$ oldest sets, and then it only removes these. See the example below.

Examples

julia> x1 = BallInf(ones(3), 3.); x2 = Ball1(ones(3), 5.);

julia> cms1 = CacheMinkowskiSum(2); cms2 = CacheMinkowskiSum(2);

julia> d = ones(3);

julia> a1 = array(cms1); a2 = array(cms2);

julia> push!(a1, x1); push!(a2, x1);

julia> σ(d, cms1); σ(d, cms2);

julia> push!(a1, x2); push!(a2, x2);

julia> σ(d, cms1);

julia> idx1 = forget_sets!(cms1) # support vector was computed for both sets
2

julia> idx1 = forget_sets!(cms2) # support vector was only computed for first set
1
source

Inherited from LazySet:

Maps

Linear Map

LinearMap{N<:Real, S<:LazySet{N}, NM, MAT<:AbstractMatrix{NM}} <: LazySet{N}

Type that represents a linear transformation $M⋅S$ of a convex set $S$.

Fields

  • M – matrix/linear map
  • X – convex set

Notes

This type is parametric in the elements of the linear map, NM, which is independent of the numeric type of the target set (N). Typically NM = N, but there may be exceptions, e.g., if NM is an interval that holds numbers of type N, where N is a floating point number type such as Float64.

source
Base.:*Method.
    *(M::AbstractMatrix{N}, X::LazySet{N}) where {N<:Real}

Return the linear map of a convex set.

Input

  • M – matrix/linear map
  • X – convex set

Output

A lazy linear map, i.e. a LinearMap instance.

source
Base.:*Method.
    *(a::N, X::LazySet{N}) where {N<:Real}

Return a linear map of a convex set by a scalar value.

Input

  • a – scalar
  • X – convex set

Output

The linear map of the convex set.

source
Base.:*Method.
    *(a::N, lm::LM)::LM where {N<:Real, LM<:LinearMap{N}}

Return a linear map scaled by a scalar value.

Input

  • a – scalar
  • lm – linear map

Output

The scaled linear map.

source
Base.:*Method.
    *(M::AbstractMatrix{N}, Z::ZeroSet{N})::ZeroSet{N} where {N<:Real}

A linear map of a zero set, which is simplified to a zero set (the absorbing element).

Input

  • M – abstract matrix
  • Z – zero set

Output

The zero set with the output dimension of the linear map.

source
LazySets.dimMethod.
dim(lm::LinearMap)::Int

Return the dimension of a linear map.

Input

  • lm – linear map

Output

The ambient dimension of the linear map.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, lm::LinearMap{N}; kwargs...) where {N<:Real}

Return the support function of the linear map.

Input

  • d – direction
  • lm – linear map
  • kwargs – additional arguments that are passed to the support function algorithm

Output

The support function in the given direction. If the direction has norm zero, the result depends on the wrapped set.

Notes

If $L = M⋅S$, where $M$ is a matrix and $S$ is a convex set, it follows that $ρ(d, L) = ρ(M^T d, S)$ for any direction $d$.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, lm::LinearMap{N}) where {N<:Real}

Return the support vector of the linear map.

Input

  • d – direction
  • lm – linear map

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the wrapped set.

Notes

If $L = M⋅S$, where $M$ is a matrix and $S$ is a convex set, it follows that $σ(d, L) = M⋅σ(M^T d, S)$ for any direction $d$.

source
Base.:∈Method.
∈(x::AbstractVector{N}, lm::LinearMap{N})::Bool where {N<:Real}

Check whether a given point is contained in a linear map of a convex set.

Input

  • x – point/vector
  • lm – linear map of a convex set

Output

true iff $x ∈ lm$.

Algorithm

Note that $x ∈ M⋅S$ iff $M^{-1}⋅x ∈ S$. This implementation does not explicitly invert the matrix, which is why it also works for non-square matrices.

Examples

julia> lm = LinearMap([2.0 0.0; 0.0 1.0], BallInf([1., 1.], 1.));

julia> ∈([5.0, 1.0], lm)
false
julia> ∈([3.0, 1.0], lm)
true

An example with non-square matrix:

julia> B = BallInf(zeros(4), 1.);

julia> M = [1. 0 0 0; 0 1 0 0]/2;

julia> ∈([0.5, 0.5], M*B)
true
source
an_element(lm::LinearMap{N})::Vector{N} where {N<:Real}

Return some element of a linear map.

Input

  • lm – linear map

Output

An element in the linear map. It relies on the an_element function of the wrapped set.

source
LazySets.isboundedMethod.
isbounded(lm::LinearMap; cond_tol::Number=DEFAULT_COND_TOL)::Bool

Determine whether a linear map is bounded.

Input

  • lm – linear map
  • cond_tol – (optional) tolerance of matrix condition (used to check whether the matrix is invertible)

Output

true iff the linear map is bounded.

Algorithm

We first check if the matrix is zero or the wrapped set is bounded. If not, we perform a sufficient check whether the matrix is invertible. If the matrix is invertible, then the map being bounded is equivalent to the wrapped set being bounded, and hence the map is unbounded. Otherwise, we check boundedness via isbounded_unit_dimensions.

source
Base.isemptyMethod.
isempty(lm::LinearMap)::Bool

Return if a linear map is empty or not.

Input

  • lm – linear map

Output

true iff the wrapped set is empty.

source
vertices_list(lm::LinearMap{N})::Vector{Vector{N}} where {N<:Real}

Return the list of vertices of a (polyhedral) linear map.

Input

  • lm – linear map

Output

A list of vertices.

Algorithm

We assume that the underlying set X is polyhedral. Then the result is just the linear map applied to the vertices of X.

source
constraints_list(lm::LinearMap{N}) where {N<:Real}

Return the list of constraints of a (polyhedral) linear map.

Input

  • lm – linear map

Output

The list of constraints of the linear map.

Notes

We assume that the underlying set X is polyhedral, i.e., offers a method constraints_list(X).

Algorithm

We fall back to a concrete set representation and apply linear_map.

source
linear_map(M::AbstractMatrix{N}, lm::LinearMap{N}) where {N}

Return the linear map of a lazy linear map.

Input

  • M – matrix
  • lm – linear map

Output

The polytope representing the linear map of the lazy linear map of a set.

source
intersection(L::LinearMap{N}, S::LazySet{N}) where {N}

Return the intersection of a lazy linear map and a convex set.

Input

  • L – linear map
  • S – convex set

Output

The polytope obtained by the intersection of l.M * L.X and S.

source

Inherited from LazySet:

Exponential Map

ExponentialMap
dim(::ExponentialMap)
ρ(::AbstractVector{N}, ::ExponentialMap{N}) where {N<:Real}
σ(::AbstractVector{N}, ::ExponentialMap{N}) where {N<:Real}
∈(::AbstractVector{N}, ::ExponentialMap{N}) where {N<:Real}
isbounded(::ExponentialMap)
isempty(::ExponentialMap)
vertices_list(::ExponentialMap{N}) where {N<:Real}

Inherited from LazySet:

ExponentialProjectionMap
dim(::ExponentialProjectionMap)
σ(::AbstractVector{N}, ::ExponentialProjectionMap{N}) where {N<:Real}
isbounded(::ExponentialProjectionMap)
isempty(::ExponentialProjectionMap)

Inherited from LazySet:

SparseMatrixExp
*(::SparseMatrixExp{N}, ::LazySet{N}) where {N<:Real}
get_row(::SparseMatrixExp, ::Int)
ProjectionSparseMatrixExp{N<:Real}

Type that represents the projection of a sparse matrix exponential, i.e., $L⋅\exp(M)⋅R$ for a given sparse matrix $M$.

Fields

  • L – left multiplication matrix
  • E – sparse matrix exponential
  • R – right multiplication matrix
source
Base.:*Method.
    *(projspmexp::ProjectionSparseMatrixExp,
      X::LazySet)::ExponentialProjectionMap

Return the application of a projection of a sparse matrix exponential to a convex set.

Input

  • projspmexp – projection of a sparse matrix exponential
  • X – convex set

Output

The application of the projection of a sparse matrix exponential to the convex set.

source

Reset Map

ResetMap{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents a lazy reset map. A reset map is a special case of an affine map $A x + b, x ∈ X$ where the linear map $A$ is the identity matrix with zero entries in all reset dimensions, and the translation vector $b$ is zero in all other dimensions.

Fields

  • X – convex set
  • resets – resets (a mapping from an index to a new value)

Example

julia> X = BallInf([2.0, 2.0, 2.0], 1.0);

julia> r = Dict(1 => 4.0, 3 => 0.0);

julia> rm = ResetMap(X, r);

Here rm modifies the set X such that x1 is reset to 4 and x3 is reset to 0, while x2 is not modified. Hence rm is equivalent to the set Hyperrectangle([4.0, 2.0, 0.0], [0.0, 1.0, 0.0]), i.e., an axis-aligned line segment embedded in 3D.

The corresponding affine map $A x + b$ would be:

\[ egin{pmatrix} 0 & 0 & 0 \ 0 & 1 & 0 \ 0 & 0 & 0 nd{pmatrix} x + egin{pmatrix} 4 & 0 & 0 nd{pmatrix}\]

Use the function get_A (resp. get_b) to create the matrix A (resp. vector b) corresponding to a given reset map.

The (in this case unique) support vector of rm in direction ones(3) is:

julia> σ(ones(3), rm)
3-element Array{Float64,1}:
 4.0
 3.0
 0.0
source
LazySets.dimMethod.
dim(rm::ResetMap)

Return the dimension of a reset map.

Input

  • rm – reset map

Output

The dimension of a reset map.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, rm::ResetMap{N}) where {N<:Real}

Return the support function of a reset map.

Input

  • d – direction
  • rm – reset map

Output

The support function in the given direction.

Notes

We use the usual dot-product definition, but for unbounded sets we redefine the product between $0$ and $±∞$ as $0$; Julia returns NaN here.

julia> Inf * 0.0
NaN

See the discussion here.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, rm::ResetMap{N}) where {N<:Real}

Return the support vector of a reset map.

Input

  • d – direction
  • rm – reset map

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the wrapped set.

source
an_element(rm::ResetMap)

Return some element of a reset map.

Input

  • rm – reset map

Output

An element in the reset map. It relies on the an_element function of the wrapped set.

source
Base.isemptyMethod.
isempty(rm::ResetMap)::Bool

Return if a reset map is empty or not.

Input

  • rm – reset map

Output

true iff the wrapped set is empty.

source
LazySets.get_AMethod.
get_A(rm::ResetMap{N}) where {N<:Real}

Return the $A$ matrix of the affine map $A x + b, x ∈ X$ represented by a reset map.

Input

  • rm – reset map

Output

The (sparse) square matrix for the affine map $A x + b, x ∈ X$ represented by the reset map.

Algorithm

We construct the identity matrix and set all entries in the reset dimensions to zero.

source
LazySets.get_bMethod.
get_b(rm::ResetMap{N}) where {N<:Real}

Return the $b$ vector of the affine map $A x + b, x ∈ X$ represented by a reset map.

Input

  • rm – reset map

Output

The (sparse) vector for the affine map $A x + b, x ∈ X$ represented by the reset map. The vector contains the reset value for all reset dimensions, and is zero for all other dimensions.

source
constraints_list(rm::ResetMap{N}) where {N<:Real}

Return the list of constraints of a polytopic reset map.

Input

  • rm – reset map of a polytope

Output

The list of constraints of the reset map.

Notes

We assume that the underlying set X is a polytope, i.e., is bounded and offers a method constraints_list(X).

Algorithm

We fall back to constraints_list of a LinearMap of the A-matrix in the affine-map view of a reset map. Each reset dimension $i$ is projected to zero, expressed by two constraints for each reset dimension. Then it remains to shift these constraints to the new value.

For instance, if the dimension $5$ was reset to $4$, then there will be constraints $x₅ ≤ 0$ and $-x₅ ≤ 0$. We then modify the right-hand side of these constraints to $x₅ ≤ 4$ and $-x₅ ≤ -4$, respectively.

source
constraints_list(rm::ResetMap{N, S}) where
    {N<:Real, S<:AbstractHyperrectangle}

Return the list of constraints of a hyperrectangular reset map.

Input

  • rm – reset map of a hyperrectangular set

Output

The list of constraints of the reset map.

Algorithm

We iterate through all dimensions. If there is a reset, we construct the corresponding (flat) constraints. Otherwise, we construct the corresponding constraints of the underlying set.

source

Inherited from LazySet:

Symmetric Interval Hull

SymmetricIntervalHull{N<:Real, S<:LazySet{N}} <: AbstractHyperrectangle{N}

Type that represents the symmetric interval hull of a convex set.

Fields

  • X – convex set
  • cache – partial storage of already computed bounds, organized as mapping from dimension to tuples (bound, valid), where valid is a flag indicating if the bound entry has been computed

Notes

The symmetric interval hull can be computed with $2n$ support vector queries of unit vectors, where $n$ is the dimension of the wrapped set (i.e., two queries per dimension). When asking for the support vector for a direction $d$, one needs $2k$ such queries, where $k$ is the number of non-zero entries in $d$.

However, if one asks for many support vectors in a loop, the number of computations may exceed $2n$. To be most efficient in such cases, this type stores the intermediately computed bounds in the cache field.

source
LazySets.dimMethod.
dim(sih::SymmetricIntervalHull)::Int

Return the dimension of a symmetric interval hull of a convex set.

Input

  • sih – symmetric interval hull of a convex set

Output

The ambient dimension of the symmetric interval hull of a convex set.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, sih::SymmetricIntervalHull{N}) where {N<:Real}

Return the support vector of a symmetric interval hull of a convex set in a given direction.

Input

  • d – direction
  • sih – symmetric interval hull of a convex set

Output

The support vector of the symmetric interval hull of a convex set in the given direction. If the direction has norm zero, the origin is returned.

Algorithm

For each non-zero entry in d we need to either look up the bound (if it has been computed before) or compute it, in which case we store it for future queries. One such computation just asks for the support vector of the underlying set for both the positive and negative unit vector in the respective dimension.

source
LazySets.centerMethod.
center(sih::SymmetricIntervalHull{N})::Vector{N} where {N<:Real}

Return the center of a symmetric interval hull of a convex set.

Input

  • sih – symmetric interval hull of a convex set

Output

The origin.

source
radius_hyperrectangle(sih::SymmetricIntervalHull{N}
                     )::Vector{N} where {N<:Real}

Return the box radius of a symmetric interval hull of a convex set in every dimension.

Input

  • sih – symmetric interval hull of a convex set

Output

The box radius of the symmetric interval hull of a convex set.

Notes

This function computes the symmetric interval hull explicitly.

source
radius_hyperrectangle(sih::SymmetricIntervalHull{N},
                      i::Int)::N where {N<:Real}

Return the box radius of a symmetric interval hull of a convex set in a given dimension.

Input

  • sih – symmetric interval hull of a convex set
  • i – dimension of interest

Output

The radius in the given dimension. If it was computed before, this is just a look-up, otherwise it requires two support vector computations.

source

Inherited from LazySet:

Inherited from AbstractPolytope:

Inherited from AbstractPolyhedron:

Inherited from AbstractCentrallySymmetricPolytope:

Inherited from AbstractHyperrectangle:

Translation

Translation{N<:Real, VN<:AbstractVector{N}, S<:LazySet{N}} <: LazySet{N}

Type that represents a lazy translation.

The translation of set X along vector v is the map:

\[x ↦ x + v,\qquad x ∈ X\]

A translation is a special case of an affine map $A x + b, x ∈ X$ where the linear map $A$ is the identity matrix and the translation vector $b = v$.

Fields

  • X – convex set
  • v – vector that defines the translation

Example

julia> X = BallInf([2.0, 2.0, 2.0], 1.0);

julia> v = [1.0, 0.0, 0.0]; # translation along dimension 1

julia> tr = Translation(X, v);

julia> typeof(tr)
Translation{Float64,Array{Float64,1},BallInf{Float64}}

julia> tr.X
BallInf{Float64}([2.0, 2.0, 2.0], 1.0)

julia> tr.v
3-element Array{Float64,1}:
 1.0
 0.0
 0.0

The sum operator + is overloaded to create translations:

julia> X + v == Translation(X, v)
true

And so does the Minkowski sum operator, :

julia> X ⊕ v == Translation(X, v)
true

The translation of a translation is performed immediately:

julia> tr = (X+v)+v
Translation{Float64,Array{Float64,1},BallInf{Float64}}(BallInf{Float64}([2.0, 2.0, 2.0], 1.0), [2.0, 0.0, 0.0])

julia> tr.v
3-element Array{Float64,1}:
 2.0
 0.0
 0.0

The dimension of a translation is obtained with the dim function:

julia> dim(tr)
3

For the support vector (resp. support function) along vector d, use σ and ρ respectively:

julia> σ([1.0, 0.0, 0.0], tr)
3-element Array{Float64,1}:
 5.0
 3.0
 3.0

julia> ρ([1.0, 0.0, 0.0], tr)
5.0

See the docstring of each of these functions for details.

The an_element function is useful to obtain an element of a translation:

julia> e = an_element(tr)
3-element Array{Float64,1}:
 4.0
 2.0
 2.0

The lazy linear map of a translation is again a translation, since the following simplification rule applies: $M * (X⊕v) = (M*X) ⊕ (M*v)$:

julia> Q = Matrix(2.0I, 3, 3) * tr;

julia> Q isa Translation && Q.v == 2 * tr.v
true

Use the isempty method to query if the translation is empty; it falls back to the isempty method of the wrapped set:

julia> isempty(tr)
false

The list of constraints of the translation of a polyhedron (in general, a set whose constraints_list is available) can be computed from a lazy translation:

julia> constraints_list(tr)
6-element Array{HalfSpace{Float64},1}:
 HalfSpace{Float64}([1.0, 0.0, 0.0], 5.0)
 HalfSpace{Float64}([0.0, 1.0, 0.0], 3.0)
 HalfSpace{Float64}([0.0, 0.0, 1.0], 3.0)
 HalfSpace{Float64}([-1.0, -0.0, -0.0], -3.0)
 HalfSpace{Float64}([-0.0, -1.0, -0.0], -1.0)
 HalfSpace{Float64}([-0.0, -0.0, -1.0], -1.0)
source
Base.:+Method.
+(X::LazySet, v::AbstractVector)

Convenience constructor for a translation.

Input

  • X – convex set
  • v – vector

Output

The symbolic translation of $X$ along vector $v$.

source
LazySets.:⊕Method.
⊕(X::LazySet, v::AbstractVector)

Unicode alias constructor ⊕ (oplus) for the lazy translation operator.

source
LazySets.dimMethod.
dim(tr::Translation)::Int

Return the dimension of a translation.

Input

  • tr – translation

Output

The dimension of a translation.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, tr::Translation{N}) where {N<:Real}

Return the support function of a translation.

Input

  • d – direction
  • tr – translation

Output

The support function in the given direction.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, tr::Translation{N}) where {N<:Real}

Return the support vector of a translation.

Input

  • d – direction
  • tr – translation

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the wrapped set.

source
an_element(tr::Translation)

Return some element of a translation.

Input

  • tr – translation

Output

An element in the translation.

Notes

This function first asks for an_element function of the wrapped set, then translates this element according to the given translation vector.

source
Base.isemptyMethod.
isempty(tr::Translation)::Bool

Return if a translation is empty or not.

Input

  • tr – translation

Output

true iff the wrapped set is empty.

source
constraints_list(tr::Translation{N}, ::Val{true}) where {N<:Real}

Return the list of constraints of the translation of a set.

Input

  • tr – lazy translation of a polyhedron

Output

The list of constraints of the translation.

Notes

We assume that the set wrapped by the lazy translation X offers a method constraints_list(⋅).

Algorithm

Let the translation be defined by the set of points y such that y = x + v for all x ∈ X. Then, each defining halfspace a⋅x ≤ b is transformed to a⋅y ≤ b + a⋅v.

source
LazySets.LinearMapMethod.
LinearMap(M::AbstractMatrix{N}, tr::Translation{N}) where {N<:Real}

Return the lazy linear map of a translation.

Input

  • M – matrix
  • tr – translation

Output

The translation defined by the linear map.

Notes

This method defines the simplification rule: $M * (X⊕v) = (M*X) ⊕ (M*v)$, returning a new translation.

source

Union

Note that a union of convex sets is generally not convex. Hence these set types are not part of the convex-set family LazySet.

Binary Set Union

UnionSet{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}}

Type that represents the set union of two convex sets.

Fields

  • X – convex set
  • Y – convex set
source
Base.:∪Method.

Alias for UnionSet.

source
LazySets.dimMethod.
dim(cup::UnionSet)::Int

Return the dimension of the set union of two convex sets.

Input

  • cup – union of two convex sets

Output

The ambient dimension of the union of two convex sets.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cup::UnionSet{N}; [algorithm]="support_vector") where {N<:Real}

Return the support vector of the union of two convex sets in a given direction.

Input

  • d – direction
  • cup – union of two convex sets
  • algorithm – (optional, default: "supportvector"): the algorithm to compute the support vector; if "supportvector", use the support vector of each argument; if "support_function" use the support function of each argument and evaluate the support vector of only one of them

Output

The support vector in the given direction.

Algorithm

The support vector of the union of two convex sets $X$ and $Y$ can be obtained as the vector that maximizes the support function of either $X$ or $Y$, i.e. it is sufficient to find the $\argmax(ρ(d, X), ρ(d, Y)])$ and evaluate its support vector.

The default implementation, with option algorithm="support_vector", computes the support vector of $X$ and $Y$ and then compares the support function using a dot product. If it happens that the support function can be more efficiently computed (without passing through the support vector), consider using the alternative algorithm="support_function" implementation, which evaluates the support function of each set directly and then calls only the support vector of either $X$ or $Y$.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, cup::UnionSet{N}) where {N<:Real}

Return the support function of the union of two convex sets in a given direction.

Input

  • d – direction
  • cup – union of two convex sets

Output

The support function in the given direction.

Algorithm

The support function of the union of two convex sets $X$ and $Y$ is the maximum of the support functions of $X$ and $Y$.

source
an_element(cup::UnionSet{N})::Vector{N} where {N<:Real}

Return some element of a union of two convex sets.

Input

  • cup – union of two convex sets

Output

An element in the union of two convex sets.

Algorithm

We use an_element on the first wrapped set.

source
Base.:∈Method.
∈(x::AbstractVector{N}, cup::UnionSet{N})::Bool where {N<:Real}

Check whether a given point is contained in a union of two convex sets.

Input

  • x – point/vector
  • cup – union of two convex sets

Output

true iff $x ∈ cup$.

source
Base.isemptyMethod.
isempty(cup::UnionSet)::Bool

Check whether a union of two convex sets is empty.

Input

  • cup – union of two convex sets

Output

true iff the union is empty.

source
LazySets.isboundedMethod.
isbounded(cup::UnionSet)::Bool

Determine whether a union of two convex sets is bounded.

Input

  • cup – union of two convex sets

Output

true iff the union is bounded.

source

$n$-ary Set Union

UnionSetArray{N<:Real, S<:LazySet{N}}

Type that represents the set union of a finite number of convex sets.

Fields

  • array – array of convex sets
source
LazySets.dimMethod.
dim(cup::UnionSetArray)::Int

Return the dimension of the set union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

The ambient dimension of the union of a finite number of convex sets.

source
LazySets.arrayMethod.
array(cup::UnionSetArray{N, S})::Vector{S} where {N<:Real, S<:LazySet{N}}

Return the array of a union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

The array that holds the union of a finite number of convex sets.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cup::UnionSetArray{N}; [algorithm]="support_vector") where {N<:Real}

Return the support vector of the union of a finite number of convex sets in a given direction.

Input

  • d – direction
  • cup – union of a finite number of convex sets
  • algorithm – (optional, default: "supportvector"): the algorithm to compute the support vector; if "supportvector", use the support vector of each argument; if "support_function" use the support function of each argument and evaluate the support vector of only one of them

Output

The support vector in the given direction.

Algorithm

The support vector of the union of a finite number of convex sets $X₁, X₂, ...$ can be obtained as the vector that maximizes the support function, i.e. it is sufficient to find the $\argmax(ρ(d, X₂), ρ(d, X₂), ...])$ and evaluate its support vector.

The default implementation, with option algorithm="support_vector", computes the support vector of all $X₁, X₂, ...$ and then compares the support function using a dot product. If it happens that the support function can be more efficiently computed (without passing through the support vector), consider using the alternative algorithm="support_function" implementation, which evaluates the support function of each set directly and then calls only the support vector of one of the $Xᵢ$.

source
LazySets.ρMethod.
ρ(d::AbstractVector{N}, cup::UnionSetArray{N}) where {N<:Real}

Return the support function of the union of a finite number of convex sets in a given direction.

Input

  • d – direction
  • cup – union of a finite number of convex sets

Output

The support function in the given direction.

Algorithm

The support function of the union of a finite number of convex sets $X₁, X₂, ...$ can be obtained as the maximum of $ρ(d, X₂), ρ(d, X₂), ...$.

source
an_element(cup::UnionSetArray{N})::Vector{N} where {N<:Real}

Return some element of a union of a finite number of convex sets.

Input

  • cup – union of a finite number of convex sets

Output

An element in the union of a finite number of convex sets.

Algorithm

We use an_element on the first wrapped set.

source
Base.:∈Method.
∈(x::AbstractVector{N}, cup::UnionSetArray{N})::Bool where {N<:Real}

Check whether a given point is contained in a union of a finite number of convex sets.

Input

  • x – point/vector
  • cup – union of a finite number of convex sets

Output

true iff $x ∈ cup$.

source
Base.isemptyMethod.
isempty(cup::UnionSetArray)::Bool

Check whether a union of a finite number of convex sets is empty.

Input

  • cup – union of a finite number of convex sets

Output

true iff the union is empty.

source
LazySets.isboundedMethod.
isbounded(cup::UnionSetArray)::Bool

Determine whether a union of a finite number of convex sets is bounded.

Input

  • cup – union of a finite number of convex sets

Output

true iff the union is bounded.

source

Complement

Note that a complement of a convex set is generally not convex. Hence this set type is not part of the convex-set family LazySet.

Binary Set Union

Complement{N<:Real, S<:LazySet{N}}

Type that represents the complement of a convex set.

Fields

  • X – convex set

Notes

Since X is assumed to be closed, unless X is empty or the universe, its complement is open (i.e., not closed). Since X is assumed to be closed, unless X is empty, the universe, or a half-space, its complement is not convex.

The complement of the complement is the original set again.

Examples

julia> B = BallInf(zeros(2), 1.); C = Complement(B)
Complement{Float64,BallInf{Float64}}(BallInf{Float64}([0.0, 0.0], 1.0))

julia> Complement(C)
BallInf{Float64}([0.0, 0.0], 1.0)
source
LazySets.dimMethod.
dim(C::Complement)

Return the dimension of the complement of a convex set.

Input

  • C – complement of a convex set

Output

The dimension of the complement of a convex set.

source
Base.:∈Method.
∈(x::AbstractVector{N}, C::Complement{N})::Bool where {N<:Real}

Check whether a given point is contained in the complement of a convex set.

Input

  • x – point/vector
  • C – complement of a convex set

Output

true iff the vector is contained in the complement.

Algorithm

\[ x ∈ X^C ⟺ x ∉ X\]
source
Base.isemptyMethod.
isempty(C::Complement)::Bool

Return if the complement of a convex set is empty or not.

Input

  • C – complement of a convex set

Output

false unless the original set is universal.

Algorithm

We use the isuniversal method.

source