Intersection
Binary intersection (Intersection)
LazySets.Intersection
— TypeIntersection{N, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}
Type that represents the intersection of two sets.
Fields
X
– setY
– setcache
– internal cache for avoiding recomputation; seeIntersectionCache
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,Array{Float64,1}},BallInf{Float64,Array{Float64,1}}}
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,Array{Float64,1},Array{Float64,1}}([0.375, 0.0], [0.125, 0.5])
Base.:∩
— Method∩
Alias for Intersection
.
LazySets.dim
— Methoddim(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.
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
– directioncap
– 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 functions of $X$ and $Y$.
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}}}
Return the support function of the intersection of a compact set and a half-space/hyperplane/line in a given direction.
Input
d
– directioncap
– lazy intersection of a compact set and a half-space/hyperplane/ linealgorithm
– (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 directiond
and the hyperplane's normal vectorn
"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:
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
– directioncap
– intersection of a compact set and a polyhedronkwargs
– 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 the constraints_list
of the polyhedron.
LazySets.ρ
— Methodρ(d::AbstractVector, cap::Intersection{N, S1, S2}; kwargs...
) where {N<:Real, S1<:AbstractPolyhedron{N}, S2<:AbstractPolyhedron{N}}
Return an upper bound on the support function of the intersection between two polyhedral sets.
Input
d
– directioncap
– intersection of two polyhedral setskwargs
– additional arguments that are passed to the support-function algorithm
Output
The support function for the given direction.
Algorithm
We combine the constraints of the two polyhedra to a new HPolyhedron
, for which we then evaluate the support function.
LazySets.σ
— Methodσ(d::AbstractVector, cap::Intersection)
Return the support vector of an intersection of two sets in a given direction.
Input
d
– directioncap
– intersection of two sets
Output
The support vector in the given direction.
LazySets.isbounded
— Methodisbounded(cap::Intersection)
Determine 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
.
Base.isempty
— Methodisempty(cap::Intersection)
Return if the intersection of two sets is empty or not.
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.
Base.:∈
— Method∈(x::AbstractVector, cap::Intersection)
Check whether a given point is contained in an intersection of two sets.
Input
x
– point/vectorcap
– intersection of two sets
Output
true
iff $x ∈ cap$.
LazySets.constraints_list
— Methodconstraints_list(cap::Intersection)
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_list
s 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.
LazySets.vertices_list
— Methodvertices_list(cap::Intersection)
Return the list of vertices of a lazy intersection of two 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.
LazySets.isempty_known
— Methodisempty_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.
LazySets.set_isempty!
— Methodset_isempty!(cap::Intersection, isempty::Bool)
Set the status of emptiness in the cache.
Input
cap
– intersection of two setsisempty
– new status of emptiness
LazySets.swap
— Methodswap(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.
LazySets.use_precise_ρ
— Functionuse_precise_ρ(cap::Intersection)
Determine 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.
LazySets._line_search
— Function_line_search(ℓ, X, H::Union{<:HalfSpace, <:Hyperplane, <:Line2D}; [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
ℓ
– directionX
– setH
– 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 >= 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
LazySets._projection
— Function_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 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
ℓ
– directionX
– setH
– hyperplanelazy_linear_map
– (optional, default:false
) to perform the projection lazily or concretelylazy_2d_intersection
– (optional, default:true
) to perform the 2D intersection between the projected set and the line lazily or concretelyalgorithm_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
.
LazySets.linear_map
— Methodlinear_map(M::AbstractMatrix, cap::Intersection)
Return the concrete linear map of a lazy intersection.
Input
M
– matrixcap
– lazy intersection
Output
The set obtained by applying the given linear map to the lazy intersection.
Algorithm
This method computes the concrete intersection.
LazySets.plot_recipe
— Methodplot_recipe(cap::Intersection{N}, [ε]::N=-one(N),
[Nφ]::Int=PLOT_POLAR_DIRECTIONS) where {N<:Real}
Convert a lazy intersection to a pair (x, y)
of points for plotting.
Input
cap
– lazy intersectionε
– (optional, default0
) ignored, used for dispatchNφ
– (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.
RecipesBase.apply_recipe
— Methodplot_intersection(cap::Intersection{N}, [ε]::N=zero(N),
[Nφ]::Int=PLOT_POLAR_DIRECTIONS) where {N}
Plot a lazy intersection.
Input
cap
– lazy intersectionε
– (optional, default0
) ignored, used for dispatchNφ
– (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 LazySets.Approximations.overapproximate
for details).
Examples
julia> using LazySets.Approximations
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 a higher value for Nφ
, 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, -1., 100) # equivalent to the above line
Inherited from LazySet
:
Intersection cache
LazySets.IntersectionCache
— TypeIntersectionCache
Container for information cached by a lazy Intersection
object.
Fields
isempty
– is the intersection empty? There are three possible states, encoded asInt8
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
$n$-ary intersection (IntersectionArray)
LazySets.IntersectionArray
— TypeIntersectionArray{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.
Constructors:
IntersectionArray(array::Vector{<:LazySet})
– default constructorIntersectionArray([n]::Int=0, [N]::Type=Float64)
– constructor for an empty sum with optional size hint and numeric type
LazySets.dim
— Methoddim(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.
LazySets.σ
— Methodσ(d::AbstractVector, ia::IntersectionArray)
Return the support vector of an intersection of a finite number of sets in a given direction.
Input
d
– directionia
– intersection of a finite number of sets
Output
The support vector in the given direction. If the direction has norm zero, the result depends on the individual sets.
LazySets.isbounded
— Methodisbounded(ia::IntersectionArray)
Determine 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
.
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/vectoria
– intersection of a finite number of sets
Output
true
iff $x ∈ ia$.
LazySets.array
— Methodarray(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.
LazySets.constraints_list
— Methodconstraints_list(H::AbstractHyperrectangle{N}) where {N}
Return the list of constraints of an axis-aligned hyperrectangular set.
Input
H
– hyperrectangular set
Output
A list of linear constraints.
constraints_list(P::Ball1{N}) where {N}
Return the list of constraints defining a ball in the 1-norm.
Input
B
– ball in the 1-norm
Output
The list of constraints of the ball.
Algorithm
The constraints can be defined as $d_i^T (x-c) ≤ r$ for all $d_i$, where $d_i$ is a vector with elements $1$ or $-1$ in $n$ dimensions. To span all possible $d_i$, the function Iterators.product
is used.
constraints_list(x::Interval{N}) where {N}
Return the list of constraints of the given interval.
Input
x
– interval
Output
The list of constraints of the interval represented as two one-dimensional half-spaces.
constraints_list(L::Line{N, VN}) where {N, VN}
Return the list of constraints of a line.
Input
L
– line
Output
A list containing 2n-2
half-spaces whose intersection is L
, where n
is the ambient dimension of L
.
constraints_list(U::Universe{N}) where {N}
Return the list of constraints defining a universe.
Input
U
– universe
Output
The empty list of constraints, as the universe is unconstrained.
constraints_list(P::HParallelotope{N, VN}) where {N, VN}
Return the list of constraints of the given parallelotope.
Input
P
– parallelotope in constraint representation
Output
The list of constraints of P
.
constraints_list(cpa::CartesianProductArray{N}) where {N}
Return the list of constraints of a (polyhedral) Cartesian product of a finite number of sets.
Input
cpa
– Cartesian product array
Output
A list of constraints.
constraints_list(ia::IntersectionArray{N}) where {N}
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_list
s of the sets and remove redundant constraints.
constraints_list(rm::ResetMap{N}) where {N}
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.
constraints_list(rm::ResetMap{N, S}) where {N, 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.
Inherited from LazySet
: