Approximations

Approximations

This section of the manual describes the Cartesian decomposition algorithms and the approximation of high-dimensional convex sets using projections.

Module Approximations.jl – polygonal approximation of convex sets through support vectors.

source

Cartesian Decomposition

decompose(S::LazySet{N},
          partition::AbstractVector{<:AbstractVector{Int}},
          block_options
         ) where {N<:Real}

Decompose a high-dimensional set into a Cartesian product of overapproximations of the projections over the specified subspaces.

Input

  • S – set
  • partition – vector of blocks (i.e., of vectors of integers) (see the Notes below)
  • block_options – mapping from block indices in partition to a corresponding overapproximation option; we only require access via [⋅] (but see also the Notes below)

Output

A CartesianProductArray containing the low-dimensional approximated projections.

Algorithm

For each block a specific project method is called, dispatching on the corresponding overapproximation option.

Notes

The argument partition requires some discussion. Typically, the list of blocks should form a partition of the set $\{1, \dots, n\}$ represented as a list of consecutive blocks, where $n$ is the ambient dimension of set S.

However, technically there is no problem if the blocks are not consecutive, blocks are missing, blocks occur more than once, or blocks are overlapping. This function will, however, stick to the order of blocks, so the resulting set must be interpreted with care in such cases. One use case is the need of a projection consisting of several blocks.

For convenience, the argument block_options can also be given as a single option instead of a mapping, which is then interpreted as the option for all blocks.

Examples

This function supports different options: one can specify the target set, the degree of accuracy, and template directions. These options are exemplified below, where we use the following example.

julia> using LazySets.Approximations: decompose

julia> S = Ball2(zeros(4), 1.);  # set to be decomposed (4D 2-norm unit ball)

julia> P2d = [1:2, 3:4];  # a partition with two blocks of size two

julia> P1d = [[1], [2], [3], [4]];  # a partition with four blocks of size one

Different set types

We can decompose using polygons in constraint representation:

julia> all([ai isa HPolygon for ai in array(decompose(S, P2d, HPolygon))])
true

For decomposition into 1D subspaces, we can use Interval:

julia> all([ai isa Interval for ai in array(decompose(S, P1d, Interval))])
true

However, if you need to specify different set types for different blocks, the interface presented so far does not apply. See the paragraph Advanced input for different block approximations below for how to do that.

Refining the decomposition I: $ε$-close approximation

The $ε$ option can be used to refine a decomposition, i.e., obtain a more accurate result. We use the Iterative refinement algorithm from the Approximations module.

To illustrate this, consider again the set S from above. We decompose into two 2D polygons. Using smaller $ε$ implies a better precision, thus more constraints in each 2D decomposition. In the following example, we look at the number of constraints in the first block.

julia> d(ε, bi) = array(decompose(S, P2d, (HPolygon => ε)))[bi]
d (generic function with 1 method)

julia> [length(constraints_list(d(ε, 1))) for ε in [Inf, 0.1, 0.01]]
3-element Array{Int64,1}:
  4
  8
 32

Refining the decomposition II: template polyhedra

Another way to refine a decomposition is by using template polyhedra. The idea is to specify a set of template directions and then to compute on each block the polytopic overapproximation obtained by evaluating the support function of the given input set over the template directions.

For example, octagonal 2D approximations of the set S are obtained with:

julia> using LazySets.Approximations: OctDirections

julia> B = decompose(S, P2d, OctDirections);

julia> length(B.array) == 2 && all(dim(bi) == 2 for bi in B.array)
true

See Template directions for the available template directions. Note that, in contrast to the polygonal $ε$-close approximation from above, this method can be applied to blocks of any size.

julia> B = decompose(S, [1:4], OctDirections);

julia> length(B.array) == 1 && dim(B.array[1]) == 4
true

Advanced input for different block approximations

Instead of defining the approximation option uniformly for each block, we can define different approximations for different blocks. The third argument has to be a mapping from block index (in the partition) to the corresponding approximation option.

For example:

julia> res = array(decompose(S, P2d, Dict(1 => Hyperrectangle, 2 => 0.1)));

julia> typeof(res[1]), typeof(res[2])
(Hyperrectangle{Float64,Array{Float64,1},Array{Float64,1}}, HPolygon{Float64})
source
project(S::LazySet{N},
        block::AbstractVector{Int},
        [::Nothing=nothing],
        [n]::Int=dim(S)
       ) where {N<:Real}

Project a high-dimensional set to a given block by using a concrete linear map.

Input

  • S – set
  • block – block structure - a vector with the dimensions of interest
  • nothing – (default: nothing) used for dispatch
  • n – (optional, default: dim(S)) ambient dimension of the set S

Output

A set representing the projection of the set S to block block.

Algorithm

We apply the function linear_map.

source
project(S::LazySet{N},
        block::AbstractVector{Int},
        set_type::Type{<:LinearMap},
        [n]::Int=dim(S)
       ) where {N<:Real}

Project a high-dimensional set to a given block by using a lazy linear map.

Input

  • S – set
  • block – block structure - a vector with the dimensions of interest
  • LinearMap – used for dispatch
  • n – (optional, default: dim(S)) ambient dimension of the set S

Output

A lazy LinearMap representing the projection of the set S to block block.

source
project(S::LazySet{N},
        block::AbstractVector{Int},
        set_type::Type{<:LazySet},
        [n]::Int=dim(S)
       ) where {N<:Real}

Project a high-dimensional set to a given block and set type, possibly involving an overapproximation.

Input

  • S – set
  • block – block structure - a vector with the dimensions of interest
  • set_type – target set type
  • n – (optional, default: dim(S)) ambient dimension of the set S

Output

A set of type set_type representing an overapproximation of the projection of S.

Algorithm

  1. Project the set S with M⋅S, where M is the identity matrix in the block

