Check for Disjointness of Sets

The function isdisjoint checks whether the intersection of two sets is empty. It can optionally produce a witness if the intersection is nonempty.

Examples

We use the following four sets for illustration.

using LazySets, LazySets.Approximations, Plots
B1 = Ball1(-ones(2), 1.)
B2 = Ball2(ones(2), 1.)
BI = BallInf(zeros(2), 1.)
H = Hyperrectangle(ones(2), ones(2))
sets = [B1, B2, BI, H]

function plot_sets(sets)
    for S in sets
        println(S)
        plot!(S, 1e-2, fillalpha=0.1)
    end
end

function plot_points(points, prefix)
    for i in eachindex(points)
        p = points[i]
        num_occur = length(findfirst(x -> x == p, points[1:i]))
        x = p[1]
        y = p[2]
        if num_occur == 1
            x += 0.15
        elseif num_occur == 2
            y += 0.15
        elseif num_occur == 3
            x -= 0.15
        else
            y -= 0.15
        end
        plot!(Singleton(p))
        plot!(annotations=(x, y, text("$(prefix)$(i)")))
    end
end

plot1 = plot()
plot_sets(sets)
plot1
Example block output
println(isdisjoint(BI, H))
w1 = isdisjoint(BI, H, true)[2]
println(isdisjoint(B1, BI))
w2 = isdisjoint(B1, BI, true)[2]
println(isdisjoint(B1, H))
false
false
true
witnesses = [w1, w2]

plot1 = plot()
plot_sets(sets)
plot_points(witnesses, "w")
plot1
Example block output

Methods

Base.isdisjointFunction
isdisjoint(X::LazySet, Y::LazySet, [witness]::Bool=false)

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

Input

  • X – set
  • Y – 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
Base.isdisjointFunction
isdisjoint(H1::AbstractHyperrectangle, H2::AbstractHyperrectangle,
           [witness]::Bool=false)

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

Input

  • H1 – hyperrectangular set
  • H2 – 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 ∩ 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
Base.isdisjointFunction
isdisjoint(X::LazySet, S::AbstractSingleton, [witness]::Bool=false)

Check whether a set and a set with a single value do not intersect, and otherwise optionally compute a witness.

Input

  • X – set
  • S – set with a single value
  • 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
Base.isdisjointFunction
isdisjoint(S1::AbstractSingleton, S2::AbstractSingleton,
           [witness]::Bool=false)

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

Input

  • S1 – set with a single value
  • S2 – set with a single value
  • 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
Base.isdisjointFunction
isdisjoint(Z::AbstractZonotope, H::Hyperplane, [witness]::Bool=false)

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

Input

  • Z – zonotopic set
  • 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
Base.isdisjointFunction
isdisjoint(L1::LineSegment, L2::LineSegment, [witness]::Bool=false)

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

Input

  • L1 – line segment
  • L2 – line segment
  • witness – (optional, default: false) compute a witness if activated

Output

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

Algorithm

The algorithm is inspired from here, which again is the special 2D case of a 3D algorithm from [1].

We first check if the two line segments are parallel, and if so, if they are collinear. In the latter case, we check membership 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.

[1] Ronald Goldman. Intersection of two lines in three-space. Graphics Gems

source
Base.isdisjointFunction
isdisjoint(X::LazySet, hp::Hyperplane, [witness]::Bool=false)

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

Input

  • X – convex 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$

Algorithm

A 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, X) ≤ b ≤ ρ(a, X)\]

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
Base.isdisjointFunction
isdisjoint(X::LazySet, hs::HalfSpace, [witness]::Bool=false)

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

Input

  • X – 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$

Algorithm

A set intersects with a half-space iff the support function in the negative direction of the half-space's normal vector $a$ is less than the constraint $b$ of the half-space: $-ρ(-a, X) ≤ b$.

For compact set X, we equivalently have that the support vector in the negative direction $-a$ is contained in the half-space: $σ(-a) ∈ hs$. The support vector is thus also a witness if the sets are not disjoint.

source
Base.isdisjointFunction

Extended help

isdisjoint(H1::HalfSpace, H2::HalfSpace, [witness]::Bool=false)

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 $H1 : a_1⋅x = b_1$ and $H2 : 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 $H1$. 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
Base.isdisjointFunction
isdisjoint(P::AbstractPolyhedron, X::LazySet, [witness]::Bool=false;
           [solver]=nothing, [algorithm]="exact")

Check whether a polyhedral set and another set do not intersect, and otherwise optionally compute a witness.

Input

  • P – polyhedral set
  • X – set (see the Notes section below)
  • witness – (optional, default: false) compute a witness if activated
  • solver – (optional, default: nothing) the backend used to solve the linear program
  • algorithm – (optional, default: "exact") algorithm keyword, one of: * "exact" (exact, uses a feasibility LP) *"sufficient" (sufficient, uses half-space checks)

Output

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

Notes

For algorithm == "exact", we assume that constraints_list(X) is defined. For algorithm == "sufficient", witness production is not supported.

For solver == nothing, we fall back to default_lp_solver(N).

Algorithm

For algorithm == "exact", see isempty(P::HPoly, ::Bool).

For algorithm == "sufficient", we rely on the intersection check between the set X and each constraint in P. This requires one support-function evaluation of X for each constraint of P. With this algorithm, the method may return false even in the case where the intersection is empty. On the other hand, if the algorithm returns true, then it is guaranteed that the intersection is empty.

source
Base.isdisjointFunction
isdisjoint(U::UnionSet, X::LazySet, [witness]::Bool=false)

