Intersection

Binary intersection (Intersection)

LazySets.IntersectionType
Intersection{N, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents the intersection of two sets.

Fields

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

Notes

If the arguments of the lazy intersection are half-spaces, the set is simplified to a polyhedron in constraint representation (HPolyhedron).

The intersection preserves convexity: if the set arguments are convex, then their intersection is convex as well.

Examples

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

julia> X, Y = BallInf([0.0, 0.0], 0.5), BallInf([1.0, 0.0], 0.75);

julia> Z = X ∩ Y;

julia> typeof(Z)
Intersection{Float64, BallInf{Float64, Vector{Float64}}, BallInf{Float64, Vector{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, Vector{Float64}, Vector{Float64}}([0.375, 0.0], [0.125, 0.5])
source
Base.:∩Method
∩(X::LazySet, Y::LazySet)

Alias for the lazy intersection.

Notes

The function symbol can be typed via \cap[TAB].

source
LazySets.dimMethod
dim(cap::Intersection)

Return the dimension of an intersection of two sets.

Input

  • cap – intersection of two sets

Output

The ambient dimension of the intersection of two sets.

source
LazySets.ρMethod
ρ(d::AbstractVector, cap::Intersection)

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

Input

  • d – direction
  • cap – intersection of two sets

Output

An upper 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-function evaluations for $X$ and $Y$.

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

Evaluate 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/line; evaluate 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 first set of the intersection (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:

[1] G. Frehse, R. Ray. Flowpipe-Guard Intersection for Reachability Computations with Support Functions. [2] C. Le Guernic. Reachability Analysis of Hybrid Systems with Linear Continuous Dynamics, PhD thesis [3] T. Rockafellar, R. Wets. Variational Analysis.

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

Return an upper bound on the support function of the intersection between a compact set and a polyhedron along a given direction.

Input

  • d – direction
  • cap – intersection of a compact set and a polyhedron
  • 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 polyhedron and then take the minimum. This gives an overapproximation of the exact support value.

This algorithm is inspired from [1].

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

Notes

This method relies on the constraints_list of the polyhedron.

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

Evaluate the support function of the intersection between two polyhedral sets.

Input

  • d – direction
  • cap – intersection of two polyhedral sets
  • kwargs – additional arguments that are passed to the support-function algorithm

Output

The evaluation of the support function in the given direction.

Algorithm

We combine the constraints of the two polyhedra to a new HPolyhedron, for which we then evaluate the support function.

source
LazySets.σMethod
σ(d::AbstractVector, cap::Intersection)

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

Input

  • d – direction
  • cap – intersection of two sets

Output

A support vector in the given direction.

Algorithm

We compute the concrete intersection, which may be expensive.

source
LazySets.isboundedMethod
isbounded(cap::Intersection)

Check whether an intersection of two sets is bounded.

Input

  • cap – intersection of two 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 LazySets._isbounded_unit_dimensions.

source
Base.isemptyMethod
isempty(cap::Intersection)

Check whether the intersection of two sets is empty.

Input

  • cap – intersection of two 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, cap::Intersection)

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

Input

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

Output

true iff $x ∈ cap$.

Algorithm

A point $x$ is in the intersection iff it is in each set.

source
LazySets.constraints_listMethod
constraints_list(cap::Intersection)

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

Input

  • cap – intersection of two (polyhedral) sets

Output

A 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
LazySets.vertices_listMethod
vertices_list(cap::Intersection)

Return a list of vertices of a lazy intersection of two (polyhedral) sets.

Input

  • cap – intersection of two (polyhedral) sets

Output

A list containing the vertices of the lazy intersection of two sets.

Notes

We assume that the underlying sets are polyhedral and that the intersection is bounded.

Algorithm

We compute the concrete intersection using intersection and then take the vertices of that representation.

source
LazySets.isempty_knownMethod
isempty_known(cap::Intersection)

Ask whether the status of emptiness is known.

Input

  • cap – intersection of two sets

Output

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

source
LazySets.set_isempty!Method
set_isempty!(cap::Intersection, isempty::Bool)

Set the status of emptiness in the cache.

Input

  • cap – intersection of two sets
  • isempty – new status of emptiness
source
LazySets.swapMethod
swap(cap::Intersection)

Return a new Intersection object with the arguments swapped.

Input

  • cap – intersection of two 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
LazySets.use_precise_ρFunction
use_precise_ρ(cap::Intersection)

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

Input

  • cap – intersection of two 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::Union{<:HalfSpace, <:Hyperplane, <:Line2D}; [kwargs...])

Given a convex set $X$ and a half-space $H = \{x: a^T x ≤ b \}$ or a hyperplane/line $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 – convex set
  • H – half-space or hyperplane or line

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 >= 1

julia> using Optim

julia> using LazySets: _line_search

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

julia> v[1]
1.0

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

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

julia> v[1]
1.0

Instead of Brent's method we can use the Golden Section method:

julia> v = _line_search([1.0, 0.0], X, H, upper=1e3, method=GoldenSection());

julia> v[1]
1.0
source
LazySets._projectionFunction
_projection(ℓ, X, H::Union{Hyperplane{N}, Line2D{N}};
            [lazy_linear_map]=false,
            [lazy_2d_intersection]=true,
            [algorithm_2d_intersection]=nothing,
            [kwargs...]) where {N}

Given a 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 – convex set
  • H – hyperplane
  • lazy_linear_map – (optional, default: false) flag to perform the projection lazily or concretely
  • lazy_2d_intersection – (optional, default: true) flag 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 used

Output

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

Algorithm

This projection method is based on Prop. 8.2, [1, page 103].

In the original algorithm, Section 8.2 of [1], 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 [1]), if one uses dispatch on ρ(y_dir, Sℓ⋂Lγ; kwargs...) given that Sℓ is itself a zonotope.

Notes

This function depends on the calculation of the support function of another set in two dimensions. Obviously one does not want to use algorithm="projection" again for this second calculation. The option algorithm_2d_intersection is used for that: if 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 arguments.

[1] C. Le Guernic. Reachability Analysis of Hybrid Systems with Linear Continuous Dynamics, PhD thesis.

source
LazySets.linear_mapMethod
linear_map(M::AbstractMatrix, cap::Intersection)

Return the concrete linear map of an intersection of two sets.

Input

  • M – matrix
  • cap – intersection of two sets

Output

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

Algorithm

This method computes the concrete intersection.

source
LazySets.plot_recipeMethod
plot_recipe(cap::Intersection{N}, [ε]::N=-one(N),
            [Nφ]::Int=PLOT_POLAR_DIRECTIONS) where {N}

Convert an intersection of two sets to a pair (x, y) of points for plotting.

Input

  • cap – intersection of two sets
  • ε – (optional, default 0) ignored, used for dispatch
  • – (optional, default: PLOT_POLAR_DIRECTIONS) number of polar directions used in the template overapproximation

Output

A pair (x, y) of points that can be plotted.

source
RecipesBase.apply_recipeMethod
plot_intersection(cap::Intersection{N}, [ε]::Real=zero(N),
                  [Nφ]::Int=PLOT_POLAR_DIRECTIONS) where {N}

Plot a lazy intersection.

Input

  • cap – lazy intersection
  • ε – (optional, default 0) ignored, used for dispatch
  • – (optional, default: PLOT_POLAR_DIRECTIONS) number of polar directions used in the template overapproximation

Notes

This function is separated from the main LazySet plot recipe because iterative refinement is not available for lazy intersections (since it uses the support vector (but see #1187)).

Also note that if the set is a nested intersection, you may have to manually overapproximate this set before plotting (see overapproximate for details).

Examples

julia> X = Ball2(zeros(2), 1.) ∩ Ball2(ones(2), 1.5);  # lazy intersection

julia> plot(X)

You can specify the accuracy of the overapproximation of the lazy intersection by passing an explicit value for , which stands for the number of polar directions used in the overapproximation. This number can also be passed to the plot function directly.

julia> plot(overapproximate(X, PolarDirections(100)))

julia> plot(X, 0.0, 100)  # equivalent to the above line
source

Inherited from LazySet:

Intersection cache

LazySets.IntersectionCacheType
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)

LazySets.IntersectionArrayType
IntersectionArray{N, S<:LazySet{N}} <: LazySet{N}

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

Fields

  • array – array of sets

Notes

This type assumes that the dimensions of all elements match.

The EmptySet is the absorbing element for IntersectionArray.

The intersection preserves convexity: if the set arguments are convex, then their intersection is convex as well.

source
Base.:∩Method
∩(X::LazySet, Xs::LazySet...)
∩(Xs::Vector{<:LazySet})

Alias for the n-ary lazy intersection.

source
LazySets.dimMethod
dim(ia::IntersectionArray)

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

Input

  • ia – intersection of a finite number of sets

Output

The ambient dimension of the intersection of a finite number of sets, or 0 if there is no set in the array.

source
LazySets.σMethod
σ(d::AbstractVector, ia::IntersectionArray)

Return a 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 sets

Output

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

Algorithm

This implementation computes the concrete intersection, which can be expensive.

source
LazySets.isboundedMethod
isbounded(ia::IntersectionArray)

Check whether an intersection of a finite number of sets is bounded.

Input

  • ia – intersection of a finite number of 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 LazySets._isbounded_unit_dimensions.

source
Base.:∈Method
∈(x::AbstractVector, ia::IntersectionArray)

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

Input

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

Output

true iff $x ∈ ia$.

Algorithm

A point $x$ is in the intersection iff it is in each set.

source
LazySets.arrayMethod
array(ia::IntersectionArray)

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

Input

  • ia – intersection of a finite number of sets

Output

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

source
LazySets.constraints_listMethod
constraints_list(ia::IntersectionArray)

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: