Overapproximation

LazySets.Approximations.overapproximateFunction
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)
  • kwargs – further keyword arguments (ignored)

Output

The input set.

source
overapproximate(S::LazySet, T::Type{<:LazySet}, [args]...; [kwargs]...)

Default overapproximation method that falls back to convert.

Input

  • X – set
  • Type{S} – target set type
  • args – further arguments
  • kwargs – further keyword arguments

Output

The result of convert, or a MethodError if no such method exists.

source
overapproximate(S::LazySet)

Alias for overapproximate(S, Hyperrectangle) resp. box_approximation(S).

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

Alias for box_approximation(S).

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

Alias for ballinf_approximation(S).

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

Overapproximate a given 2D set using iterative refinement.

Input

  • S – two-dimensional bounded set
  • HPolygon – target set type
  • ε – (optional, default: Inf) error tolerance
  • prune – (optional, default: true) flag for removing redundant constraints in the end

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 polygonal overapproximation with respect to the Hausdorff distance.

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

Alias for overapproximate(S, HPolygon, ε).

source
overapproximate(X::LazySet{N}, dirs::AbstractDirections;
                [prune]::Bool=true) where {N}

Overapproximate a (possibly unbounded) set with template directions.

Input

  • X – set
  • dirs – directions
  • prune – (optional, default: true) flag for removing redundant constraints

Output

A polyhedron overapproximating the set X with the directions from dirs. The overapproximation is computed using the support function. The result is an HPolytope if it is bounded and otherwise an HPolyhedron.

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

Overapproximate a set with template directions.

Input

  • X – set
  • dirs – type of direction representation

Output

A polyhedron overapproximating the set X with the directions from dirs. The result is an HPolytope if it is bounded and otherwise an HPolyhedron.

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

Overapproximate the intersection between a set and a polyhedron given a set of template directions.

Input

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

Output

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

Algorithm

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

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

This algorithm is inspired from [1].

Notes

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

This method may return a non-empty set even if the original set is empty.

[1] G. Frehse, R. Ray. Flowpipe-Guard Intersection for Reachability Computations with Support Functions. ADHS 2012.

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

Overapproximate 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
  • dirs – 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 dirs.

source
overapproximate(Z::AbstractZonotope, ::Type{<:HParallelotope},
                [indices]=1:dim(Z))

Overapproximate a zonotopic set with a parallelotope in constraint representation.

Input

  • Z – zonotopic set
  • HParallelotope – target set type
  • indices – (optional; default: 1:dim(Z)) generator indices selected when constructing the parallelotope

Output

An overapproximation of the given zonotopic set using a parallelotope.

Algorithm

The algorithm is based on Proposition 8 discussed in Section 5 of [1].

[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::Intersection{N, <:AbstractZonotope, <:Hyperplane},
                dirs::AbstractDirections) where {N}

Overapproximate the intersection between a zonotopic set and a hyperplane with a polyhedron or polytope using the given directions.

Input

  • X – intersection between a zonotopic set and a hyperplane
  • dirs – type of direction representation

Output

An overapproximation of the intersection between a zonotopic set and a hyperplane. If the directions are bounding, the result is an HPolytope, otherwise the result is an HPolyhedron.

Algorithm

This function implements [Algorithm 8.1, 1].

[1] Colas Le Guernic. Reachability Analysis of Hybrid Systems with Linear continuous dynamics (Doctoral dissertation). 2009.

source
overapproximate(QM::QuadraticMap{N, <:SparsePolynomialZonotope},
                ::Type{<:SparsePolynomialZonotope}) where {N}

Overapproximate a quadratic map of a sparse polynomial zonotope with a sparse polynomial zonotope.

Input

  • QM – quadratic map of a sparse polynomial zonotope
  • SparsePolynomialZonotope – target type

Output

A sparse polynomial zonotope overapproximating the quadratic map of a sparse polynomial zonotope.

Algorithm

This method implements Proposition 13 of [1].

[1] N. Kochdumper and M. Althoff. Sparse Polynomial Zonotopes: A Novel Set Representation for Reachability Analysis. Transactions on Automatic Control

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

Return the overapproximation of a set with an interval.