Check whether a union of two sets and another set do not intersect, and otherwise optionally compute a witness.

Input

  • U – union of two sets
  • X – set
  • witness – (optional, default: false) compute a witness if activated

Output

true iff $\text{U} ∩ X = ∅$.

source
Base.isdisjointFunction
isdisjoint(U::UnionSetArray, X::LazySet, [witness]::Bool=false)

Check whether a union of a finite number of sets and another set do not intersect, and otherwise optionally compute a witness.

Input

  • U – union of a finite number of sets
  • X – set
  • witness – (optional, default: false) compute a witness if activated

Output

true iff $\text{U} ∩ X = ∅$.

source
Base.isdisjointFunction
isdisjoint(U::Universe, X::LazySet, [witness]::Bool=false)

Check whether a universe and another set do not intersect, and otherwise optionally compute a witness.

Input

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

Output

true iff $X ≠ ∅$.

source
Base.isdisjointFunction
isdisjoint(C::Complement, X::LazySet, [witness]::Bool=false)

Check whether the complement of a set and another set do not intersect, and otherwise optionally compute a witness.

Input

  • C – complement of a set
  • X – set
  • witness – (optional, default: false) compute a witness if activated

Output

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

Algorithm

We fall back to X ⊆ C.X, which can be justified as follows:

\[ X ∩ Y^C = ∅ ⟺ X ⊆ Y\]

source
Base.isdisjointFunction
isdisjoint(Z1::AbstractZonotope, Z2::AbstractZonotope,
           [witness]::Bool=false; [solver]=nothing)

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

Input

  • Z1 – zonotopic set
  • Z2 – zonotopic set
  • witness – (optional, default: false) compute a witness if activated
  • solver – (optional, default: nothing) the backend used to solve the linear program

Output

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

Algorithm

The algorithm is taken from [1].

$Z1 ∩ Z2 = ∅$ iff $c_1 - c_2 ∉ Z(0, (g_1, g_2))$ where $c_i$ and $g_i$ are the center and generators of zonotope Zi and $Z(c, g)$ represents the zonotope with center $c$ and generators $g$.

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

source
Base.isdisjointFunction
isdisjoint(cpa::CartesianProductArray, P::AbstractPolyhedron,
           [witness]::Bool=false)

Check whether a polytopic Cartesian product array and a polyhedral set do not intersect, and otherwise optionally compute a witness.

Input

  • cpa – Cartesian products of a finite number of polytopes
  • P – polyhedral set
  • witness – (optional, default: false) compute a witness if activated

Output

  • If witness option is deactivated: true iff $\text{cpa} ∩ P = ∅$
  • If witness option is activated:
    • (true, []) iff $\text{cpa} ∩ P = ∅$
    • (false, v) iff $\text{cpa} ∩ P ≠ ∅$ and $v ∈ \text{cpa} ∩ P$

Algorithm

We first identify the blocks of cpa in which P is constrained. Then we project cpa to those blocks and convert the result to an HPolytope (or HPolyhedron if the set type is not known to be bounded) Q. Finally we determine whether Q and the projected P intersect.

source
Base.isdisjointFunction
isdisjoint(X::CartesianProductArray, Y::CartesianProductArray,
           [witness]::Bool=false)

Check whether two Cartesian products of a finite number of sets with the same block structure do not intersect, and otherwise optionally compute a witness.

Input

  • X – Cartesian products of a finite number of sets
  • Y – Cartesian products of a finite number of sets
  • 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$

Notes

The implementation requires (and checks) that the Cartesian products have the same block structure.

Witness production is currently not supported.

source
Base.isdisjointFunction
isdisjoint(cpa::CartesianProductArray, H::AbstractHyperrectangle,
           [witness]::Bool=false)

Check whether a Cartesian product of a finite number of sets and a hyperrectangular set do not intersect, and otherwise optionally compute a witness.

Input

  • cpa – Cartesian product of a finite number of sets
  • H – hyperrectangular set
  • witness – (optional, default: false) compute a witness if activated

Output

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

Algorithm

The sets cpa and H are disjoint if and only if at least one block of cpa and the corresponding projection of H are disjoint. We perform these checks sequentially.

source
Base.isdisjointFunction
isdisjoint(X::LazySet, Y::LazySet, [witness]::Bool=false)

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

Input

  • X – set
  • Y – 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
isdisjoint(P::AbstractPolyhedron, X::LazySet, [witness]::Bool=false;
           [solver]=nothing, [algorithm]="exact")

Check whether a polyhedral set and another set do not intersect, and otherwise optionally compute a witness.

Input

  • P – polyhedral set
  • X – set (see the Notes section below)
  • witness – (optional, default: false) compute a witness if activated
  • solver – (optional, default: nothing) the backend used to solve the linear program
  • algorithm – (optional, default: "exact") algorithm keyword, one of: * "exact" (exact, uses a feasibility LP) *"sufficient" (sufficient, uses half-space checks)

Output

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

Notes

For algorithm == "exact", we assume that constraints_list(X) is defined. For algorithm == "sufficient", witness production is not supported.

For solver == nothing, we fall back to default_lp_solver(N).

Algorithm

For algorithm == "exact", see isempty(P::HPoly, ::Bool).

For algorithm == "sufficient", we rely on the intersection check between the set X and each constraint in P. This requires one support-function evaluation of X for each constraint of P. With this algorithm, the method may return false even in the case where the intersection is empty. On the other hand, if the algorithm returns true, then it is guaranteed that the intersection is empty.

source