Binary Functions on Sets

Binary Functions on Sets

This section of the manual describes the binary functions for set types.

Check for emptiness of intersection

LazySets.isdisjointFunction.
isdisjoint(X, Y)

An alternative name for is_intersection_empty(X, Y).

source
is_intersection_empty(X::LazySet{N},
                      Y::LazySet{N},
                      witness::Bool=false
                      )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether two sets do not intersect, and otherwise optionally compute a witness.

Input

  • X – set
  • Y – another set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $X ∩ Y = ∅$
  • If witness option is activated:
    • (true, []) iff $X ∩ Y = ∅$
    • (false, v) iff $X ∩ Y ≠ ∅$ and $v ∈ X ∩ Y$

Algorithm

This is a fallback implementation that computes the concrete intersection, intersection, of the given sets.

A witness is constructed using the an_element implementation of the result.

source
is_intersection_empty(H1::AbstractHyperrectangle{N},
                      H2::AbstractHyperrectangle{N},
                      witness::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether two hyperrectangles do not intersect, and otherwise optionally compute a witness.

Input

  • H1 – first hyperrectangle
  • H2 – second hyperrectangle
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $H1 ∩ H2 = ∅$
  • If witness option is activated:
    • (true, []) iff $H1 ∩ H2 = ∅$
    • (false, v) iff $H1 ∩ H2 ≠ ∅$ and $v ∈ H1 ∩ H2$

Algorithm

$H1 ∩ H2 ≠ ∅$ iff $|c_2 - c_1| ≤ r_1 + r_2$, where $≤$ is taken component-wise.

A witness is computed by starting in one center and moving toward the other center for as long as the minimum of the radius and the center distance. In other words, the witness is the point in H1 that is closest to the center of H2.

source
is_intersection_empty(X::LazySet{N},
                      S::AbstractSingleton{N},
                      witness::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a convex set and a singleton do not intersect, and otherwise optionally compute a witness.

Input

  • X – convex set
  • S – singleton
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $S ∩ X = ∅$
  • If witness option is activated:
    • (true, []) iff $S ∩ X = ∅$
    • (false, v) iff $S ∩ X ≠ ∅$ and v = element(S) $∈ S ∩ X$

Algorithm

$S ∩ X = ∅$ iff element(S) $∉ X$.

source
is_intersection_empty(H::AbstractHyperrectangle{N},
                      S::AbstractSingleton{N},
                      witness::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a hyperrectangle and a singleton do not intersect, and otherwise optionally compute a witness.

Input

  • H – hyperrectangle
  • S – singleton
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $H ∩ S = ∅$
  • If witness option is activated:
    • (true, []) iff $H ∩ S = ∅$
    • (false, v) iff $H ∩ S ≠ ∅$ and v = element(S) $∈ H ∩ S$

Algorithm

$H ∩ S = ∅$ iff element(S) $∉ H$.

source
is_intersection_empty(S1::AbstractSingleton{N},
                      S2::AbstractSingleton{N},
                      witness::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether two singletons do not intersect, and otherwise optionally compute a witness.

Input

  • S1 – first singleton
  • S2 – second singleton
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $S1 ∩ S2 = ∅$
  • If witness option is activated:
    • (true, []) iff $S1 ∩ S2 = ∅$
    • (false, v) iff $S1 ∩ S2 ≠ ∅$ and v = element(S1) $∈ S1 ∩ S2$

Algorithm

$S1 ∩ S2 = ∅$ iff $S1 ≠ S2$.

source
is_intersection_empty(Z::Zonotope{N}, H::Hyperplane{N}, witness::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a zonotope and a hyperplane do not intersect, and otherwise optionally compute a witness.

Input

  • Z – zonotope
  • H – hyperplane
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $Z ∩ H = ∅$
  • If witness option is activated:
    • (true, []) iff $Z ∩ H = ∅$
    • (false, v) iff $Z ∩ H ≠ ∅$ and $v ∈ Z ∩ H$

Algorithm

$Z ∩ H = ∅$ iff $(b - a⋅c) ∉ \left[ ± ∑_{i=1}^p |a⋅g_i| \right]$, where $a$, $b$ are the hyperplane coefficients, $c$ is the zonotope's center, and $g_i$ are the zonotope's generators.

For witness production we fall back to a less efficient implementation for general sets as the first argument.

source
is_intersection_empty(B1::Ball2{N},
                      B2::Ball2{N},
                      witness::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:AbstractFloat}

Check whether two balls in the 2-norm do not intersect, and otherwise optionally compute a witness.

Input

  • B1 – first ball in the 2-norm
  • B2 – second ball in the 2-norm
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $B1 ∩ B2 = ∅$
  • If witness option is activated:
    • (true, []) iff $B1 ∩ B2 = ∅$
    • (false, v) iff $B1 ∩ B2 ≠ ∅$ and $v ∈ B1 ∩ B2$

Algorithm

$B1 ∩ B2 = ∅$ iff $‖ c_2 - c_1 ‖_2 > r_1 + r_2$.

A witness is computed depending on the smaller/bigger ball (to break ties, choose B1 for the smaller ball) as follows.

  • If the smaller ball's center is contained in the bigger ball, we return it.
  • Otherwise start in the smaller ball's center and move toward the other center until hitting the smaller ball's border. In other words, the witness is the point in the smaller ball that is closest to the center of the bigger ball.
source
is_intersection_empty(ls1::LineSegment{N},
                      ls2::LineSegment{N},
                      witness::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether two line segments do not intersect, and otherwise optionally compute a witness.

Input

  • ls1 – first line segment
  • ls2 – second line segment
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $ls1 ∩ ls2 = ∅$
  • If witness option is activated:
    • (true, []) iff $ls1 ∩ ls2 = ∅$
    • (false, v) iff $ls1 ∩ ls2 ≠ ∅$ and $v ∈ ls1 ∩ ls2$

Algorithm

The algorithm is inspired from here, which again is the special 2D case of a 3D algorithm by Ronald Goldman's article on the Intersection of two lines in three-space in Graphics Gems, Andrew S. (ed.), 1990.

We first check if the two line segments are parallel, and if so, if they are collinear. In the latter case, we check containment of any of the end points in the other line segment. Otherwise the lines are not parallel, so we can solve an equation of the intersection point, if it exists.

source
is_intersection_empty(X::LazySet{N},
                      hp::Union{Hyperplane{N}, Line{N}},
                      [witness]::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a compact set an a hyperplane do not intersect, and otherwise optionally compute a witness.

Input

  • X – compact set
  • hp – hyperplane
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $X ∩ hp = ∅$
  • If witness option is activated:
    • (true, []) iff $X ∩ hp = ∅$
    • (false, v) iff $X ∩ hp ≠ ∅$ and $v ∈ X ∩ hp$

Notes

We assume that X is compact. Otherwise, the support vector queries may fail.

Algorithm

A compact convex set intersects with a hyperplane iff the support function in the negative resp. positive direction of the hyperplane's normal vector $a$ is to the left resp. right of the hyperplane's constraint $b$:

\[-ρ(-a) ≤ b ≤ ρ(a)\]

For witness generation, we compute a line connecting the support vectors to the left and right, and then take the intersection of the line with the hyperplane. We follow this algorithm for the line-hyperplane intersection.

source
is_intersection_empty(X::LazySet{N},
                      hs::HalfSpace{N},
                      [witness]::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a compact set an a half-space do not intersect, and otherwise optionally compute a witness.

Input

  • X – compact set
  • hs – half-space
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $X ∩ hs = ∅$
  • If witness option is activated:
    • (true, []) iff $X ∩ hs = ∅$
    • (false, v) iff $X ∩ hs ≠ ∅$ and $v ∈ X ∩ hs$

Notes

We assume that X is compact. Otherwise, the support vector queries may fail.

Algorithm

A compact convex set intersects with a half-space iff the support vector in the negative direction of the half-space's normal vector $a$ is contained in the half-space: $σ(-a) ∈ hs$. The support vector is thus also a witness.

Optional keyword arguments can be passed to the ρ function. In particular, if X is a lazy intersection, options can be passed to the line search algorithm.

source
is_intersection_empty(hs1::HalfSpace{N},
                      hs2::HalfSpace{N},
                      [witness]::Bool=false
                     )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether two half-spaces do not intersect, and otherwise optionally compute a witness.

Input

  • hs1 – half-space
  • hs2 – half-space
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $hs1 ∩ hs2 = ∅$
  • If witness option is activated:
    • (true, []) iff $hs1 ∩ hs2 = ∅$
    • (false, v) iff $hs1 ∩ hs2 ≠ ∅$ and $v ∈ hs1 ∩ hs2$

Algorithm

Two half-spaces do not intersect if and only if their normal vectors point in the opposite direction and there is a gap between the two defining hyperplanes.

The latter can be checked as follows: Let $hs_1 : a_1⋅x = b_1$ and $hs2 : a_2⋅x = b_2$. Then we already know that $a_2 = -k⋅a_1$ for some positive scaling factor $k$. Let $x_1$ be a point on the defining hyperplane of $hs_1$. We construct a line segment from $x_1$ to the point $x_2$ on the defining hyperplane of $hs_2$ by shooting a ray from $x_1$ with direction $a_1$. Thus we look for a factor $s$ such that $(x_1 + s⋅a_1)⋅a_2 = b_2$. This gives us $s = (b_2 - x_1⋅a_2) / (-k a_1⋅a_1)$. The gap exists if and only if $s$ is positive.

If the normal vectors do not point in opposite directions, then the defining hyperplanes intersect and we can produce a witness as follows. All points $x$ in this intersection satisfy $a_1⋅x = b_1$ and $a_2⋅x = b_2$. Thus we have $(a_1 + a_2)⋅x = b_1+b_2$. We now find a dimension where $a_1 + a_2$ is non-zero, say, $i$. Then the result is a vector with one non-zero entry in dimension $i$, defined as $[0, …, 0, (b_1 + b_2)/(a_1[i] + a_2[i]), 0, …, 0]$. Such a dimension $i$ always exists.

source
is_intersection_empty(X::LazySet{N},
                      P::Union{HPolyhedron{N}, HPolytope{N}, AbstractHPolygon{N}},
                      witness::Bool=false
                     )::Bool where {N<:Real}

Check whether a compact set and a polytope do not intersect.

Input

  • X – compact set
  • P – polytope or polygon in constraint-representation

Output

true iff $X ∩ P = ∅$.

Notes

We assume that X is compact. Otherwise, the support vector queries may fail. Witness production is not supported.

Algorithm

The algorithm relies on the intersection check between the set $X$ and each constraint in $P$. It costs $m$ support vector evaluations of $X$, where $m$ is the number of constraints in $P$.

Note that this method can be used with any set $P$ whose constraints are known.

source
is_intersection_empty(X::LazySet{N},
                      Y::LazySet{N},
                      witness::Bool=false
                      )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether two sets do not intersect, and otherwise optionally compute a witness.

Input

  • X – set
  • Y – another set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $X ∩ Y = ∅$
  • If witness option is activated:
    • (true, []) iff $X ∩ Y = ∅$
    • (false, v) iff $X ∩ Y ≠ ∅$ and $v ∈ X ∩ Y$

Algorithm

This is a fallback implementation that computes the concrete intersection, intersection, of the given sets.

A witness is constructed using the an_element implementation of the result.

source
is_intersection_empty(X::LazySet{N},
                      P::Union{HPolyhedron{N}, HPolytope{N}, AbstractHPolygon{N}},
                      witness::Bool=false
                     )::Bool where {N<:Real}

Check whether a compact set and a polytope do not intersect.

Input

  • X – compact set
  • P – polytope or polygon in constraint-representation

Output

true iff $X ∩ P = ∅$.

Notes

We assume that X is compact. Otherwise, the support vector queries may fail. Witness production is not supported.

Algorithm

The algorithm relies on the intersection check between the set $X$ and each constraint in $P$. It costs $m$ support vector evaluations of $X$, where $m$ is the number of constraints in $P$.

Note that this method can be used with any set $P$ whose constraints are known.

source

Convex hull

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

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

Input

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

Output

The HPolyhedron (resp. HPolytope) obtained by the concrete convex hull of P1 and P2.

Notes

For performance reasons, it is suggested to use the CDDLib.Library() backend for the convex_hull.

For further information on the supported backends see Polyhedra's documentation.

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

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

Input

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

Output

The VPolytope obtained by the concrete convex hull of P1 and P2.

Notes

For performance reasons, it is suggested to use the CDDLib.Library() backend for the convex_hull.

source
convex_hull(P::VPolygon{N}, Q::VPolygon{N};
            [algorithm]::String="monotone_chain")::VPolygon{N} where {N<:Real}

Return the convex hull of two polygons in vertex representation.

Input

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

Output

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

Algorithm

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

source

Intersection of two sets

intersection(L1::Line{N}, L2::Line{N}
            )::Union{Singleton{N}, Line{N}, EmptySet{N}} where {N<:Real}

Return the intersection of two 2D lines.

Input

  • L1 – first line
  • L2 – second line

Output

If the lines are identical, the result is the first line. If the lines are parallel and not identical, the result is the empty set. Otherwise the result is the only intersection point.

Examples

The line $y = -x + 1$ intersected with the line $y = x$:

julia> intersection(Line([-1., 1.], 0.), Line([1., 1.], 1.))
Singleton{Float64,Array{Float64,1}}([0.5, 0.5])

julia> intersection(Line([1., 1.], 1.), Line([1., 1.], 1.))
Line{Float64,Array{Float64,1}}([1.0, 1.0], 1.0)
source
intersection(H1::AbstractHyperrectangle{N},
             H2::AbstractHyperrectangle{N}
            )::Union{<:Hyperrectangle{N}, EmptySet{N}} where {N<:Real}

Return the intersection of two hyperrectangles.

Input

  • H1 – first hyperrectangle
  • H2 – second hyperrectangle

Output

If the hyperrectangles do not intersect, the result is the empty set. Otherwise the result is the hyperrectangle that describes the intersection.

Algorithm

In each isolated direction i we compute the rightmost left border and the leftmost right border of the hyperrectangles. If these borders contradict, then the intersection is empty. Otherwise the result uses these borders in each dimension.

source
intersection(x::Interval{N},
             y::Interval{N}
             )::Union{Interval{N}, EmptySet{N}} where {N<:Real}

Return the intersection of two intervals.

Input

  • x – first interval
  • y – second interval

Output

If the intervals do not intersect, the result is the empty set. Otherwise the result is the interval that describes the intersection.

source
intersection(P1::AbstractHPolygon{N},
             P2::AbstractHPolygon{N}
            )::Union{HPolygon{N}, EmptySet{N}} where {N<:Real}

Return the intersection of two polygons in constraint representation.

Input

  • P1 – first polygon
  • P2 – second polygon

Output

If the polygons do not intersect, the result is the empty set. Otherwise the result is the polygon that describes the intersection.

Algorithm

We just combine the constraints of both polygons. To obtain a linear-time algorithm, we interleave the constraints. If there are two constraints with the same normal vector, we choose the tighter one.

source
intersection(P::HPoly{N},
             hs::HalfSpace{N};
             backend=default_polyhedra_backend(P, N),
             prunefunc=removehredundancy!) where {N<:Real}

Compute the intersection of a polytope in H-representation and a half-space.

Input

  • P – polytope
  • hs – half-space
  • backend – (optional, default: default_polyhedra_backend(P, N)) the polyhedral computations backend, see Polyhedra's documentation for further information
  • prunefunc – (optional, default: removehredundancy!) function to post-process the polytope after adding the additional constraint

Output

The same polytope in H-representation with just one more constraint.

source
intersection(P1::HPoly{N},
             P2::HPoly{N};
             backend=default_polyhedra_backend(P1, N),
             prunefunc=removehredundancy!) where {N<:Real}

Compute the intersection of two polyhedra in H-representation.

Input

  • P1 – polytope
  • P2 – polytope
  • backend – (optional, default: default_polyhedra_backend(P1, N)) the polyhedral computations backend, see Polyhedra's documentation for further information
  • prunefunc – (optional, default: removehredundancy!) function to post-process the polytope after adding the additional constraint

Output

A new same polytope in H-representation with just one more constraint.

source
intersection(P1::HPoly{N},
             P2::VPolytope{N};
             backend=default_polyhedra_backend(P1, N),
             prunefunc=removehredundancy!) where {N<:Real}

Compute the intersection of two polytopes in either H-representation or V-representation.

Input

  • P1 – polytope
  • P2 – polytope
  • backend – (optional, default: default_polyhedra_backend(P1, N)) the polyhedral computations backend, see Polyhedra's documentation for further information
  • prunefunc – (optional, default: removehredundancy!) function to post-process the output of intersect

Output

The polytope obtained by the intersection of P1 and P2.

source
intersection(P1::HPoly{N},
             P2::AbstractPolytope{N};
             backend=default_polyhedra_backend(P1, N),
             prunefunc=removehredundancy!) where {N<:Real}

Compute the intersection of a polyhedron and a polytope.

Input

  • P1 – polyhedron
  • P2 – polytope
  • backend – (optional, default: default_polyhedra_backend(P1, N)) the polyhedral computations backend, see Polyhedra's documentation for further information
  • prunefunc – (optional, default: removehredundancy!) function to post-process the output of intersect

Output

The polytope in H-representation obtained by the intersection of P1 and P2.

source
intersection(P1::S1, P2::S2) where {S1<:AbstractPolytope{N},
                                    S2<:AbstractPolytope{N}} where {N<:Real}

Compute the intersection of two polytopic sets.

Input

  • P1 – polytope
  • P2 – another polytope

Output

The polytope obtained by the intersection of P1 and P2. Usually the V-representation is used.

Notes

This fallback implementation requires Polyhedra to evaluate the concrete intersection. Inputs that are not of type HPolytope or VPolytope are converted to an HPolytope through the constraints_list function.

source

Subset check

Base.:⊆Method.
⊆(S::LazySet{N}, H::AbstractHyperrectangle{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a convex set is contained in a hyperrectangular set, and if not, optionally compute a witness.

Input

  • S – inner convex set
  • H – outer hyperrectangular set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $S ⊆ H$
  • If witness option is activated:
    • (true, []) iff $S ⊆ H$
    • (false, v) iff $S \not\subseteq H$ and $v ∈ S \setminus H$

Algorithm

$S ⊆ H$ iff $\operatorname{ihull}(S) ⊆ H$, where $\operatorname{ihull}$ is the interval hull operator.

source
Base.:⊆Method.
⊆(P::AbstractPolytope{N}, S::LazySet{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a polytope is contained in a convex set, and if not, optionally compute a witness.

Input

  • P – inner polytope
  • S – outer convex set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $P ⊆ S$
  • If witness option is activated:
    • (true, []) iff $P ⊆ S$
    • (false, v) iff $P \not\subseteq S$ and $v ∈ P \setminus S$

Algorithm

Since $S$ is convex, $P ⊆ S$ iff $v_i ∈ S$ for all vertices $v_i$ of $P$.

source
Base.:⊆Method.
⊆(P::AbstractPolytope{N}, H::AbstractHyperrectangle, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a polytope is contained in a hyperrectangular set, and if not, optionally compute a witness.

Input

  • P – inner polytope
  • H – outer hyperrectangular set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $P ⊆ H$
  • If witness option is activated:
    • (true, []) iff $P ⊆ H$
    • (false, v) iff $P \not\subseteq H$ and $v ∈ P \setminus H$

Notes

This copy-pasted method just exists to avoid method ambiguities.

Algorithm

Since $H$ is convex, $P ⊆ H$ iff $v_i ∈ H$ for all vertices $v_i$ of $P$.

source
Base.:⊆Method.
⊆(H1::AbstractHyperrectangle{N},
  H2::AbstractHyperrectangle{N},
  [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a given hyperrectangular set is contained in another hyperrectangular set, and if not, optionally compute a witness.

Input

  • H1 – inner hyperrectangular set
  • H2 – outer hyperrectangular set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $H1 ⊆ H2$
  • If witness option is activated:
    • (true, []) iff $H1 ⊆ H2$
    • (false, v) iff $H1 \not\subseteq H2$ and $v ∈ H1 \setminus H2$

Algorithm

$H1 ⊆ H2$ iff $c_1 + r_1 ≤ c_2 + r_2 ∧ c_1 - r_1 ≥ c_2 - r_2$ iff $r_1 - r_2 ≤ c_1 - c_2 ≤ -(r_1 - r_2)$, where $≤$ is taken component-wise.

source
Base.:⊆Method.
⊆(S::AbstractSingleton{N}, set::LazySet{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a given set with a single value is contained in a convex set, and if not, optionally compute a witness.

Input

  • S – inner set with a single value
  • set – outer convex set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $S ⊆ \text{set}$
  • If witness option is activated:
    • (true, []) iff $S ⊆ \text{set}$
    • (false, v) iff $S \not\subseteq \text{set}$ and $v ∈ S \setminus \text{set}$
source
Base.:⊆Method.
⊆(S::AbstractSingleton{N},
  H::AbstractHyperrectangle{N},
  [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a given set with a single value is contained in a hyperrectangular set, and if not, optionally compute a witness.

Input

  • S – inner set with a single value
  • H – outer hyperrectangular set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $S ⊆ H$
  • If witness option is activated:
    • (true, []) iff $S ⊆ H$
    • (false, v) iff $S \not\subseteq H$ and $v ∈ S \setminus H$

Notes

This copy-pasted method just exists to avoid method ambiguities.

source
Base.:⊆Method.
⊆(S1::AbstractSingleton{N},
  S2::AbstractSingleton{N},
  [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a given set with a single value is contained in another set with a single value, and if not, optionally compute a witness.

Input

  • S1 – inner set with a single value
  • S2 – outer set with a single value
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $S1 ⊆ S2$ iff $S1 == S2$
  • If witness option is activated:
    • (true, []) iff $S1 ⊆ S2$
    • (false, v) iff $S1 \not\subseteq S2$ and $v ∈ S1 \setminus S2$
source
Base.:⊆Method.
⊆(B1::Ball2{N}, B2::Ball2{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:AbstractFloat}

Check whether a ball in the 2-norm is contained in another ball in the 2-norm, and if not, optionally compute a witness.

Input

  • B1 – inner ball in the 2-norm
  • B2 – outer ball in the 2-norm
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $B1 ⊆ B2$
  • If witness option is activated:
    • (true, []) iff $B1 ⊆ B2$
    • (false, v) iff $B1 \not\subseteq B2$ and $v ∈ B1 \setminus B2$

Algorithm

$B1 ⊆ B2$ iff $‖ c_1 - c_2 ‖_2 + r_1 ≤ r_2$

source
Base.:⊆Method.
⊆(B::Union{Ball2{N}, Ballp{N}},
  S::AbstractSingleton{N},
  [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:AbstractFloat}

Check whether a ball in the 2-norm or p-norm is contained in a set with a single value, and if not, optionally compute a witness.

Input

  • B – inner ball in the 2-norm or p-norm
  • S – outer set with a single value
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $B ⊆ S$
  • If witness option is activated:
    • (true, []) iff $B ⊆ S$
    • (false, v) iff $B \not\subseteq S$ and $v ∈ B \setminus S$
source
Base.:⊆Method.
⊆(L::LineSegment{N}, S::LazySet{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a line segment is contained in a convex set, and if not, optionally compute a witness.

Input

  • L – inner line segment
  • S – outer convex set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $L ⊆ S$
  • If witness option is activated:
    • (true, []) iff $L ⊆ S$
    • (false, v) iff $L \not\subseteq S$ and $v ∈ L \setminus S$

Algorithm

Since $S$ is convex, $L ⊆ S$ iff $p ∈ S$ and $q ∈ S$, where $p, q$ are the end points of $L$.

source
Base.:⊆Method.
⊆(L::LineSegment{N}, H::AbstractHyperrectangle{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a line segment is contained in a hyperrectangular set, and if not, optionally compute a witness.

Input

  • L – inner line segment
  • H – outer hyperrectangular set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $L ⊆ H$
  • If witness option is activated:
    • (true, []) iff $L ⊆ H$
    • (false, v) iff $L \not\subseteq H$ and $v ∈ L \setminus H$

Notes

This copy-pasted method just exists to avoid method ambiguities.

Algorithm

Since $H$ is convex, $L ⊆ H$ iff $p ∈ H$ and $q ∈ H$, where $p, q$ are the end points of $L$.

source
Base.:⊆Method.
⊆(x::Interval, y::Interval)

Check whether an interval is contained in another interval.

Input

  • x – interval
  • y – interval

Output

true iff $x ⊆ y$.

source
Base.:⊆Method.
⊆(∅::EmptySet{N}, X::LazySet{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether an empty set is contained in another set.

Input

  • – empty set
  • X – another set
  • witness – (optional, default: false) compute a witness if activated (ignored, just kept for interface reasons)

Output

true.

source
Base.:⊆Method.
⊆(X::LazySet{N}, ∅::EmptySet{N}, [witness]::Bool=false
 )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}

Check whether a set is contained in an empty set.

Input

  • X – another set
  • – empty set
  • witness – (optional, default: false) compute a witness if activated

Output

true iff X is empty.

Algorithm

We rely on isempty(X) for the emptiness check and on an_element(X) for witness production.

source