coordinates and zero otherwise.

  1. Overapproximate the projected lazy set using overapproximate and

set_type.

source
project(S::LazySet{N},
        block::AbstractVector{Int},
        set_type_and_precision::Pair{<:UnionAll, <:Real},
        [n]::Int=dim(S)
       ) where {N<:Real}

Project a high-dimensional set to a given block and set type with a certified error bound.

Input

  • S – set
  • block – block structure - a vector with the dimensions of interest
  • set_type_and_precision – pair (T, ε) of a target set type T and an error bound ε for approximation
  • n – (optional, default: dim(S)) ambient dimension of the set S

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

  1. Project the set S with M⋅S, where M is the identity matrix in the block

coordinates and zero otherwise.

  1. Overapproximate the projected lazy set with the given error bound ε.
source
project(S::LazySet{N},
        block::AbstractVector{Int},
        ε::Real,
        [n]::Int=dim(S)
       ) where {N<:Real}

Project a high-dimensional set to a given block and set type with a certified error bound.

Input

  • S – set
  • block – block structure - a vector with the dimensions of interest
  • ε – error bound for approximation
  • n – (optional, default: dim(S)) ambient dimension of the set S

Output

A set representing the epsilon-close approximation of the projection of S.

Algorithm

  1. Project the set S with M⋅S, where M is the identity matrix in the block

coordinates and zero otherwise.

  1. Overapproximate the projected lazy set with the given error bound ε.

The target set type is chosen automatically.

source
project(S::LazySet{N},
        block::AbstractVector{Int},
        directions::Type{<:AbstractDirections},
        [n]::Int
       ) where {N<:Real}

Project a high-dimensional set to a given block using template directions.

Input

  • S – set
  • block – block structure - a vector with the dimensions of interest
  • directions – template directions
  • n – (optional, default: dim(S)) ambient dimension of the set S

Output

The template direction approximation of the projection of S.

source
project(H::HalfSpace{N}, block::AbstractVector{Int})

Concrete projection of a half-space.

Input

  • H – set
  • block – block structure, a vector with the dimensions of interest

Output

A set representing the projection of the half-space H on the dimensions specified by block.

Algorithm

If the unconstrained dimensions of H are a subset of the block variables, the projection is applied to the normal direction of H. Otherwise, the projection results in the universal set.

The latter can be seen as follows. Without loss of generality consider a projection onto a single and constrained dimension $xₖ$ (projections in multiple dimensions can be modeled as repeated one-dimensional projections). We can write the projection as an existentially quantified linear constraint:

\[ ∃xₖ: a₁x₁ + … + aₖxₖ + … + aₙxₙ ≤ b\]

Since $aₖ ≠ 0$, there is always a value for $xₖ$ that satisfies the constraint for any valuation of the other variables.

Examples

Consider the half-space $x + y + 0⋅z ≤ 1$, whose ambient dimension is 3. The (trivial) projection in the three dimensions is achieved letting the block of variables to be [1, 2, 3]:

julia> H = HalfSpace([1.0, 1.0, 0.0], 1.0)
HalfSpace{Float64,Array{Float64,1}}([1.0, 1.0, 0.0], 1.0)

julia> project(H, [1, 2, 3])
HalfSpace{Float64,Array{Float64,1}}([1.0, 1.0, 0.0], 1.0)

Projecting along dimensions 1 and 2 only:

julia> project(H, [1, 2])
HalfSpace{Float64,Array{Float64,1}}([1.0, 1.0], 1.0)

In general, use the call syntax project(H, constrained_dimensions(H)) to return the half-space projected on the dimensions where it is constrained only:

julia> project(H, constrained_dimensions(H))
HalfSpace{Float64,Array{Float64,1}}([1.0, 1.0], 1.0)

If a constrained dimension is projected, we get the universal set of the dimension corresponding to the projection.

julia> project(H, [1, 3])
Universe{Float64}(2)

julia> project(H, [1])
Universe{Float64}(1)
source
project(P::AbstractPolyhedron{N}, block::AbstractVector{Int}) where {N}

Concrete projection of a polyhedral set.

Input

  • P – set
  • block – block structure, a vector with the dimensions of interest

Output

An HPolyhedron representing the projection of P on the dimensions specified by block.

Algorithm

  • If the unconstrained dimensions of P are a subset of the block variables, each half-sace c of P is transformed to HalfSpace(c.a[block], c.b).
  • In the general case, we compute the concrete linear map of the projection matrix associated to the given block structure.

Examples

Consider the four-dimensional cross-polytope (unit ball in the 1-norm):

julia> P = Ball1(zeros(4), 1.0);

All dimensions are constrained, and computing the (trivial) projection on the whole space behaves as expected:

julia> constrained_dimensions(P)
4-element Array{Int64,1}:
 1
 2
 3
 4

julia> P_1234 = project(P, [1, 2, 3, 4]);

julia> P_1234 == convert(HPolyhedron, P)
true

Each constraint of the cross polytope is constrained in all dimensions.

Now let's take a ball in the infinity norm and remove some constraints:

julia> B = BallInf(zeros(4), 1.0);

julia> c = constraints_list(B)[1:2]
2-element Array{HalfSpace{Float64,LazySets.Arrays.SingleEntryVector{Float64}},1}:
 HalfSpace{Float64,LazySets.Arrays.SingleEntryVector{Float64}}([1.0, 0.0, 0.0, 0.0], 1.0)
 HalfSpace{Float64,LazySets.Arrays.SingleEntryVector{Float64}}([0.0, 1.0, 0.0, 0.0], 1.0)

julia> P = HPolyhedron(c);

julia> constrained_dimensions(P)
2-element Array{Int64,1}:
 1
 2

Finally we take the concrete projection onto variables 1 and 2:

julia> project(P, [1, 2]) |> constraints_list
2-element Array{HalfSpace{Float64,VN} where VN<:AbstractArray{Float64,1},1}:
 HalfSpace{Float64,Array{Float64,1}}([1.0, 0.0], 1.0)
 HalfSpace{Float64,Array{Float64,1}}([0.0, 1.0], 1.0)
source

Convenience functions

 uniform_partition(n::Int, block_size::Int)

Compute a uniform block partition of the given size.

Input

  • n – number of dimensions of the partition
  • block_size – size of each block

Output

A vector of ranges, Vector{UnitRange{Int}}, such that the size of each block is the same, if possible.

Examples

If the number of dimensions n is 2, we have two options: either two blocks of size 1 or one block of size 2:

julia> LazySets.Approximations.uniform_partition(2, 1)
2-element Array{UnitRange{Int64},1}:
 1:1
 2:2

julia> LazySets.Approximations.uniform_partition(2, 2)
1-element Array{UnitRange{Int64},1}:
 1:2

If the block size argument is not compatible with (i.e. does not divide) n, the output is filled with one block of the size needed to reach n:

julia> LazySets.Approximations.uniform_partition(3, 1)
3-element Array{UnitRange{Int64},1}:
 1:1
 2:2
 3:3

julia> LazySets.Approximations.uniform_partition(3, 2)
2-element Array{UnitRange{Int64},1}:
 1:2
 3:3

julia> LazySets.Approximations.uniform_partition(10, 6)
2-element Array{UnitRange{Int64},1}:
 1:6
 7:10
source

Overapproximations

overapproximate(X::S, ::Type{S}, args...) where {S<:LazySet}

Overapproximating a set of type S with type S is a no-op.

Input

  • X – set
  • Type{S} – target set type
  • args – further arguments (ignored)

Output

The input set.

source
overapproximate(S::LazySet{N},
                ::Type{<:HPolygon},
                [ε]::Real=Inf) where {N<:Real}

Return an approximation of a given 2D set using iterative refinement.

Input

  • S – convex set, assumed to be two-dimensional
  • HPolygon – type for dispatch
  • ε – (optional, default: Inf) error tolerance

Output

A polygon in constraint representation.

Notes

The result is always a convex overapproximation of the input set.

If no error tolerance ε is given, or is Inf, the result is a box-shaped polygon. For convex input sets, the result is an ε-close approximation as a polygon, with respect to the Hausdorff distance.

source
overapproximate(S::LazySet, ε::Real)

Alias for overapproximate(S, HPolygon, ε).

source
overapproximate(S::LazySet, Type{<:Hyperrectangle})

Return an approximation of a given set as a hyperrectangle.

Input

  • S – set
  • Hyperrectangle – type for dispatch

Output

A hyperrectangle.

Algorithm

The center of the hyperrectangle is obtained by averaging the support function of the given set in the canonical directions, and the lengths of the sides can be recovered from the distance among support functions in the same directions.

source
overapproximate(S::CartesianProductArray{N, <:AbstractHyperrectangle{N}},
                ::Type{<:Hyperrectangle}) where {N<:Real}

Return a tight overapproximation of the Cartesian product array of a finite number of convex sets with and hyperrectangle.

Input

  • S – Cartesian product array of a finite number of convex set
  • Hyperrectangle – type for dispatch

Output

A hyperrectangle.

Algorithm

This method falls back to the corresponding convert method. Since the sets wrapped by the Cartesian product array are hyperrectangles, it can be done efficiently without overapproximation.

source
overapproximate(S::CartesianProduct{N, <:AbstractHyperrectangle{N}, <:AbstractHyperrectangle{N}},
                ::Type{<:Hyperrectangle}) where {N<:Real}

Return a tight overapproximation of the Cartesian product of two hyperrectangles by a new hyperrectangle.

Input

  • S – Cartesian product of two hyperrectangular sets
  • Hyperrectangle – type for dispatch

Output

A hyperrectangle.

Algorithm

This method falls back to the corresponding convert method. Since the sets wrapped by the Cartesian product are hyperrectangles, it can be done efficiently without overapproximation.

source
overapproximate(lm::LinearMap{N, <:AbstractHyperrectangle{N}},
                ::Type{<:Hyperrectangle}) where {N}

Return a tight overapproximation of the linear map of a hyperrectangular set using a hyperrectangle.

Input

  • S – linear map of a hyperrectangular set
  • Hyperrectangle – type for dispatch

Output

A hyperrectangle.

Algorithm

If c and r denote the center and vector radius of a hyperrectangle H, a tight hyperrectangular overapproximation of M * H is obtained by transforming c ↦ M*c and r ↦ abs.(M) * c, where abs.(⋅) denotes the element-wise absolute value operator.

source
overapproximate(r::Rectification{N}, ::Type{<:Hyperrectangle}
               ) where {N<:Real}

Overapproximate the rectification of a convex set by a tight hyperrectangle.

Input

  • r – rectification of a convex set
  • Hyperrectangle – type for dispatch

Output

A hyperrectangle.

Algorithm

Box approximation and rectification distribute. Hence we first check whether the wrapped set is empty. If so, we return the empty set. Otherwise, we compute the box approximation of the wrapped set, rectify the resulting box (which is simple), and finally convert the resulting set to a box.

source
overapproximate(S::LazySet)

Alias for overapproximate(S, Hyperrectangle).

source
overapproximate(S::LazySet{N}, ::Type{<:BallInf}) where {N<:Real}

Overapproximate a convex set by a tight ball in the infinity norm.

Input

  • S – convex set
  • BallInf – type for dispatch

Output

A tight ball in the infinity norm.

Algorithm

The center and radius of the box are obtained by evaluating the support function of the given convex set along the canonical directions.

source
overapproximate(X::ConvexHull{N, <:AbstractZonotope{N}, <:AbstractZonotope{N}},
                ::Type{<:Zonotope}) where {N<:Real}

Overapproximate the convex hull of two zonotopes.

Input

  • X – convex hull of two zonotopes
  • Zonotope – type for dispatch
  • algorithm – (optional; default: "mean") choice of algorithm; possible values are "mean" and "join"

Output

A zonotope $Z$ such that $X ⊆ Z$.

Algorithm

The algorithm can be controlled by the parameter algorithm. Note that the results of the two implemented algorithms are generally incomparable.

'mean' method

If algorithm == "mean", we choose the method proposed in [1]. The convex hull of two zonotopes $Z₁$ and $Z₂$ of the same order, which we write

\[Z_j = ⟨c^{(j)}, g^{(j)}_1, …, g^{(j)}_p⟩\]

for $j = 1, 2$, can be overapproximated as follows:

\[CH(Z_1, Z_2) ⊆ \frac{1}{2}⟨c^{(1)}+c^{(2)}, g^{(1)}_1+g^{(2)}_1, …, g^{(1)}_p+g^{(2)}_p, c^{(1)}-c^{(2)}, g^{(1)}_1-g^{(2)}_1, …, g^{(1)}_p-g^{(2)}_p⟩.\]

If the zonotope order is not the same, this algorithm calls reduce_order to reduce the order to the minimum of the arguments.

It should be noted that the output zonotope is not necessarily the minimal enclosing zonotope, which is in general expensive in high dimensions. This is further investigated in [2].

'join' method

If algorithm == "join", we choose the method proposed in [3, Definition 1]. The convex hull $X$ of two zonotopes $Z₁$ and $Z₂$ is overapproximated by a zonotope $Z₃$ such that the box approximation of $X$ is identical with the box approximation of $Z₃$. Let $□(X)$ denote the box approximation of $X$. The center of $Z₃$ is the center of $□(X)$.

The generator construction consists of two phases. In the first phase, we construct generators $g$ as a combination of one generator from $Z₁$, say, $g₁$, with another generator from $Z₂$, say, $g₂$. The entry of $g$ in the $i$-th dimension is given as

\[ g[i] = \arg\min_{\min(g₁[i], g₂[i]) ≤ x ≤ \max(g₁[i], g₂[i])} |x|.\]

If $g$ is the zero vector, it can be omitted.

In the second phase, we construct another generator for each dimension. These generators are scaled unit vectors. The following formula defines the sum of all those generators.

\[ \sup(□(X)) - c - ∑_g |g|\]

where $c$ is the center of the new zonotope and the $g$s are the generators constructed in the first phase.

References

[1] Reachability of Uncertain Linear Systems Using Zonotopes, A. Girard. HSCC 2005.

[2] Zonotopes as bounding volumes, L. J. Guibas et al, Proc. of Symposium on Discrete Algorithms, pp. 803-812.

[3] The zonotope abstract domain Taylor1+. K. Ghorbal, E. Goubault, S. Putot. CAV 2009.

source
overapproximate(lm::LinearMap{N, <:AbstractZonotope{N}},
                ::Type{<:Zonotope}) where {N<:Real}

Overapproximate a lazy linear map of a zonotopic set with a zonotope.

Input

  • lm – lazy linear map of a zonotopic set
  • Zonotope – type for dispatch

Output

The tight zonotope corresponding to lm.

source
overapproximate(Z::AbstractZonotope, ::Type{<:Hyperrectangle})

Return a tight overapproximation of a zonotope with an axis-aligned box.

Input

  • Z – zonotope
  • Hyperrectangle – type for dispatch

Output

A hyperrectangle.

Algorithm

This function implements the method in [Section 5.1.2, 1]. A zonotope $Z = ⟨c, G⟩$ can be overapproximated tightly by an axis-aligned box (i.e. a Hyperrectangle) such that its center is $c$ and the radius along dimension $i$ is the column-sum of the absolute values of the $i$-th row of $G$ for $i = 1,…, p$, where $p$ is the number of generators of $Z$.

[1] Althoff, M., Stursberg, O., & Buss, M. (2010). Computing reachable sets of hybrid systems using a combination of zonotopes and polytopes. Nonlinear analysis: hybrid systems, 4(2), 233-249.

source
overapproximate(X::LazySet{N}, dir::AbstractDirections{N}) where {N}

Overapproximating a set with template directions.

Input

  • X – set
  • dir – (concrete) direction representation

Output

A polyhedron overapproximating the set X with the directions from dir. If the directions are known to be bounded, the result is an HPolytope, otherwise the result is an HPolyhedron.

source
overapproximate(X::LazySet{N}, dir::Type{<:AbstractDirections}) where {N}

Overapproximating a set with template directions.

Input

  • X – set
  • dir – type of direction representation

Output

A polyhedron overapproximating the set X with the directions from dir. If the directions are known to be bounded, the result is an HPolytope, otherwise the result is an HPolyhedron.

source
overapproximate(S::LazySet{N}, ::Type{<:Interval}) where {N<:Real}

Return the overapproximation of a unidimensional set with an interval.

Input

  • S – one-dimensional set
  • Interval – type for dispatch

Output

An interval.

Algorithm

We use two support-function evaluations.

source
overapproximate(cap::Intersection{N}, ::Type{<:Interval}) where {N<:Real}

Return the overapproximation of a unidimensional intersection with an interval.

Input

  • cap – one-dimensional lazy intersection
  • Interval – type for dispatch

Output

An interval.

Algorithm

The algorithm recursively overapproximates the two intersected sets with intervals and then intersects these.

source
overapproximate(cap::IntersectionArray{N}, ::Type{<:Interval}) where {N<:Real}

Return the overapproximation of a unidimensional intersection with an interval.