Input

  • S – one-dimensional set
  • Interval – target type

Output

An interval.

Algorithm

We use the extrema function.

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

Return the overapproximation of a lazy intersection with an interval.

Input

  • cap – one-dimensional intersection
  • Interval – target type

Output

An interval.

Algorithm

The algorithm recursively overapproximates the two intersected sets with intervals and then intersects these. (Note that this fails if the sets are unbounded.)

For convex sets this algorithm is precise.

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

Return the overapproximation of an intersection array with an interval.

Input

  • cap – one-dimensional intersection array
  • Interval – target type

Output

An interval.

Algorithm

The algorithm recursively overapproximates the intersected sets with intervals and then intersects these. (Note that this fails if the sets are unbounded.)

For convex sets this algorithm is precise.

source
overapproximate(Z::AbstractZonotope, ::Type{<:Zonotope}, r::Real)

Reduce the order of a zonotopic set.

Input

  • Z – zonotopic set
  • Zonotope – target set type
  • r – desired order

Output

A new zonotope with r generators, if possible.

Algorithm

This method falls back to reduce_order with the default algorithm.

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

Overapproximate the convex hull of two zonotopic sets.

Input

  • X – convex hull of two zonotopic sets
  • Zonotope – target set type
  • 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 zonotopic sets $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 to obtain 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 zonotopic sets $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] = \argmin_{\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, A. T. Nguyen, L. Zhang. SODA 2003.

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

source
overapproximate(P::SimpleSparsePolynomialZonotope, ::Type{<:Zonotope})

Overapproximate a simple sparse polynomial zonotope with a zonotope.

Input

  • P – simple sparse polynomial zonotope
  • Zonotope – target set type

Output

A zonotope.

source
overapproximate(P::SimpleSparsePolynomialZonotope, ::Type{<:Zonotope},
                dom::IntervalBox)

Overapproximate a simple sparse polynomial zonotope over the parameter domain dom with a zonotope.

Input

  • P – simple sparse polynomial zonotope
  • Zonotope – target set type
  • dom – parameter domain, which should be a subset of [-1, 1]^q, where q = nparams(P)

Output

A zonotope.

source
overapproximate(P::SparsePolynomialZonotope, ::Type{<:Zonotope})

Overapproximate a sparse polynomial zonotope with a zonotope.

Input

  • P – sparse polynomial zonotope
  • Zonotope – target set type

Output

A zonotope.

source
overapproximate(P::DensePolynomialZonotope, ::Type{<:Zonotope})

Overapproximate a polynomial zonotope with a zonotope.

Input

  • P – polynomial zonotope
  • Zonotope – target set type

Output

A zonotope.

Algorithm

This method implements Proposition 1 in [1].

[1] M. Althoff in Reachability analysis of nonlinear systems using conservative polynomialization and non-convex sets, Hybrid Systems: Computation and Control, 2013, pp. 173-182.

source
overapproximate(X::LazySet, ZT::Type{<:Zonotope},
                dir::Union{AbstractDirections, Type{<:AbstractDirections}};
                [algorithm]="vrep", kwargs...)

Overapproximate a set with a zonotope.

Input

  • X – set
  • Zonotope – target set type
  • dir – directions used for the generators
  • algorithm – (optional, default: "vrep") algorithm used to compute the overapproximation
  • kwargs – further algorithm-specific options

Output

A zonotope that overapproximates X and uses at most the generator directions provided in dir (redundant directions will be ignored).

Notes

Two algorithms are available:

  • "vrep" – Overapproximate a polytopic set with a zonotope of minimal total generator sum using only generators in the given directions. Under this constraint, the zonotope has the minimal sum of generator vectors. See the docstring of _overapproximate_zonotope_vrep for further details.

  • "cpa" – Overapproximate a polytopic set with a zonotope using a Cartesian decomposition into two-dimensional blocks. See the docstring of _overapproximate_zonotope_cpa for further details.

source
overapproximate(r::Rectification{N, <:AbstractZonotope},
                ::Type{<:Zonotope}) where {N}

Overapproximate the rectification of a zonotopic set with a zonotope.

Input

  • r – lazy rectification of a zonotopic set
  • Zonotope – target set type

Output

A zonotope overapproximation of the set obtained by rectifying Z.

Algorithm

This method implements [Theorem 3.1, 1].

[1] Singh, G., Gehr, T., Mirman, M., Püschel, M., & Vechev, M. Fast and effective robustness certification. NeurIPS 2018.

source
overapproximate(CHA::ConvexHullArray{N, <:AbstractZonotope},
                ::Type{<:Zonotope}) where {N}

Overapproximate the convex hull of a list of zonotopic sets with a zonotope.

Input

  • CHA – convex hull array of zonotopic sets
  • Zonotope – target set type

Output

A zonotope overapproximation of the convex hull array of zonotopic sets.

Algorithm

This method iteratively applies the overapproximation algorithm to the convex hull of two zonotopic sets from the given array of zonotopic sets.

source
overapproximate(QM::QuadraticMap{N, <:AbstractZonotope},
                ::Type{<:Zonotope}) where {N}

Overapproximate a quadratic map of a zonotopic set with a zonotope.

Input

  • QM – quadratic map of a zonotopic set
  • Zonotope – target set type

Output

A zonotope overapproximating the quadratic map of a zonotopic set.

Notes

Mathematically, a quadratic map of a zonotope with matrices $Q$ is defined as:

\[ Z_Q = \right\{ λ \mid λ_i = x^T Q\^{(i)} x,~i = 1, …, n,~x ∈ Z \left\}\]

Algorithm

This method implements [Lemma 1, 1].

[1] Matthias Althoff and Bruce H. Krogh. Avoiding geometric intersection operations in reachability analysis of hybrid systems. HSCC 2012.

source
overapproximate(X::Intersection{N, <:AbstractZonotope, <:Hyperplane},
                ::Type{<:Zonotope})

Overapproximate the intersection of a zonotopic set and a hyperplane with a zonotope.

Input

  • X – intersection of a zonotopic set and a hyperplane
  • Zonotope – target set type

Output

A zonotope overapproximating the intersection.

Algorithm

This method implements Algorithm 3 in [1].

[1] Moussa Maïga, Nacim Ramdani, Louise Travé-Massuyès, Christophe Combastel: A CSP versus a zonotope-based method for solving guard set intersection in nonlinear hybrid reachability. Mathematics in Computer Science (8) 2014.

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

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 – target set type

Output

A CartesianProductArray representing the decomposed linear map.

source
overapproximate(lm::LinearMap{N, <:CartesianProductArray},
                ::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 – target set type
  • dir – template directions for overapproximation

Output

A CartesianProductArray representing the decomposed linear map.

source
overapproximate(lm::LinearMap{N, <:CartesianProductArray},
                ::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 – target set type
  • set_type – set type for overapproximation

Output

A CartesianProductArray representing the decomposed linear map.

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

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

Input

  • rm – reset map
  • CartesianProductArray – target set type
  • 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 overapproximate method.

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

Overapproximate the intersection of a Cartesian product of a finite number of sets and a polyhedron with a new Cartesian product.

Input

  • cap – lazy intersection of a Cartesian product array and a polyhedron
  • CartesianProductArray – target set type
  • 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(vTM::Vector{TaylorModel1{T, S}}, ::Type{<:Zonotope};
                [remove_redundant_generators]::Bool=true
                [normalize]::Bool=true) where {T, S}

Overapproximate a Taylor model in one variable with a zonotope.

Input

  • vTM – vector of TaylorModel1
  • Zonotope – target set type
  • remove_redundant_generators – (optional; default: true) flag to remove redundant generators of the resulting zonotope
  • normalize – (optional; default: true) flag to skip the normalization of the Taylor models

Output

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

Examples

If the polynomials are linear, this method exactly transforms to a zonotope. 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 Vector{TaylorModel1{Float64, Float64}}:
  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> Z = overapproximate(vTM, Zonotope);

julia> center(Z)
2-element Vector{Float64}:
  1.0
 -2.1

julia> Matrix(genmat(Z))
2×3 Matrix{Float64}:
 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 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, Vector{Float64}, Vector{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 IntervalBox is the same as X
Hyperrectangle{Float64, StaticArraysCore.SVector{2, Float64}, StaticArraysCore.SVector{2, Float64}}([1.0, -2.1000000000000005], [2.5, 6.500000000000001])

However, the zonotope returns better results if we want to approximate the Taylor model because it is not axis-aligned:

julia> d = [-0.35, 0.93];

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

This method 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)
 0.9 + 3.0 t + 1.0 t² + 𝒪(t⁴)

julia> vTM = [TaylorModel1(pi, I, x₀, D) for pi in [p1, p2, p3]]
3-element Vector{TaylorModel1{Float64, Float64}}:
           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 Vector{Float64}:
  1.0
 -2.1
  2.4

julia> Matrix(genmat(Z))
3×4 Matrix{Float64}:
 2.0  0.5  0.0  0.0
 6.0  0.0  0.5  0.0
 6.0  0.0  0.0  5.0

The last generator corresponds to the addition of 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 $ℝ^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};
                [remove_redundant_generators]::Bool=true
                [normalize]::Bool=true) where {N, T, S}

Overapproximate a multivariate Taylor model with a zonotope.

Input

  • vTM – vector of TaylorModelN
  • Zonotope – target set type
  • remove_redundant_generators – (optional; default: true) flag to remove redundant generators of the resulting zonotope
  • normalize – (optional; default: true) flag to skip the normalization of the Taylor models

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, TaylorModels

julia> const IA = IntervalArithmetic;

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

julia> x₀ = IA.IntervalBox(0..0, 2) # expansion point
[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 Vector{TaylorModelN{2, Float64, Float64}}:
            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 Vector{Float64}:
   5.5
 124.0

julia> Matrix(genmat(Z))
2×2 Matrix{Float64}:
   0.0  -6.0
 124.5   0.0

Algorithm

We refer to the algorithm description for the univariate case.

source
LazySets.Approximations._overapproximate_zonotope_vrepFunction
_overapproximate_zonotope_vrep(X::LazySet{N},
                               dir::AbstractDirections;
                               solver=default_lp_solver(N)) where {N}

Overapproximate a polytopic set with a zonotope of minimal total generator sum using only generators in the given directions.

Input

  • X – polytopic set
  • dir – directions used for the generators
  • solver – (optional, default: default_lp_solver(N)) the backend used to solve the linear program

Output

A zonotope that overapproximates X and uses at most the directions provided in dir (redundant directions will be ignored). Under this constraint, the zonotope has the minimal sum of generator vectors.

Notes

The algorithm only requires one representative of each generator direction and their additive inverse (e.g., only one of [1, 0] and [-1, 0]) and assumes that the directions are normalized. We preprocess the directions in that respect.

Algorithm

We solve a linear program parametric in the vertices $v_j$ of X and the directions $d_k$ in dir presented in Section 4.2 in [1], adapting the notation to the one used in this library.

\[ \min ∑_{k=1}^l α_k \ s.t. \ c + ∑_{k=1}^l b_{kj} * d_k = v_j \quad ∀ j \ -α_k ≤ b_{kj} ≤ α_k \quad ∀ k, j \ α_k ≥ 0 \quad ∀ k\]

The resulting zonotope has center c and generators α_k · d_k.

Note that the first type of side constraints is vector-based and that the nonnegativity constraints (last type) are not stated explicitly in [1].

[1] Zonotopes as bounding volumes. L. J. Guibas, A. T. Nguyen, L. Zhang. SODA 2003.

source
LazySets.Approximations._overapproximate_zonotope_cpaFunction
_overapproximate_zonotope_cpa(X::LazySet, dir::Type{<:AbstractDirections})

Overapproximate a polytopic set with a zonotope using Cartesian decomposition.

Input

  • X – polytopic set
  • dir – directions used for the generators

Output

A zonotope that overapproximates X.

Notes

The algorithm decomposes X into 2D sets and overapproximates those sets with zonotopes, and finally converts the Cartesian product of the sets to a zonotope.

Algorithm

The implementation is based on Section 8.2.4 in [1].

[1] Le Guernic, C. Reachability analysis of hybrid systems with linear continuous dynamics (Doctoral dissertation). 2009.

source