Input

  • cap – one-dimensional lazy intersection
  • Interval – type for dispatch

Output

An interval.

Algorithm

The algorithm recursively overapproximates the two intersected sets with intervals and then intersects these.

source
overapproximate(cap::Intersection{N, <:LazySet, <:AbstractPolyhedron{N}},
                dir::AbstractDirections{N};
                kwargs...
               ) where {N<:Real}

Return the overapproximation of the intersection between a compact set and a polytope given a set of template directions.

Input

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

Output

A polytope in H-representation such that the normal direction of each half-space is given by an element of dir.

Algorithm

Let di be a direction drawn from the set of template directions dir. Let X be the compact set and let P be the polytope. We overapproximate the set X ∩ H with a polytope in constraint representation using a given set of template directions dir.

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 overapproximations can return a non-empty set even if the original intersection is empty.

source
overapproximate(cap::Intersection{N, <:HalfSpace{N}, <:AbstractPolytope{N}},
                dir::AbstractDirections{N};
                [kwargs]...
               ) where {N<:Real}

Return the overapproximation of the intersection between a half-space and a polytope given a set of template directions.

Input

  • cap – intersection of a half-space and a polytope
  • dir – template directions
  • kwargs – additional arguments that are passed to the support function algorithm

Output

A polytope in H-representation such that the normal direction of each half-space is given by an element of dir.

source
overapproximate(lm::LinearMap{N, <:CartesianProductArray{N}},
                ::Type{CartesianProductArray{N, S}}
               ) where {N, S<:LazySet{N}}

Decompose a lazy linear map of a Cartesian product array while keeping the original block structure.

Input

  • lm – lazy linear map of Cartesian product array
  • CartesianProductArray – type for dispatch

Output

A CartesianProductArray representing the decomposed linear map.

source
overapproximate(lm::LinearMap{N, <:CartesianProductArray{N}},
                ::Type{<:CartesianProductArray},
                dir::Type{<:AbstractDirections}) where {N}

Decompose a lazy linear map of a Cartesian product array with template directions while keeping the original block structure.

Input

  • lm – lazy linear map of a Cartesian product array
  • CartesianProductArray – type for dispatch
  • dir – template directions for overapproximation

Output

A CartesianProductArray representing the decomposed linear map.

source
overapproximate(lm::LinearMap{N, <:CartesianProductArray{N}},
                ::Type{<:CartesianProductArray},
                set_type::Type{<:LazySet}) where {N}

Decompose a lazy linear map of a Cartesian product array with a given set type while keeping the original block structure.

Input

  • lm – lazy linear map of a Cartesian product array
  • CartesianProductArray – type for dispatch
  • set_type – set type for overapproximation

Output

A CartesianProductArray representing the decomposed linear map.

source
overapproximate(rm::ResetMap{N, <:CartesianProductArray{N}},
                ::Type{<:CartesianProductArray}, oa) where {N}

Overapproximate a reset map (that only resets to zero) of a Cartesian product by a new Cartesian product.

Input

  • rm – reset map
  • CartesianProductArray – type for dispatch
  • oa – overapproximation option

Output

A Cartesian product with the same block structure.

Notes

This implementation currently only supports resets to zero.

Algorithm

We convert the ResetMap into a LinearMap and then call the corresponding overapproximation method.

source
overapproximate(cap::Intersection{N,
                                  <:CartesianProductArray{N},
                                  <:AbstractPolyhedron{N}},
                ::Type{CartesianProductArray}, oa) where {N}

Return the intersection of the Cartesian product of a finite number of convex sets and a polyhedron.

Input

  • cap – lazy intersection of a Cartesian product array and a polyhedron
  • CartesianProductArray – type for dispatch
  • oa – overapproximation option

Output

A CartesianProductArray that overapproximates the intersection of cpa and P.

Algorithm

The intersection only needs to be computed in the blocks of cpa that are constrained in P. Hence we first collect those constrained blocks in a lower-dimensional Cartesian product array and then convert to an HPolytope X. Then we take the intersection of X and the projection of Y onto the corresponding dimensions. (This projection is purely syntactic and exact.) Finally we decompose the result again and plug together the unaffected old blocks and the newly computed blocks. The result is a CartesianProductArray with the same block structure as in X.

source
overapproximate(Z::Zonotope{N}, ::Type{<:Zonotope}, r::Union{Integer, Rational}) where {N<:Real}

Reduce the order of a zonotope by overapproximating with a zonotope with less generators.

Input

  • Z – zonotope
  • Zonotope – desired type for dispatch
  • r – desired order

Output

A new zonotope with less generators, if possible.

Algorithm

This function implements the algorithm described in A. Girard's Reachability of Uncertain Linear Systems Using Zonotopes, HSCC. Vol. 5. 2005.

If the desired order is smaller than one, the zonotope is not reduced.

source
overapproximate(vTM::Vector{TaylorModel1{T, S}},
                ::Type{<:Zonotope}) where {T, S}

Overapproximate a taylor model in one variable with a zonotope.

Input

  • vTMTaylorModel1
  • Zonotope – type for dispatch

Output

A zonotope that overapproximates the range of the given taylor model.

Examples

If the polynomials are linear, this functions exactly transforms to a zonotope. However, the nonlinear case necessarily introduces overapproximation error. Consider the linear case first:

julia> using LazySets, TaylorModels

julia> const IA = IntervalArithmetic;

julia> I = IA.Interval(-0.5, 0.5) # interval remainder
[-0.5, 0.5]

julia> x₀ = IA.Interval(0.0) # expansion point
[0, 0]

julia> D = IA.Interval(-3.0, 1.0)
[-3, 1]

julia> p1 = Taylor1([2.0, 1.0], 2) # define a linear polynomial
 2.0 + 1.0 t + 𝒪(t³)

julia> p2 = Taylor1([0.9, 3.0], 2) # define another linear polynomial
 0.9 + 3.0 t + 𝒪(t³)

julia> vTM = [TaylorModel1(pi, I, x₀, D) for pi in [p1, p2]]
2-element Array{TaylorModel1{Float64,Float64},1}:
 2.0 + 1.0 t + [-0.5, 0.5]
 0.9 + 3.0 t + [-0.5, 0.5]

Here, vTM is a taylor model vector, since each component is a taylor model in one variable (TaylorModel1). Using overapproximate(vTM, Zonotope) we can compute its associated zonotope in generator representation:

julia> using LazySets.Approximations

julia> Z = overapproximate(vTM, Zonotope);

julia> center(Z)
2-element Array{Float64,1}:
  1.0
 -2.1

julia> Matrix(genmat(Z))
2×3 Array{Float64,2}:
 2.0  0.5  0.0
 6.0  0.0  0.5

Note how the generators of this zonotope mainly consist of two pieces: one comes from the linear part of the polynomials, and another one that corresponds to the interval remainder. This conversion gives the same upper and lower bounds as the range evaluation using interval arithmetic:

julia> X = box_approximation(Z)
Hyperrectangle{Float64}([1.0, -2.1], [2.5, 6.5])

julia> Y = evaluate(vTM[1], vTM[1].dom) × evaluate(vTM[2], vTM[2].dom)
[-1.5, 3.5] × [-8.60001, 4.40001]

julia> H = convert(Hyperrectangle, Y) # this IntevalBox is the same as X
Hyperrectangle{Float64}([1.0, -2.1], [2.5, 6.5])

However, the zonotope returns better results if we want to approximate the TM, since it is not axis-aligned:

julia> d = [-0.35, 0.93];

julia> ρ(d, Z) < ρ(d, X)
true

This function also works if the polynomials are non-linear; for example suppose that we add a third polynomial with a quadratic term:

julia> p3 = Taylor1([0.9, 3.0, 1.0], 3);

julia> vTM = [TaylorModel1(pi, I, x₀, D) for pi in [p1, p2, p3]]
3-element Array{TaylorModel1{Float64,Float64},1}:
           2.0 + 1.0 t + [-0.5, 0.5]
           0.9 + 3.0 t + [-0.5, 0.5]
  0.9 + 3.0 t + 1.0 t² + [-0.5, 0.5]

julia> Z = overapproximate(vTM, Zonotope);

julia> center(Z)
3-element Array{Float64,1}:
  1.0
 -2.1
  0.8999999999999999

julia> Matrix(genmat(Z))
3×4 Array{Float64,2}:
 2.0  0.5  0.0  0.0
 6.0  0.0  0.5  0.0
 6.0  0.0  0.0  6.5

The fourth and last generator corresponds to the addition between the interval remainder and the box overapproximation of the nonlinear part of p3 over the domain.

Algorithm

Let $\text{vTM} = (p, I)$ be a vector of $m$ taylor models, where $I$ is the interval remainder in $\mathbb{R}^m$. Let $p_{lin}$ (resp. $p_{nonlin}$) correspond to the linear (resp. nonlinear) part of each scalar polynomial.

The range of $\text{vTM}$ can be enclosed by a zonotope with center $c$ and matrix of generators $G$, $Z = ⟨c, G⟩$, by performing a conservative linearization of $\text{vTM}$:

\[ vTM' = (p', I') := (p_{lin} − p_{nonlin} , I + \text{Int}(p_{nonlin})).\]

This algorithm proceeds in two steps:

1- Conservatively linearize $\text{vTM}$ as above and compute a box overapproximation of the nonlinear part. 2- Transform the linear taylor model to a zonotope exactly through variable normalization onto the symmetric intervals $[-1, 1]$.

source
overapproximate(vTM::Vector{TaylorModelN{N, T, S}},
                ::Type{<:Zonotope}) where {N,T, S}

Overapproximate a multivariate taylor model with a zonotope.

Input

  • vTMTaylorModelN
  • Zonotope – type for dispatch

Output

A zonotope that overapproximates the range of the given taylor model.

Examples

Consider a vector of two 2-dimensional taylor models of order 2 and 4 respectively.

julia> using LazySets, LazySets.Approximations, TaylorModels

julia> const IA = IntervalArithmetic;

julia> x₁, x₂ = set_variables(Float64, ["x₁", "x₂"], order=8)
2-element Array{TaylorN{Float64},1}:
  1.0 x₁ + 𝒪(‖x‖⁹)
  1.0 x₂ + 𝒪(‖x‖⁹)

julia> x₀ = IntervalBox(0..0, 2) # expansion point
[0, 0] × [0, 0]

julia> Dx₁ = IA.Interval(0.0, 3.0) # domain for x₁
[0, 3]

julia> Dx₂ = IA.Interval(-1.0, 1.0) # domain for x₂
[-1, 1]

julia> D = Dx₁ × Dx₂ # take the Cartesian product of the domain on each variable
[0, 3] × [-1, 1]

julia> r = IA.Interval(-0.5, 0.5) # interval remainder
[-0.5, 0.5]

julia> p1 = 1 + x₁^2 - x₂
 1.0 - 1.0 x₂ + 1.0 x₁² + 𝒪(‖x‖⁹)

julia> p2 = x₂^3 + 3x₁^4 + x₁ + 1
 1.0 + 1.0 x₁ + 1.0 x₂³ + 3.0 x₁⁴ + 𝒪(‖x‖⁹)

julia> vTM = [TaylorModelN(pi, r, x₀, D) for pi in [p1, p2]]
2-element Array{TaylorModelN{2,Float64,Float64},1}:
             1.0 - 1.0 x₂ + 1.0 x₁² + [-0.5, 0.5]
   1.0 + 1.0 x₁ + 1.0 x₂³ + 3.0 x₁⁴ + [-0.5, 0.5]

julia> Z = overapproximate(vTM, Zonotope);

julia> center(Z)
2-element Array{Float64,1}:
   5.5
 124.0

julia> Matrix(genmat(Z))
2×4 Array{Float64,2}:
 0.0  -1.0  5.0    0.0
 1.5   0.0  0.0  123.0

Algorithm

We refer to the algorithm description for the univariate case.

source

Underapproximations

underapproximate(X::LazySet{N}, dirs::AbstractDirections;
                [apply_convex_hull]::Bool=false) where {N<:Real}

Compute the underapproximation of a convex set by sampling support vectors.

Input

  • X – set
  • dirs – directions
  • apply_convex_hull – (optional, default: false) if true, post-process the support vectors with a convex hull operation

Output

The VPolytope obtained by taking the convex hull of the support vectors of X along the directions determined by dirs.

source

Box Approximations

ballinf_approximation(S::LazySet{N}) where {N<:Real}

Overapproximate a convex set by a tight ball in the infinity norm.

Input

  • S – convex set

Output

A tight ball in the infinity norm.

source
box_approximation(S::LazySet{N}) where {N<:Real}

Overapproximate a convex set by a tight hyperrectangle.

Input

  • S – convex set

Output

A tight hyperrectangle.

source
interval_hull

Alias for box_approximation.

source
box_approximation_symmetric(S::LazySet{N}) where {N<:Real}

Overapproximate a convex set by a tight hyperrectangle centered in the origin.

Input

  • S – convex set

Output

A tight hyperrectangle centered in the origin.

Algorithm

The center of the box is the origin, and the radius is obtained by computing the maximum value of the support function evaluated at the canonical directions.

source
symmetric_interval_hull

Alias for box_approximation_symmetric.

source
box_approximation_helper(S::LazySet{N}) where {N<:Real}

Common code of box_approximation and box_approximation_symmetric.

Input

  • S – convex set

Output

A tuple containing the data that is needed to construct a tightly overapproximating hyperrectangle.

  • c – center
  • r – radius

Algorithm

The center of the hyperrectangle is obtained by averaging the support function of the given convex set in the canonical directions. The lengths of the sides can be recovered from the distance among support functions in the same directions.

source

Iterative refinement

LocalApproximation{N<:Real, VN<:AbstractVector{N}}

Type that represents a local approximation in 2D.

Fields

  • p1 – first inner point
  • d1 – first direction
  • p2 – second inner point
  • d2 – second direction
  • q – intersection of the lines l1 ⟂ d1 at p1 and l2 ⟂ d2 at p2
  • refinable – states if this approximation is refinable
  • err – error upper bound

Notes

The criteria for being refinable are determined in the method new_approx.

source
PolygonalOverapproximation{N<:Real, SN<:LazySet{N}, VN<:AbstractVector{N}}

Type that represents the polygonal approximation of a convex set.

Fields

  • S – convex set
  • approx_stack – stack of local approximations that still need to be examined
  • constraints – vector of linear constraints that are already finalized (i.e., they satisfy the given error bound)
source
new_approx(S::LazySet, p1::VN, d1::VN,
           p2::VN, d2::VN) where {N<:AbstractFloat, VN<:AbstractVector{N}}

Create a LocalApproximation instance for the given excerpt of a polygonal approximation.

Input

  • S – convex set
  • p1 – first inner point
  • d1 – first direction
  • p2 – second inner point
  • d2 – second direction

Output

A local approximation of S in the given directions.

source
addapproximation!(Ω::PolygonalOverapproximation, p1::VN, d1::VN,
                  p2::VN, d2::VN) where {N<:Real, VN<:AbstractVector{N}}

Input

  • Ω – polygonal overapproximation of a convex set
  • p1 – first inner point
  • d1 – first direction
  • p2 – second inner point
  • d2 – second direction

Output

The list of local approximations in Ω of the set Ω.S is updated in-place and the new approximation is returned by this function.

source
refine(approx::LocalApproximation, S::LazySet)

Refine a given local approximation of the polygonal approximation of a convex set by splitting along the normal direction of the approximation.

Input

  • approx – local approximation to be refined
  • S – 2D convex set

Output

The tuple consisting of the refined right and left local approximations.

source
tohrep(Ω::PolygonalOverapproximation)

Convert a polygonal overapproximation into a concrete polygon.

Input

  • Ω – polygonal overapproximation of a convex set

Output

A polygon in constraint representation.

Algorithm

Internally we keep the constraints sorted. Hence we do not need to use addconstraint! when creating the HPolygon.

source
approximate(S::LazySet{N}, ε::N) where {N<:AbstractFloat}

Return an ε-close approximation of the given 2D convex set (in terms of Hausdorff distance) as an inner and an outer approximation composed by sorted local Approximation2D.

Input

  • S – 2D convex set
  • ε – error bound

Output

An ε-close approximation of the given 2D convex set.

source
constraint(approx::LocalApproximation)

Convert a local approximation to a linear constraint.

Input

  • approx – local approximation

Output

A linear constraint.

source

Template directions

AbstractDirections{N}

Abstract type for template direction representations.

Notes

All subtypes should implement the standard iterator methods from Base and the function dim(d<:AbstractDirections).

source
isbounding(ad::AbstractDirections)

Checks if an overapproximation with a list of template directions results in a bounded set, given a bounded input set.

Input

  • ad – template directions

Output

Given a bounded set $X$, we can construct an outer approximation of $X$ by using the template directions ad as normal vectors of the facets. If this function returns true, then the result is again a bounded set (i.e., a polytope). Note that the result does not depend on the specific shape of $X$, as long as $X$ is bounded.

Notes

By default, this function returns false in order to be conservative. Custom subtypes of AbstractDirections should hence add a method for this function.

source
BoxDirections{N} <: AbstractDirections{N}

Box direction representation.

Fields

  • n – dimension
source
OctDirections{N} <: AbstractDirections{N}

Octagon direction representation.

Fields

  • n – dimension

Notes

Octagon directions consist of all vectors that are zero almost everywhere except in two dimensions $i$, $j$ (possibly $i = j$) where it is $±1$.

source
BoxDiagDirections{N} <: AbstractDirections{N}

Box-diagonal direction representation.

Fields

  • n – dimension

Notes

Box-diagonal directions can be seen as the union of diagonal directions (all entries are ±1) and box directions (one entry is ±1, all other entries are 0). The iterator first enumerates all diagonal directions, and then all box directions.

source
PolarDirections{N<:AbstractFloat} <: AbstractDirections{N}

Polar directions representation.

Fields

  • – length of the partition of the polar angle

Notes

The PolarDirections constructor provides a sample of the unit sphere in $\mathbb{R}^2$, which is parameterized by the polar angles $φ ∈ Dφ := [0, 2π]$ respectively; see the wikipedia entry Polar coordinate system. The domain $Dφ$ is discretized in $Nφ$ pieces. Then the Cartesian components of each direction are obtained with

\[[cos(φᵢ), sin(φᵢ)].\]

Examples

The integer passed as an argument is used to discretize $φ$:

julia> using LazySets.Approximations: PolarDirections

julia> pd = PolarDirections(2)
PolarDirections{Float64}(2, Array{Float64,1}[[1.0, 0.0], [-1.0, 1.22465e-16]])

julia> pd.Nφ
2
source
SphericalDirections{N<:AbstractFloat} <: AbstractDirections{N}

Spherical directions representation.

Fields

  • – length of the partition of the azimuthal angle
  • – length of the partition of the polar angle
  • stack – list of computed directions

Notes

The SphericalDirections constructor provides a sample of the unit sphere in $\mathbb{R}^3$, which is parameterized by the azimuthal and polar angles $θ ∈ Dθ := [0, π]$ and $φ ∈ Dφ := [0, 2π]$ respectively, see the wikipedia entry Spherical coordinate system. The domains $Dθ$ and $Dφ$ are discretized in $Nθ$ and $Nφ$ respectively. Then the Cartesian components of each direction are obtained with

\[[sin(θᵢ)*cos(φᵢ), sin(θᵢ)*sin(φᵢ), cos(θᵢ)].\]

The north and south poles are treated separately so that those points are not considered more than once.

Examples

A SphericalDirections can be built in different ways. If you pass only one integer, it is used to discretize both $θ$ and $φ$:

julia> using LazySets.Approximations: SphericalDirections

julia> sd = SphericalDirections(3)
SphericalDirections{Float64}(3, 3, Array{Float64,1}[[0.0, 0.0, 1.0], [0.0, 0.0, -1.0], [1.0, 0.0, 6.12323e-17], [-1.0, 1.22465e-16, 6.12323e-17]])

julia> sd.Nθ, sd.Nφ 
(3, 3)

Pass two integers to control the discretization in $θ$ and in $φ$ separately:

julia> sd_4_5 = SphericalDirections(4, 5);

julia> length(sd_4_5)
10

julia> sd_4_8 = SphericalDirections(4, 8);

julia> length(sd_4_8)
16
source
CustomDirections{N, VN<:AbstractVector{N}} <: AbstractDirections{N}

User-defined template directions.

Fields

  • directions – list of template directions
  • n – dimension
  • isbounding – boundedness status

Notes

Since custom directions may not be bounded, boundedness can either be asserted in the constructor or it will be determined automatically by performing an experiment (more precisely, we overapproximate the unit ball in the infinity norm using the given directions).

The dimension will also be determined automatically, unless the empty vector is passed (in which case the optional argument n needs to be specified).

source

See also overapproximate(X::LazySet, dir::AbstractDirections)::HPolytope.

Hausdorff distance

hausdorff_distance(X::LazySet{N}, Y::LazySet{N}; [p]::N=N(Inf),
                   [ε]::N=N(1e-3)) where {N<:Real}

Compute the Hausdorff distance between two convex sets up to a given threshold.

Input

  • X – convex set
  • Y – convex set
  • p – (optional, default: Inf) norm parameter of the Hausdorff distance
  • ε – (optional, default: 1e-3) precision threshold; the true Hausdorff distance may diverge from the result by at most this value

Output

A value from the $ε$-neighborhood of the Hausdorff distance between $X$ and $Y$.

Notes

Given a $p$-norm, the Hausdorff distance $d_H^p(X, Y)$ between sets $X$ and $Y$ is defined as follows:

\[ d_H^p(X, Y) = \inf\{δ ≥ 0 \mid Y ⊆ X ⊕ δ 𝐵_p^n \text{ and } X ⊆ Y ⊕ δ 𝐵_p^n\}\]

Here $𝐵_p^n$ is the $n$-dimensional unit ball in the $p$-norm.

The implementation may internally rely on the support function of $X$ and $Y$; hence any imprecision in the implementation of the support function may affect the result. At the time of writing, the only set type with imprecise support function is the lazy Intersection.

Algorithm

We perform binary search for bounding the Hausdorff distance in an interval $[l, u]$, where initially $l$ is $0$ and $u$ is described below. The binary search terminates when $u - l ≤ ε$, i.e., the interval becomes sufficiently small.

To find an upper bound $u$, we start with the heuristics of taking the biggest distance in the axis-parallel directions. As long as this bound does not work, we increase the bound by $2$.

Given a value $δ$, to check whether the sets are within Hausdorff distance $δ$, we simply check the inclusions given above, where on the right-hand side we use a lazy MinkowskiSum with a Ballp centered in the origin.

source