Common Set Representations
This section of the manual describes the basic set representation types.
- Common Set Representations
Balls
Euclidean norm ball
LazySets.Ball2
— Type.Ball2{N<:AbstractFloat} <: AbstractCentrallySymmetric{N}
Type that represents a ball in the 2-norm.
Fields
center
– center of the ball as a real vectorradius
– radius of the ball as a real scalar ($≥ 0$)
Notes
Mathematically, a ball in the 2-norm is defined as the set
where $c ∈ \mathbb{R}^n$ is its center and $r ∈ \mathbb{R}_+$ its radius. Here $‖ ⋅ ‖_2$ denotes the Euclidean norm (also known as 2-norm), defined as $‖ x ‖_2 = \left( \sum\limits_{i=1}^n |x_i|^2 \right)^{1/2}$ for any $x ∈ \mathbb{R}^n$.
Examples
Create a five-dimensional ball B
in the 2-norm centered at the origin with radius 0.5:
julia> B = Ball2(zeros(5), 0.5)
Ball2{Float64}([0.0, 0.0, 0.0, 0.0, 0.0], 0.5)
julia> dim(B)
5
Evaluate B
's support vector in the direction $[1,2,3,4,5]$:
julia> σ([1.,2.,3.,4.,5.], B)
5-element Array{Float64,1}:
0.06741998624632421
0.13483997249264842
0.20225995873897262
0.26967994498529685
0.3370999312316211
LazySets.σ
— Method.σ(d::AbstractVector{N}, B::Ball2{N}) where {N<:AbstractFloat}
Return the support vector of a 2-norm ball in a given direction.
Input
d
– directionB
– ball in the 2-norm
Output
The support vector in the given direction. If the direction has norm zero, the origin is returned.
Notes
Let $c$ and $r$ be the center and radius of a ball $B$ in the 2-norm, respectively. For nonzero direction $d$ we have
This function requires computing the 2-norm of the input direction, which is performed in the given precision of the numeric datatype of both the direction and the set. Exact inputs are not supported.
Base.:∈
— Method.∈(x::AbstractVector{N}, B::Ball2{N})::Bool where {N<:AbstractFloat}
Check whether a given point is contained in a ball in the 2-norm.
Input
x
– point/vectorB
– ball in the 2-norm
Output
true
iff $x ∈ B$.
Notes
This implementation is worst-case optimized, i.e., it is optimistic and first computes (see below) the whole sum before comparing to the radius. In applications where the point is typically far away from the ball, a fail-fast implementation with interleaved comparisons could be more efficient.
Algorithm
Let $B$ be an $n$-dimensional ball in the 2-norm with radius $r$ and let $c_i$ and $x_i$ be the ball's center and the vector $x$ in dimension $i$, respectively. Then $x ∈ B$ iff $\left( ∑_{i=1}^n |c_i - x_i|^2 \right)^{1/2} ≤ r$.
Examples
julia> B = Ball2([1., 1.], sqrt(0.5))
Ball2{Float64}([1.0, 1.0], 0.7071067811865476)
julia> ∈([.5, 1.6], B)
false
julia> ∈([.5, 1.5], B)
true
LazySets.center
— Method.center(B::Ball2{N})::Vector{N} where {N<:AbstractFloat}
Return the center of a ball in the 2-norm.
Input
B
– ball in the 2-norm
Output
The center of the ball in the 2-norm.
Base.rand
— Method.rand(::Type{Ball2}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Ball2{N}
Create a random ball in the 2-norm.
Input
Ball2
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random ball in the 2-norm.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1. Additionally, the radius is nonnegative.
LazySets.translate
— Method.translate(B::Ball2{N}, v::AbstractVector{N}) where {N<:AbstractFloat}
Translate (i.e., shift) a ball in the 2-norm by a given vector.
Input
B
– ball in the 2-normv
– translation vector
Output
A translated ball in the 2-norm.
Algorithm
We add the vector to the center of the ball.
Inherited from LazySet
:
Inherited from AbstractCentrallySymmetric
:
Infinity norm ball
LazySets.BallInf
— Type.BallInf{N<:Real} <: AbstractHyperrectangle{N}
Type that represents a ball in the infinity norm.
Fields
center
– center of the ball as a real vectorradius
– radius of the ball as a real scalar ($≥ 0$)
Notes
Mathematically, a ball in the infinity norm is defined as the set
where $c ∈ \mathbb{R}^n$ is its center and $r ∈ \mathbb{R}_+$ its radius. Here $‖ ⋅ ‖_∞$ denotes the infinity norm, defined as $‖ x ‖_∞ = \max\limits_{i=1,…,n} \vert x_i \vert$ for any $x ∈ \mathbb{R}^n$.
Examples
Create the two-dimensional unit ball and compute its support function along the positive $x=y$ direction:
julia> B = BallInf(zeros(2), 1.0)
BallInf{Float64}([0.0, 0.0], 1.0)
julia> dim(B)
2
julia> ρ([1., 1.], B)
2.0
LazySets.center
— Method.center(B::BallInf{N})::Vector{N} where {N<:Real}
Return the center of a ball in the infinity norm.
Input
B
– ball in the infinity norm
Output
The center of the ball in the infinity norm.
LazySets.radius
— Function.radius(B::BallInf, [p]::Real=Inf)::Real
Return the radius of a ball in the infinity norm.
Input
B
– ball in the infinity normp
– (optional, default:Inf
) norm
Output
A real number representing the radius.
Notes
The radius is defined as the radius of the enclosing ball of the given $p$-norm of minimal volume with the same center.
LazySets.radius_hyperrectangle
— Method.radius_hyperrectangle(B::BallInf{N})::Vector{N} where {N<:Real}
Return the box radius of a infinity norm ball, which is the same in every dimension.
Input
B
– infinity norm ball
Output
The box radius of the ball in the infinity norm.
LazySets.radius_hyperrectangle
— Method.radius_hyperrectangle(B::BallInf{N}, i::Int)::N where {N<:Real}
Return the box radius of a infinity norm ball in a given dimension.
Input
B
– infinity norm balli
– dimension of interest
Output
The box radius of the ball in the infinity norm in the given dimension.
Base.rand
— Method.rand(::Type{BallInf}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::BallInf{N}
Create a random ball in the infinity norm.
Input
BallInf
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random ball in the infinity norm.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1. Additionally, the radius is nonnegative.
LazySets.translate
— Method.translate(B::BallInf{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a ball in the infinity norm by a given vector.
Input
B
– ball in the infinity normv
– translation vector
Output
A translated ball in the infinity norm.
Algorithm
We add the vector to the center of the ball.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
:
Inherited from AbstractHyperrectangle
:
Manhattan norm ball
LazySets.Ball1
— Type.Ball1{N<:Real} <: AbstractCentrallySymmetricPolytope{N}
Type that represents a ball in the 1-norm (also known as the Manhattan norm). The ball is also known as a cross-polytope.
It is defined as the set
where $c ∈ \mathbb{R}^n$ is its center and $r ∈ \mathbb{R}_+$ its radius.
Fields
center
– center of the ball as a real vectorradius
– radius of the ball as a scalar ($≥ 0$)
Examples
Unit ball in the 1-norm in the plane:
julia> B = Ball1(zeros(2), 1.)
Ball1{Float64}([0.0, 0.0], 1.0)
julia> dim(B)
2
We evaluate the support vector in the East direction:
julia> σ([0.,1], B)
2-element Array{Float64,1}:
0.0
1.0
LazySets.σ
— Method.σ(d::AbstractVector{N}, B::Ball1{N}) where {N<:Real}
Return the support vector of a ball in the 1-norm in a given direction.
Input
d
– directionB
– ball in the 1-norm
Output
Support vector in the given direction.
Base.:∈
— Method.∈(x::AbstractVector{N}, B::Ball1{N})::Bool where {N<:Real}
Check whether a given point is contained in a ball in the 1-norm.
Input
x
– point/vectorB
– ball in the 1-norm
Output
true
iff $x ∈ B$.
Notes
This implementation is worst-case optimized, i.e., it is optimistic and first computes (see below) the whole sum before comparing to the radius. In applications where the point is typically far away from the ball, a fail-fast implementation with interleaved comparisons could be more efficient.
Algorithm
Let $B$ be an $n$-dimensional ball in the 1-norm with radius $r$ and let $c_i$ and $x_i$ be the ball's center and the vector $x$ in dimension $i$, respectively. Then $x ∈ B$ iff $∑_{i=1}^n |c_i - x_i| ≤ r$.
Examples
julia> B = Ball1([1., 1.], 1.);
julia> ∈([.5, -.5], B)
false
julia> ∈([.5, 1.5], B)
true
LazySets.vertices_list
— Method.vertices_list(B::Ball1{N})::Vector{Vector{N}} where {N<:Real}
Return the list of vertices of a ball in the 1-norm.
Input
B
– ball in the 1-norm
Output
A list containing the vertices of the ball in the 1-norm.
LazySets.center
— Method.center(B::Ball1{N})::Vector{N} where {N<:Real}
Return the center of a ball in the 1-norm.
Input
B
– ball in the 1-norm
Output
The center of the ball in the 1-norm.
Base.rand
— Method.rand(::Type{Ball1}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Ball1{N}
Create a random ball in the 1-norm.
Input
Ball1
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random ball in the 1-norm.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1. Additionally, the radius is nonnegative.
LazySets.constraints_list
— Method.constraints_list(P::Ball1{N})::Vector{LinearConstraint{N}} where {N<:Real}
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.
LazySets.translate
— Method.translate(B::Ball1{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a ball in the 1-norm by a given vector.
Input
B
– ball in the 1-normv
– translation vector
Output
A translated ball in the 1-norm.
Algorithm
We add the vector to the center of the ball.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
:
p-norm ball
LazySets.Ballp
— Type.Ballp{N<:AbstractFloat} <: AbstractCentrallySymmetric{N}
Type that represents a ball in the p-norm, for $1 ≤ p ≤ ∞$.
It is defined as the set
where $c ∈ \mathbb{R}^n$ is its center and $r ∈ \mathbb{R}_+$ its radius. Here $‖ ⋅ ‖_p$ for $1 ≤ p ≤ ∞$ denotes the vector $p$-norm, defined as $‖ x ‖_p = \left( \sum\limits_{i=1}^n |x_i|^p \right)^{1/p}$ for any $x ∈ \mathbb{R}^n$.
Fields
p
– norm as a real scalarcenter
– center of the ball as a real vectorradius
– radius of the ball as a scalar ($≥ 0$)
Notes
The special cases $p=1$, $p=2$ and $p=∞$ fall back to the specialized types Ball1
, Ball2
and BallInf
, respectively.
Examples
A five-dimensional ball in the $p=3/2$ norm centered at the origin of radius 0.5:
julia> B = Ballp(3/2, zeros(5), 0.5)
Ballp{Float64}(1.5, [0.0, 0.0, 0.0, 0.0, 0.0], 0.5)
julia> dim(B)
5
We evaluate the support vector in direction $[1,2,…,5]$:
julia> σ([1., 2, 3, 4, 5], B)
5-element Array{Float64,1}:
0.013516004434607558
0.05406401773843023
0.12164403991146802
0.21625607095372093
0.33790011086518895
LazySets.σ
— Method.σ(d::AbstractVector{N}, B::Ballp{N}) where {N<:AbstractFloat}
Return the support vector of a Ballp
in a given direction.
Input
d
– directionB
– ball in the p-norm
Output
The support vector in the given direction. If the direction has norm zero, the center of the ball is returned.
Algorithm
The support vector of the unit ball in the $p$-norm along direction $d$ is:
where $\tilde{v}_i = \frac{|d_i|^q}{d_i}$ if $d_i ≠ 0$ and $tilde{v}_i = 0$ otherwise, for all $i=1,…,n$, and $q$ is the conjugate number of $p$. By the affine transformation $x = r\tilde{x} + c$, one obtains that the support vector of $\mathcal{B}_p^n(c, r)$ is
where $v_i = c_i + r\frac{|d_i|^q}{d_i}$ if $d_i ≠ 0$ and $v_i = 0$ otherwise, for all $i = 1, …, n$.
Base.:∈
— Method.∈(x::AbstractVector{N}, B::Ballp{N})::Bool where {N<:AbstractFloat}
Check whether a given point is contained in a ball in the p-norm.
Input
x
– point/vectorB
– ball in the p-norm
Output
true
iff $x ∈ B$.
Notes
This implementation is worst-case optimized, i.e., it is optimistic and first computes (see below) the whole sum before comparing to the radius. In applications where the point is typically far away from the ball, a fail-fast implementation with interleaved comparisons could be more efficient.
Algorithm
Let $B$ be an $n$-dimensional ball in the p-norm with radius $r$ and let $c_i$ and $x_i$ be the ball's center and the vector $x$ in dimension $i$, respectively. Then $x ∈ B$ iff $\left( ∑_{i=1}^n |c_i - x_i|^p \right)^{1/p} ≤ r$.
Examples
julia> B = Ballp(1.5, [1., 1.], 1.)
Ballp{Float64}(1.5, [1.0, 1.0], 1.0)
julia> ∈([.5, -.5], B)
false
julia> ∈([.5, 1.5], B)
true
LazySets.center
— Method.center(B::Ballp{N})::Vector{N} where {N<:AbstractFloat}
Return the center of a ball in the p-norm.
Input
B
– ball in the p-norm
Output
The center of the ball in the p-norm.
Base.rand
— Method.rand(::Type{Ballp}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Ballp{N}
Create a random ball in the p-norm.
Input
Ballp
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random ball in the p-norm.
Algorithm
The center and radius are normally distributed with mean 0 and standard deviation 1. Additionally, the radius is nonnegative. The p-norm is a normally distributed number ≥ 1 with mean 1 and standard deviation 1.
LazySets.translate
— Method.translate(B::Ballp{N}, v::AbstractVector{N}) where {N<:AbstractFloat}
Translate (i.e., shift) a ball in the p-norm by a given vector.
Input
B
– ball in the p-normv
– translation vector
Output
A translated ball in the p- norm.
Algorithm
We add the vector to the center of the ball.
Inherited from LazySet
:
Inherited from AbstractCentrallySymmetric
:
Ellipsoid
LazySets.Ellipsoid
— Type.Ellipsoid{N<:AbstractFloat} <: AbstractCentrallySymmetric{N}
Type that represents an ellipsoid.
It is defined as the set
where $c \in \mathbb{R}^n$ is its center and $Q \in \mathbb{R}^{n×n}$ its shape matrix, which should be a positive definite matrix. An ellipsoid can also be characterized as the image of a Euclidean ball by an invertible linear transformation. It is the higher-dimensional generalization of an ellipse.
Fields
center
– center of the ellipsoidshape matrix
– real positive definite matrix, i.e. it is equal to its transpose and $x^\mathrm{T}Qx > 0$ for all nonzero $x$
Examples
If the center is not specified, it is assumed that the center is the origin. For instance, a 3D ellipsoid with center at the origin and the shape matrix being the identity can be created with:
julia> E = Ellipsoid(Matrix{Float64}(I, 3, 3))
Ellipsoid{Float64}([0.0, 0.0, 0.0], [1.0 0.0 0.0; 0.0 1.0 0.0; 0.0 0.0 1.0])
julia> dim(E)
3
julia> an_element(E)
3-element Array{Float64,1}:
0.0
0.0
0.0
This ellipsoid corresponds to the unit Euclidean ball. Let's evaluate its support vector in a given direction:
julia> σ(ones(3), E)
3-element Array{Float64,1}:
0.5773502691896258
0.5773502691896258
0.5773502691896258
A two-dimensional ellipsoid with given center and shape matrix:
julia> E = Ellipsoid(ones(2), Diagonal([2.0, 0.5]))
Ellipsoid{Float64}([1.0, 1.0], [2.0 0.0; 0.0 0.5])
LazySets.ρ
— Method.ρ(d::AbstractVector{N}, E::Ellipsoid{N}) where {N<:AbstractFloat}
Return the support function of an ellipsoid in a given direction.
Input
d
– directionE
– ellipsoid
Output
The support function of the ellipsoid in the given direction.
Algorithm
The support value is $cᵀ d + ‖Bᵀ d‖₂$ where $c$ is the center and $Q = B Bᵀ$ is the shape matrix of E
.
LazySets.σ
— Method.σ(d::AbstractVector{N}, E::Ellipsoid{N}) where {N<:AbstractFloat}
Return the support vector of an ellipsoid in a given direction.
Input
d
– directionE
– ellipsoid
Output
Support vector in the given direction.
Algorithm
Let $E$ be an ellipsoid of center $c$ and shape matrix $Q = BB^\mathrm{T}$. Its support vector along direction $d$ can be deduced from that of the unit Euclidean ball $\mathcal{B}_2$ using the algebraic relations for the support vector,
Base.:∈
— Method.∈(x::AbstractVector{N}, E::Ellipsoid{N})::Bool where {N<:AbstractFloat}
Check whether a given point is contained in an ellipsoid.
Input
x
– point/vectorE
– ellipsoid
Output
true
iff x ∈ E
.
Algorithm
The point $x$ belongs to the ellipsoid of center $c$ and shape matrix $Q$ if and only if
Base.rand
— Method.rand(::Type{Ellipsoid}; [N]::Type{<:AbstractFloat}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Ellipsoid{N}
Create a random ellipsoid.
Input
Ellipsoid
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random ellipsoid.
Algorithm
The center is a normally distributed vector with entries of mean 0 and standard deviation 1.
The idea for the shape matrix comes from here. The matrix is symmetric positive definite, but also diagonally dominant.
where $n$ = dim
(defaults to 2), and $S$ is a $n \times n$ random matrix whose coefficients are uniformly distributed in the interval $[-1, 1]$.
LazySets.center
— Method.center(E::Ellipsoid{N})::Vector{N} where {N<:AbstractFloat}
Return the center of the ellipsoid.
Input
E
– ellipsoid
Output
The center of the ellipsoid.
LazySets.translate
— Method.translate(E::Ellipsoid{N}, v::AbstractVector{N}; share::Bool=false
) where {N<:AbstractFloat}
Translate (i.e., shift) an ellipsoid by a given vector.
Input
E
– ellipsoidv
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated ellipsoid.
Notes
The shape matrix is shared with the original ellipsoid if share == true
.
Algorithm
We add the vector to the center of the ellipsoid.
Inherited from LazySet
:
Inherited from AbstractCentrallySymmetric
:
Empty set
LazySets.EmptySet
— Type.EmptySet{N<:Real} <: LazySet{N}
Type that represents the empty set, i.e., the set with no elements.
LazySets.∅
— Constant.∅
An EmptySet
instance of type Float64
.
LazySets.dim
— Method.dim(∅::EmptySet)
Return the dimension of the empty set, which is -1 by convention.
Input
∅
– an empty set
Output
-1
by convention.
LazySets.σ
— Method.σ(d::AbstractVector{N}, ∅::EmptySet{N}) where {N<:Real}
Return the support vector of an empty set.
Input
d
– direction∅
– an empty set
Output
An error.
Base.:∈
— Method.∈(x::AbstractVector{N}, ∅::EmptySet{N})::Bool where {N<:Real}
Check whether a given point is contained in an empty set.
Input
x
– point/vector∅
– empty set
Output
The output is always false
.
Examples
julia> ∈([1.0, 0.0], ∅)
false
LazySets.an_element
— Method.an_element(∅::EmptySet)
Return some element of an empty set.
Input
∅
– empty set
Output
An error.
Base.rand
— Method.rand(::Type{EmptySet}; [N]::Type{<:Real}=Float64, [dim]::Int=0,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::EmptySet{N}
Create an empty set (note that there is nothing to randomize).
Input
EmptySet
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 0) dimension (is ignored)rng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
The (only) empty set of the given numeric type.
LazySets.isbounded
— Method.isbounded(∅::EmptySet)::Bool
Determine whether an empty set is bounded.
Input
∅
– empty set
Output
true
.
Base.isempty
— Method.isempty(∅::EmptySet)::Bool
Return if the empty set is empty or not.
Input
∅
– empty set
Output
true
.
LinearAlgebra.norm
— Function.norm(S::EmptySet, [p]::Real=Inf)
Return the norm of an empty set. It is the norm of the enclosing ball (of the given $p$-norm) of minimal volume that is centered in the origin.
Input
S
– empty setp
– (optional, default:Inf
) norm
Output
An error.
LazySets.radius
— Function.radius(S::EmptySet, [p]::Real=Inf)
Return the radius of an empty set. It is the radius of the enclosing ball (of the given $p$-norm) of minimal volume with the same center.
Input
S
– empty setp
– (optional, default:Inf
) norm
Output
An error.
LazySets.diameter
— Function.diameter(S::EmptySet, [p]::Real=Inf)
Return the diameter of an empty set. It is the maximum distance between any two elements of the set, or, equivalently, the diameter of the enclosing ball (of the given $p$-norm) of minimal volume with the same center.
Input
S
– empty setp
– (optional, default:Inf
) norm
Output
An error.
LazySets.linear_map
— Method.linear_map(M::AbstractMatrix{N}, ∅::EmptySet{N}) where {N}
Return the linear map of an empty set.
Input
M
– matrix∅
– empty set
Output
The empty set.
LazySets.translate
— Method.translate(∅::EmptySet{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) an empty set by a given vector.
Input
∅
– empty setv
– translation vector
Output
The empty set.
RecipesBase.apply_recipe
— Function.plot_emptyset(∅::EmptySet, [ε::Float64=0.0]; ...)
Plot an empty set.
Input
∅
– empty setε
– (optional, default:0.0
) approximation error bound
Examples
julia> using Plots, LazySets;
julia> plot(∅);
julia> plot(∅, 1e-2);
Inherited from LazySet
:
Half-Space
LazySets.HalfSpace
— Type.HalfSpace{N<:Real} <: AbstractPolyhedron{N}
Type that represents a (closed) half-space of the form $a⋅x ≤ b$.
Fields
a
– normal directionb
– constraint
Examples
The set $y ≥ 0$ in the plane:
julia> HalfSpace([0, -1.], 0.)
HalfSpace{Float64}([0.0, -1.0], 0.0)
LazySets.LinearConstraint
— Type.LinearConstraint
Alias for HalfSpace
LazySets.dim
— Method.dim(hs::HalfSpace)::Int
Return the dimension of a half-space.
Input
hs
– half-space
Output
The ambient dimension of the half-space.
LazySets.ρ
— Method.ρ(d::AbstractVector{N}, hs::HalfSpace{N})::N where {N<:Real}
Evaluate the support function of a half-space in a given direction.
Input
d
– directionhs
– half-space
Output
The support function of the half-space. If the set is unbounded in the given direction, the result is Inf
.
LazySets.σ
— Method.σ(d::AbstractVector{N}, hs::HalfSpace{N}) where {N<:Real}
Return the support vector of a half-space.
Input
d
– directionhs
– half-space
Output
The support vector in the given direction, which is only defined in the following two cases:
- The direction has norm zero.
- The direction is the half-space's normal direction.
In both cases the result is any point on the boundary (the defining hyperplane). Otherwise this function throws an error.
Base.:∈
— Method.∈(x::AbstractVector{N}, hs::HalfSpace{N})::Bool where {N<:Real}
Check whether a given point is contained in a half-space.
Input
x
– point/vectorhs
– half-space
Output
true
iff $x ∈ hs$.
Algorithm
We just check if $x$ satisfies $a⋅x ≤ b$.
LazySets.an_element
— Method.an_element(hs::HalfSpace{N})::Vector{N} where {N<:Real}
Return some element of a half-space.
Input
hs
– half-space
Output
An element on the defining hyperplane.
Base.rand
— Method.rand(::Type{HalfSpace}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::HalfSpace{N}
Create a random half-space.
Input
HalfSpace
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random half-space.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1. Additionally, the constraint a
is nonzero.
LazySets.isbounded
— Method.isbounded(hs::HalfSpace)::Bool
Determine whether a half-space is bounded.
Input
hs
– half-space
Output
false
.
Base.isempty
— Method.isempty(hs::HalfSpace)::Bool
Return if a half-space is empty or not.
Input
hs
– half-space
Output
false
.
LazySets.constraints_list
— Method.constraints_list(hs::HalfSpace{N})::Vector{LinearConstraint{N}}
where {N<:Real}
Return the list of constraints of a half-space.
Input
hs
– half-space
Output
A singleton list containing the half-space.
LazySets.constraints_list
— Method.constraints_list(A::AbstractMatrix{N}, b::AbstractVector{N}
)::Vector{LinearConstraint{N}} where {N<:Real}
Convert a matrix-vector representation to a linear-constraint representation.
Input
A
– matrixb
– vector
Output
A list of linear constraints.
LazySets.constrained_dimensions
— Method.constrained_dimensions(hs::HalfSpace{N})::Vector{Int} where {N<:Real}
Return the indices in which a half-space is constrained.
Input
hs
– half-space
Output
A vector of ascending indices i
such that the half-space is constrained in dimension i
.
Examples
A 2D half-space with constraint $x1 ≥ 0$ is constrained in dimension 1 only.
LazySets.translate
— Method.translate(hs::HalfSpace{N}, v::AbstractVector{N}; share::Bool=false
) where {N<:Real}
Translate (i.e., shift) a half-space by a given vector.
Input
hs
– half-spacev
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated half-space.
Notes
The normal vectors of the halfspace (vector a
in a⋅x ≤ b
) is shared with the original halfspace if share == true
.
Algorithm
A half-space $a⋅x ≤ b$ is transformed to the half-space $a⋅x ≤ b + a⋅v$. In other words, we add the dot product $a⋅v$ to $b$.
LazySets.halfspace_left
— Method.halfspace_left(p::AbstractVector{N},
q::AbstractVector{N})::HalfSpace{N} where {N<:Real}
Return a half-space describing the 'left' of a two-dimensional line segment through two points.
Input
p
– first pointq
– second point
Output
The half-space whose boundary goes through the two points p
and q
and which describes the left-hand side of the directed line segment pq
.
Algorithm
The implementation is simple: the half-space $a⋅x ≤ b$ is calculated as a = [dy, -dx]
, where $d = (dx, dy)$ denotes the line segment pq
, that is, $\vec{d} = \vec{p} - \vec{q}$, and b = dot(p, a)
.
Examples
The left half-space of the "east" and "west" directions in two-dimensions are the upper and lower half-spaces:
julia> using LazySets: halfspace_left
julia> halfspace_left([0.0, 0.0], [1.0, 0.0])
HalfSpace{Float64}([0.0, -1.0], 0.0)
julia> halfspace_left([0.0, 0.0], [-1.0, 0.0])
HalfSpace{Float64}([0.0, 1.0], 0.0)
We create a box from the sequence of line segments that describe its edges:
julia> H1 = halfspace_left([-1.0, -1.0], [1.0, -1.0]);
julia> H2 = halfspace_left([1.0, -1.0], [1.0, 1.0]);
julia> H3 = halfspace_left([1.0, 1.0], [-1.0, 1.0]);
julia> H4 = halfspace_left([-1.0, 1.0], [-1.0, -1.0]);
julia> H = HPolygon([H1, H2, H3, H4]);
julia> B = BallInf([0.0, 0.0], 1.0);
julia> B ⊆ H && H ⊆ B
true
LazySets.halfspace_right
— Method.halfspace_right(p::AbstractVector{N},
q::AbstractVector{N})::HalfSpace{N} where {N<:Real}
Return a half-space describing the 'right' of a two-dimensional line segment through two points.
Input
p
– first pointq
– second point
Output
The half-space whose boundary goes through the two points p
and q
and which describes the right-hand side of the directed line segment pq
.
Algorithm
See the documentation of halfspace_left
.
LazySets.tosimplehrep
— Method.tosimplehrep(constraints::AbstractVector{LinearConstraint{N}})
where {N<:Real}
Return the simple H-representation $Ax ≤ b$ from a list of linear constraints.
Input
constraints
– a list of linear constraints
Output
The tuple (A, b)
where A
is the matrix of normal directions and b
is the vector of offsets.
LazySets.remove_redundant_constraints
— Method.remove_redundant_constraints(
constraints::AbstractVector{LinearConstraint{N}};
backend=GLPKSolverLP())::Union{AbstractVector{LinearConstraint{N}},
EmptySet{N}} where {N<:Real}
Remove the redundant constraints of a given list of linear constraints.
Input
constraints
– list of constraintsbackend
– (optional, default:GLPKSolverLP
) numeric LP solver backend
Output
The list of constraints with the redundant ones removed, or an empty set if the constraints are infeasible.
Algorithm
See remove_redundant_constraints!(::AbstractVector{LinearConstraint{<:Real}})
for details.
LazySets.remove_redundant_constraints!
— Method. remove_redundant_constraints!(
constraints::AbstractVector{LinearConstraint{N}};
[backend]=GLPKSolverLP())::Bool where {N<:Real}
Remove the redundant constraints of a given list of linear constraints; the list is updated in-place.
Input
constraints
– list of constraintsbackend
– (optional, default:GLPKSolverLP
) numeric LP solver backend
Output
true
if the function was successful and the list of constraints constraints
is modified by removing the redundant constraints, and false
only if the constraints are infeasible.
Notes
Note that the result may be true
even if the constraints are infeasible. For example, $x ≤ 0 && x ≥ 1$ will return true
without removing any constraint. To check if the constraints are infeasible, use isempty(HPolyhedron(constraints)
.
Algorithm
If there are m
constraints in n
dimensions, this function checks one by one if each of the m
constraints is implied by the remaining ones.
To check if the k
-th constraint is redundant, an LP is formulated using the constraints that have not yet being removed. If, at an intermediate step, it is detected that a subgroup of the constraints is infeasible, this function returns false
. If the calculation finished successfully, this function returns true
.
For details, see Fukuda's Polyhedra FAQ.
Inherited from LazySet
:
Hyperplane
LazySets.Hyperplane
— Type.Hyperplane{N<:Real} <: AbstractPolyhedron{N}
Type that represents a hyperplane of the form $a⋅x = b$.
Fields
a
– normal directionb
– constraint
Examples
The plane $y = 0$:
julia> Hyperplane([0, 1.], 0.)
Hyperplane{Float64}([0.0, 1.0], 0.0)
LazySets.dim
— Method.dim(hp::Hyperplane)::Int
Return the dimension of a hyperplane.
Input
hp
– hyperplane
Output
The ambient dimension of the hyperplane.
LazySets.ρ
— Method.ρ(d::AbstractVector{N}, hp::Hyperplane{N})::N where {N<:Real}
Evaluate the support function of a hyperplane in a given direction.
Input
d
– directionhp
– hyperplane
Output
The support function of the hyperplane. If the set is unbounded in the given direction, the result is Inf
.
LazySets.σ
— Method.σ(d::AbstractVector{N}, hp::Hyperplane{N}) where {N<:Real}
Return the support vector of a hyperplane.
Input
d
– directionhp
– hyperplane
Output
The support vector in the given direction, which is only defined in the following two cases:
- The direction has norm zero.
- The direction is the hyperplane's normal direction or its opposite direction.
In all cases, the result is any point on the hyperplane. Otherwise this function throws an error.
Base.:∈
— Method.∈(x::AbstractVector{N}, hp::Hyperplane{N})::Bool where {N<:Real}
Check whether a given point is contained in a hyperplane.
Input
x
– point/vectorhp
– hyperplane
Output
true
iff $x ∈ hp$.
Algorithm
We just check if $x$ satisfies $a⋅x = b$.
LazySets.an_element
— Method.an_element(hp::Hyperplane{N})::Vector{N} where {N<:Real}
Return some element of a hyperplane.
Input
hp
– hyperplane
Output
An element on the hyperplane.
Base.rand
— Method.rand(::Type{Hyperplane}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Hyperplane{N}
Create a random hyperplane.
Input
Hyperplane
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random hyperplane.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1. Additionally, the constraint a
is nonzero.
LazySets.isbounded
— Method.isbounded(hp::Hyperplane)::Bool
Determine whether a hyperplane is bounded.
Input
hp
– hyperplane
Output
false
.
Base.isempty
— Method.isempty(hp::Hyperplane)::Bool
Return if a hyperplane is empty or not.
Input
hp
– hyperplane
Output
false
.
LazySets.constrained_dimensions
— Method.constrained_dimensions(hp::Hyperplane{N})::Vector{Int} where {N<:Real}
Return the indices in which a hyperplane is constrained.
Input
hp
– hyperplane
Output
A vector of ascending indices i
such that the hyperplane is constrained in dimension i
.
Examples
A 2D hyperplane with constraint $x1 = 0$ is constrained in dimension 1 only.
LazySets.constraints_list
— Method.constraints_list(hp::Hyperplane{N})::Vector{LinearConstraint{N}}
where {N<:Real}
Return the list of constraints of a hyperplane.
Input
hp
– hyperplane
Output
A list containing two half-spaces.
LazySets.translate
— Method.translate(hp::Hyperplane{N}, v::AbstractVector{N}; share::Bool=false
) where {N<:Real}
Translate (i.e., shift) a hyperplane by a given vector.
Input
hp
– hyperplanev
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated hyperplane.
Notes
The normal vectors of the hyperplane (vector a
in a⋅x = b
) is shared with the original hyperplane if share == true
.
Algorithm
A hyperplane $a⋅x = b$ is transformed to the hyperplane $a⋅x = b + a⋅v$. In other words, we add the dot product $a⋅v$ to $b$.
Inherited from LazySet
:
Hyperrectangle
LazySets.Hyperrectangle
— Type.Hyperrectangle{N<:Real} <: AbstractHyperrectangle{N}
Type that represents a hyperrectangle.
A hyperrectangle is the Cartesian product of one-dimensional intervals.
Fields
center
– center of the hyperrectangle as a real vectorradius
– radius of the ball as a real vector, i.e., half of its width along each coordinate direction
Examples
There is also a constructor from lower and upper bounds with keyword arguments high
and low
. The following two constructions are equivalent:
julia> c = ones(2);
julia> r = [0.1, 0.2];
julia> l = [0.9, 0.8];
julia> h = [1.1, 1.2];
julia> Hyperrectangle(c, r)
Hyperrectangle{Float64}([1.0, 1.0], [0.1, 0.2])
julia> Hyperrectangle(low=l, high=h)
Hyperrectangle{Float64}([1.0, 1.0], [0.1, 0.2])
Base.rand
— Method.rand(::Type{Hyperrectangle}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Hyperrectangle{N}
Create a random hyperrectangle.
Input
Hyperrectangle
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random hyperrectangle.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1. Additionally, the radius is nonnegative.
LazySets.center
— Method.center(H::Hyperrectangle{N})::Vector{N} where {N<:Real}
Return the center of a hyperrectangle.
Input
H
– hyperrectangle
Output
The center of the hyperrectangle.
LazySets.radius_hyperrectangle
— Method.radius_hyperrectangle(H::Hyperrectangle{N})::Vector{N} where {N<:Real}
Return the box radius of a hyperrectangle in every dimension.
Input
H
– hyperrectangle
Output
The box radius of the hyperrectangle.
LazySets.radius_hyperrectangle
— Method.radius_hyperrectangle(H::Hyperrectangle{N}, i::Int)::N where {N<:Real}
Return the box radius of a hyperrectangle in a given dimension.
Input
H
– hyperrectanglei
– dimension of interest
Output
The radius in the given dimension.
LazySets.translate
— Method.translate(H::Hyperrectangle{N}, v::AbstractVector{N}; share::Bool=false
) where {N<:Real}
Translate (i.e., shift) a hyperrectangle by a given vector.
Input
H
– hyperrectanglev
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated hyperrectangle.
Notes
The radius vector is shared with the original hyperrectangle if share == true
.
Algorithm
We add the vector to the center of the hyperrectangle.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
:
Inherited from AbstractHyperrectangle
:
Interval
LazySets.Interval
— Type.Interval{N<:Real, IN<:AbstractInterval{N}} <: AbstractHyperrectangle{N}
Type representing an interval on the real line. Mathematically, it is of the form
Fields
dat
– data container for the given interval
Notes
This type relies on the IntervalArithmetic.jl library for representation of intervals and arithmetic operations.
Examples
Unidimensional intervals are symbolic representations of a real closed interval.
We can create intervals in different ways, the simpler way is to pass a pair of numbers:
julia> x = Interval(0.0, 1.0)
Interval{Float64,IntervalArithmetic.Interval{Float64}}([0, 1])
or a 2-vector:
julia> x = Interval([0.0, 1.0])
Interval{Float64,IntervalArithmetic.Interval{Float64}}([0, 1])
Note that if the package IntervalArithmetic
is loaded in the current scope, you have to prepend LazySets.
to the Interval
type since there is a name conflict otherwise.
julia> using IntervalArithmetic
WARNING: using IntervalArithmetic.Interval in module Main conflicts with an existing identifier.
julia> x = Interval(IntervalArithmetic.Interval(0.0, 1.0))
Interval{Float64,IntervalArithmetic.Interval{Float64}}([0, 1])
julia> dim(x)
1
julia> center(x)
1-element Array{Float64,1}:
0.5
This type is such that the usual pairwise arithmetic operators +
, -
, *
trigger the corresponding interval arithmetic backend method, and return a new Interval
object. For the symbolic Minkowksi sum, use MinkowskiSum
or ⊕
.
Interval of other numeric types can be created as well, eg. a rational interval:
julia> Interval(0//1, 2//1)
Interval{Rational{Int64},AbstractInterval{Rational{Int64}}}([0//1, 2//1])
LazySets.dim
— Method.dim(x::Interval)::Int
Return the ambient dimension of an interval.
Input
x
– interval
Output
The integer 1.
LazySets.σ
— Method.σ(d::AbstractVector{N}, x::Interval{N}) where {N<:Real}
Return the support vector of an ellipsoid in a given direction.
Input
d
– directionx
– interval
Output
Support vector in the given direction.
Base.:∈
— Method.∈(v::AbstractVector{N}, x::Interval{N}) where {N<:Real})
Return whether a vector is contained in the interval.
Input
v
– one-dimensional vectorx
– interval
Output
true
iff x
contains v
's first component.
Base.:∈
— Method.∈(v::N, x::Interval{N}) where {N<:Real}
Return whether a number is contained in the interval.
Input
v
– scalarx
– interval
Output
true
iff x
contains v
.
LazySets.an_element
— Method.an_element(x::Interval{N})::Vector{N} where {N<:Real}
Return some element of an interval.
Input
x
– interval
Output
The left border (min(x)
) of the interval.
LazySets.vertices_list
— Method.vertices_list(x::Interval{N})::Vector{Vector{N}} where {N<:Real}
Return the list of vertices of this interval.
Input
x
– interval
Output
The list of vertices of the interval represented as two one-dimensional vectors.
LazySets.translate
— Method.translate(x::Interval{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) an interval by a given vector.
Input
x
– intervalv
– translation vector
Output
A translated interval.
Algorithm
We add the vector to the left and right of the interval.
LazySets.center
— Method.center(x::Interval{N})::Vector{N} where {N<:Real}
Return the interval's center.
Input
x
– interval
Output
The center, or midpoint, of $x$.
Base.min
— Method.min(x::Interval{N})::N where {N<:Real}
Return the lower component of an interval.
Input
x
– interval
Output
The lower (lo
) component of the interval.
Base.max
— Method.max(x::Interval{N})::N where {N<:Real}
Return the higher or upper component of an interval.
Input
x
– interval
Output
The higher (hi
) component of the interval.
LazySets.low
— Method.low(x::Interval{N})::Vector{N} where {N<:Real}
Return the lower coordinate of an interval set.
Input
x
– interval
Output
A vector with the lower coordinate of the interval.
LazySets.high
— Method.high(x::Interval{N})::Vector{N} where {N<:Real}
Return the higher coordinate of an interval set.
Input
x
– interval
Output
A vector with the higher coordinate of the interval.
LazySets.radius_hyperrectangle
— Method.radius_hyperrectangle(x::Interval{N})::Vector{N} where {N<:Real}
Return the box radius of an interval in every dimension.
Input
x
– interval
Output
The box radius of the interval (a one-dimensional vector).
LazySets.radius_hyperrectangle
— Method.radius_hyperrectangle(x::Interval{N}, i::Int)::N where {N<:Real}
Return the box radius of an interval in a given dimension.
Input
x
– intervali
– dimension index (must be1
)
Output
The box radius in the given dimension.
Base.:+
— Method.+(x::Interval{N}, y::Interval{N}) where {N<:Real}
Return the sum of the intervals.
Input
x
– intervaly
– interval
Output
The sum of the intervals as a new Interval
set.
Base.:-
— Method.-(x::Interval{N}, y::Interval{N}) where {N<:Real}
Return the difference of the intervals.
Input
x
– intervaly
– interval
Output
The difference of the intervals as a new Interval
set.
Base.:*
— Method. *(x::Interval{N}, y::Interval{N}) where {N<:Real}
Return the product of the intervals.
Input
x
– intervaly
– interval
Output
The product of the intervals as a new Interval
set.
Base.rand
— Method.rand(::Type{Interval}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Interval{N}
Create a random interval.
Input
Interval
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 1) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random interval.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1.
RecipesBase.apply_recipe
— Method.plot_interval(I::Interval; ...)
Plot an interval.
Input
I
– interval
Examples
julia> using Plots, LazySets;
julia> I = Interval(0.0, 1.0);
julia> plot(I);
RecipesBase.apply_recipe
— Method.plot_intervals(Xk::Vector{S}; ...) where {S<:Interval}
Plot an array of intervals.
Input
Xk
– linear array of intervals
Examples
julia> using Plots, LazySets;
julia> I1 = Interval([0., 1.]);
julia> I2 = Interval([0.5, 2.]);
julia> plot([I1, I2]);
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
:
Inherited from AbstractHyperrectangle
:
Line
LazySets.Line
— Type.Line{N<:Real, VN<:AbstractVector{N}} <: AbstractPolyhedron{N}
Type that represents a line in 2D of the form $a⋅x = b$ (i.e., a special case of a Hyperplane
).
Fields
a
– normal directionb
– constraint
Examples
The line $y = -x + 1$:
julia> Line([1., 1.], 1.)
Line{Float64,Array{Float64,1}}([1.0, 1.0], 1.0)
LazySets.dim
— Method.dim(L::Line)::Int
Return the ambient dimension of a line.
Input
L
– line
Output
The ambient dimension of the line, which is 2.
LazySets.σ
— Method.σ(d::AbstractVector{N}, L::Line{N}) where {N<:Real}
Return the support vector of a line in a given direction.
Input
d
– directionL
– line
Output
The support vector in the given direction, which is defined the same way as for the more general Hyperplane
.
Base.:∈
— Method.∈(x::AbstractVector{N}, L::Line{N})::Bool where {N<:Real}
Check whether a given point is contained in a line.
Input
x
– point/vectorL
– line
Output
true
iff x ∈ L
.
Algorithm
The point $x$ belongs to the line if and only if $a⋅x = b$ holds.
LazySets.an_element
— Method.an_element(L::Line{N})::Vector{N} where {N<:Real}
Return some element of a line.
Input
L
– line
Output
An element on the line.
Algorithm
If the $b$ value of the line is zero, the result is the origin. Otherwise the result is some $x = [x1, x2]$ such that $a·[x1, x2] = b$. We first find out in which dimension $a$ is nonzero, say, dimension 1, and then choose $x1 = 1$ and accordingly $x2 = \frac{b - a1}{a2}$.
Base.rand
— Method.rand(::Type{Line}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Line{N}
Create a random line.
Input
Line
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random line.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1. Additionally, the constraint a
is nonzero.
LazySets.isbounded
— Method.isbounded(L::Line)::Bool
Determine whether a line is bounded.
Input
L
– line
Output
false
.
Base.isempty
— Method.isempty(L::Line)::Bool
Return if a line is empty or not.
Input
L
– line
Output
false
.
LazySets.constrained_dimensions
— Method.constrained_dimensions(L::Line{N})::Vector{Int} where {N<:Real}
Return the indices in which a line is constrained.
Input
L
– line
Output
A vector of ascending indices i
such that the line is constrained in dimension i
.
Examples
A line with constraint $x1 = 0$ is constrained in dimension 1 only.
LazySets.constraints_list
— Method.constraints_list(L::Line{N})::Vector{LinearConstraint{N}}
where {N<:Real}
Return the list of constraints of a line.
Input
L
– line
Output
A list containing two half-spaces.
LazySets.translate
— Method.translate(L::Line{N}, v::AbstractVector{N}; share::Bool=false
) where {N<:Real}
Translate (i.e., shift) a line by a given vector.
Input
L
– linev
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated line.
Notes
The normal vector of the line (vector a
in a⋅x = b
) is shared with the original line if share == true
.
Algorithm
A line $a⋅x = b$ is transformed to the line $a⋅x = b + a⋅v$. In other words, we add the dot product $a⋅v$ to $b$.
Inherited from LazySet
:
Line segment
LazySets.LineSegment
— Type.LineSegment{N<:Real} <: AbstractCentrallySymmetricPolytope{N}
Type that represents a line segment in 2D between two points $p$ and $q$.
Fields
p
– first pointq
– second point
Examples
A line segment along the $x = y$ diagonal:
julia> s = LineSegment([0., 0], [1., 1.])
LineSegment{Float64}([0.0, 0.0], [1.0, 1.0])
julia> dim(s)
2
Use plot(s)
to plot the extreme points of s
and the line segment joining them. Membership test is computed with ∈ (in
):
julia> [0., 0] ∈ s && [.25, .25] ∈ s && [1., 1] ∈ s && !([.5, .25] ∈ s)
true
We can check the intersection with another line segment, and optionally compute a witness (which is just the common point in this case):
julia> sn = LineSegment([1., 0], [0., 1.])
LineSegment{Float64}([1.0, 0.0], [0.0, 1.0])
julia> isempty(s ∩ sn)
false
julia> is_intersection_empty(s, sn, true)
(false, [0.5, 0.5])
LazySets.dim
— Method.dim(L::LineSegment)::Int
Return the ambient dimension of a line segment.
Input
L
– line segment
Output
The ambient dimension of the line segment, which is 2.
LazySets.σ
— Method.σ(d::AbstractVector{N}, L::LineSegment{N}) where {N<:Real}
Return the support vector of a line segment in a given direction.
Input
d
– directionL
– line segment
Output
The support vector in the given direction.
Algorithm
If the angle between the vector $q - p$ and $d$ is bigger than 90° and less than 270° (measured in counter-clockwise order), the result is $p$, otherwise it is $q$. If the angle is exactly 90° or 270°, or if the direction has norm zero, this implementation returns $q$.
Base.:∈
— Method.∈(x::AbstractVector{N}, L::LineSegment{N})::Bool where {N<:Real}
Check whether a given point is contained in a line segment.
Input
x
– point/vectorL
– line segment
Output
true
iff $x ∈ L$.
Algorithm
Let $L = (p, q)$ be the line segment with extremes $p$ and $q$, and let $x$ be the given point.
- A necessary condition for $x ∈ (p, q)$ is that the three points are aligned, thus their cross product should be zero.
- It remains to check that $x$ belongs to the box approximation of $L$. This amounts to comparing each coordinate with those of the extremes $p$ and $q$.
Notes
The algorithm is inspired from here.
LazySets.center
— Method.center(L::LineSegment{N})::Vector{N} where {N<:Real}
Return the center of a line segment.
Input
L
– line segment
Output
The center of the line segment.
LazySets.an_element
— Method.an_element(L::LineSegment{N}) where {N<:Real}
Return some element of a line segment.
Input
L
– line segment
Output
The first vertex of the line segment.
Base.rand
— Method.rand(::Type{LineSegment}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::LineSegment{N}
Create a random line segment.
Input
LineSegment
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random line segment.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1.
LazySets.halfspace_left
— Method.halfspace_left(L::LineSegment)
Return a half-space describing the 'left' of a two-dimensional line segment through two points.
Input
L
– line segment
Output
The half-space whose boundary goes through the two points p
and q
and which describes the left-hand side of the directed line segment pq
.
LazySets.halfspace_right
— Method.halfspace_right(L::LineSegment)
Return a half-space describing the 'right' of a two-dimensional line segment through two points.
Input
L
– line segment
Output
The half-space whose boundary goes through the two points p
and q
and which describes the right-hand side of the directed line segment pq
.
LazySets.vertices_list
— Method.vertices_list(L::LineSegment{N}
)::Vector{<:AbstractVector{N}} where {N<:Real}
Return the list of vertices of a line segment.
Input
L
– line segment
Output
The list of end points of the line segment.
LazySets.constraints_list
— Method.constraints_list(L::LineSegment{N})::Vector{LinearConstraint{N}} where {N<:Real}
Return the list of constraints defining a line segment in 2D.
Input
L
– line segment
Output
A vector of constraints that define the line segment.
Algorithm
$L$ is defined by 4 constraints. In this algorithm, the first two constraints are returned by $halfspace_right$ and $halfspace_left$, and the other two are obtained by considering the vector normal to the line segment that passes through each opposite vertex.
Notes
This function returns a vector of halfspaces. It does not return equality constraints.
LazySets.translate
— Method.translate(L::LineSegment{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a line segment by a given vector.
Input
L
– line segmentv
– translation vector
Output
A translated line segment.
Algorithm
We add the vector to both defining points of the line segment.
RecipesBase.apply_recipe
— Method.plot_linesegment(L::LineSegment; ...)
Plot a line segment.
Input
L
– line segment
Examples
julia> using Plots, LazySets;
julia> L = LineSegment([0., 0.], [1., 1.]);
julia> plot(L);
RecipesBase.apply_recipe
— Method.plot_linesegments(Xk::Vector{S}; ...) where {S<:LineSegment}
Plot an array of line segments.
Input
Xk
– linear array of line segments
Examples
julia> using Plots, LazySets;
julia> L1 = LineSegment([0., 0.], [1., 1.]);
julia> L2 = LineSegment([1., 0.], [0., 1.]);
julia> plot([L1, L2]);
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
:
Polygons
Constraint representation
LazySets.HPolygon
— Type.HPolygon{N<:Real} <: AbstractHPolygon{N}
Type that represents a convex polygon in constraint representation whose edges are sorted in counter-clockwise fashion with respect to their normal directions.
Fields
constraints
– list of linear constraints, sorted by the anglesort_constraints
– (optional, default:true
) flag for sorting the constraints (sortedness is a running assumption of this type)check_boundedness
– (optional, default:false
) flag for checking if the constraints make the polygon bounded; (boundedness is a running assumption of this type)
Notes
The default constructor assumes that the given list of edges is sorted. It does not perform any sorting. Use addconstraint!
to iteratively add the edges in a sorted way.
HPolygon(constraints::Vector{LinearConstraint{<:Real}})
– default constructorHPolygon()
– constructor with no constraints
LazySets.σ
— Method.σ(d::AbstractVector{N}, P::HPolygon{N};
[linear_search]::Bool=(length(P.constraints) < BINARY_SEARCH_THRESHOLD)
) where {N<:Real}
Return the support vector of a polygon in a given direction.
Input
d
– directionP
– polygon in constraint representationlinear_search
– (optional, default: see below) flag for controlling whether to perform a linear search or a binary search
Output
The support vector in the given direction. The result is always one of the vertices; in particular, if the direction has norm zero, any vertex is returned.
Algorithm
Comparison of directions is performed using polar angles; see the overload of <=
for two-dimensional vectors.
For polygons with BINARY_SEARCH_THRESHOLD = 10
or more constraints we use a binary search by default.
LazySets.translate
— Method.translate(v::AbstractVector{N}, P::HPolygon{N}; share::Bool=false
) where {N<:Real}
Translate (i.e., shift) a polygon in constraint representation by a given vector.
Input
P
– polygon in constraint representationv
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated polygon in constraint representation.
Notes
The normal vectors of the constraints (vector a
in a⋅x ≤ b
) are shared with the original constraints if share == true
.
Algorithm
We translate every constraint.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractPolygon
:
Inherited from AbstractHPolygon
:
an_element
∈
vertices_list
tohrep
tovrep
isbounded
addconstraint!
isredundant
remove_redundant_constraints!
constraints_list
Optimized constraint representation
LazySets.HPolygonOpt
— Type.HPolygonOpt{N<:Real} <: AbstractHPolygon{N}
Type that represents a convex polygon in constraint representation whose edges are sorted in counter-clockwise fashion with respect to their normal directions. This is a refined version of HPolygon
.
Fields
constraints
– list of linear constraintsind
– index in the list of constraints to begin the search to evaluate the support functionsort_constraints
– (optional, default:true
) flag for sorting the constraints (sortedness is a running assumption of this type)check_boundedness
– (optional, default:false
) flag for checking if the constraints make the polygon bounded; (boundedness is a running assumption of this type)
Notes
This structure is optimized to evaluate the support function/vector with a large sequence of directions that are close to each other. The strategy is to have an index that can be used to warm-start the search for optimal values in the support vector computation.
The default constructor assumes that the given list of edges is sorted. It does not perform any sorting. Use addconstraint!
to iteratively add the edges in a sorted way.
HPolygonOpt(constraints::Vector{LinearConstraint{<:Real}}, [ind]::Int=1)
– default constructor with optional index
LazySets.σ
— Method.σ(d::AbstractVector{N}, P::HPolygonOpt{N};
[linear_search]::Bool=(length(P.constraints) < BINARY_SEARCH_THRESHOLD)
) where {N<:Real}
Return the support vector of an optimized polygon in a given direction.
Input
d
– directionP
– optimized polygon in constraint representationlinear_search
– (optional, default: see below) flag for controlling whether to perform a linear search or a binary search
Output
The support vector in the given direction. The result is always one of the vertices; in particular, if the direction has norm zero, any vertex is returned.
Algorithm
Comparison of directions is performed using polar angles; see the overload of <=
for two-dimensional vectors.
For polygons with BINARY_SEARCH_THRESHOLD = 10
or more constraints we use a binary search by default.
LazySets.translate
— Method.translate(P::HPolygonOpt{N}, v::AbstractVector{N}; share::Bool=false
) where {N<:Real}
Translate (i.e., shift) an optimized polygon in constraint representation by a given vector.
Input
P
– optimized polygon in constraint representationv
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated optimized polygon in constraint representation.
Notes
The normal vectors of the constraints (vector a
in a⋅x ≤ b
) are shared with the original constraints if share == true
.
Algorithm
We translate every constraint.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractPolygon
:
Inherited from AbstractHPolygon
:
an_element
∈
vertices_list
tohrep
tovrep
isbounded
addconstraint!
isredundant
remove_redundant_constraints!
constraints_list
Vertex representation
LazySets.VPolygon
— Type.VPolygon{N<:Real} <: AbstractPolygon{N}
Type that represents a polygon by its vertices.
Fields
vertices
– the list of vertices
Notes
The constructor of VPolygon
runs a convex hull algorithm on its vertices by default, to remove the possibly redundant vertices. The vertices are sorted in counter-clockwise fashion. Use the flag apply_convex_hull=false
to skip the computation of the convex hull.
VPolygon(vertices::Vector{Vector{N}}; apply_convex_hull::Bool=true, algorithm::String="monotone_chain")
LazySets.σ
— Method.σ(d::AbstractVector{N}, P::VPolygon{N}) where {N<:Real}
Return the support vector of a polygon in a given direction.
Input
d
– directionP
– polygon in vertex representation
Output
The support vector in the given direction. If the direction has norm zero, the first vertex is returned.
Algorithm
This implementation performs a brute-force search, comparing the projection of each vector along the given direction. It runs in $O(n)$ where $n$ is the number of vertices.
Notes
For arbitrary points without structure this is the best one can do. However, a more efficient approach can be used if the vertices of the polygon have been sorted in counter-clockwise fashion. In that case a binary search algorithm can be used that runs in $O(\log n)$. See issue #40.
Base.:∈
— Method.∈(x::AbstractVector{N}, P::VPolygon{N})::Bool where {N<:Real}
Check whether a given point is contained in a polygon in vertex representation.
Input
x
– point/vectorP
– polygon in vertex representation
Output
true
iff $x ∈ P$.
Algorithm
This implementation exploits that the polygon's vertices are sorted in counter-clockwise fashion. Under this assumption we can just check if the vertex lies on the left of each edge, using the dot product.
Examples
julia> P = VPolygon([[2.0, 3.0], [3.0, 1.0], [5.0, 1.0], [4.0, 5.0]];
apply_convex_hull=false);
julia> ∈([4.5, 3.1], P)
false
julia> ∈([4.5, 3.0], P)
true
julia> ∈([4.4, 3.4], P) # point lies on the edge -> floating point error
false
julia> P = VPolygon([[2//1, 3//1], [3//1, 1//1], [5//1, 1//1], [4//1, 5//1]];
apply_convex_hull=false);
julia> ∈([44//10, 34//10], P) # with rational numbers the answer is correct
true
LazySets.an_element
— Method.an_element(P::VPolygon{N})::Vector{N} where {N<:Real}
Return some element of a polygon in vertex representation.
Input
P
– polygon in vertex representation
Output
The first vertex of the polygon in vertex representation.
Base.rand
— Method.rand(::Type{VPolygon}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::VPolygon{N}
Create a random polygon in vertex representation.
Input
VPolygon
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseedingnum_vertices
– (optional, default:-1
) number of vertices of the polygon (see comment below)
Output
A random polygon in vertex representation.
Algorithm
We follow the idea here based on P. Valtr. Probability that n random points are in convex position. There is also a nice video available here.
The number of vertices can be controlled with the argument num_vertices
. For a negative value we choose a random number in the range 3:10
.
LazySets.vertices_list
— Method.vertices_list(P::VPolygon{N})::Vector{Vector{N}} where {N<:Real}
Return the list of vertices of a convex polygon in vertex representation.
Input
P
– a polygon vertex representation
Output
List of vertices.
LazySets.tohrep
— Method.tohrep(P::VPolygon{N}, ::Type{HPOLYGON}=HPolygon
) where {N<:Real, HPOLYGON<:AbstractHPolygon}
Build a constraint representation of the given polygon.
Input
P
– polygon in vertex representationHPOLYGON
– (optional, default:HPolygon
) type of target polygon
Output
The same polygon but in constraint representation, an AbstractHPolygon
.
Algorithm
The algorithms consists of adding an edge for each consecutive pair of vertices. Since the vertices are already ordered in counter-clockwise fashion (CWW), the constraints will be sorted automatically (CCW) if we start with the first edge between the first and second vertex.
LazySets.tovrep
— Method.tovrep(P::VPolygon{N})::VPolygon{N} where {N<:Real}
Build a vertex representation of the given polygon.
Input
P
– polygon in vertex representation
Output
The identity, i.e., the same polygon instance.
LazySets.constraints_list
— Method.constraints_list(P::VPolygon{N})::Vector{LinearConstraint{N}} where {N<:Real}
Return the list of constraints defining a polygon in V-representation.
Input
P
– polygon in V-representation
Output
The list of constraints of the polygon.
Algorithm
First the H-representation of $P$ is computed, then its list of constraints is returned.
LazySets.translate
— Method.translate(P::VPolygon{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a polygon in vertex representation by a given vector.
Input
P
– polygon in vertex representationv
– translation vector
Output
A translated polygon in vertex representation.
Algorithm
We add the vector to each vertex of the polygon.
LazySets.remove_redundant_vertices
— Method.remove_redundant_vertices(P::VPolygon{N};
[algorithm]::String="monotone_chain")::VPolygon{N} where {N<:Real}
Return the polygon obtained by removing the redundant vertices of the given polygon.
Input
P
– polygon in vertex representationalgorithm
– (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 polygon.
Algorithm
A convex hull algorithm is used to compute the convex hull of the vertices of the given input polygon P
; see ?convex_hull
for details on the available algorithms. The vertices of the output polygon are sorted in counter-clockwise fashion.
LazySets.remove_redundant_vertices!
— Method.remove_redundant_vertices!(P::VPolygon{N};
[algorithm]::String="monotone_chain")::VPolygon{N} where {N<:Real}
Remove the redundant vertices of the given polygon.
Input
P
– polygon in vertex representationalgorithm
– (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 polygon.
Algorithm
A convex hull algorithm is used to compute the convex hull of the vertices of the given input polygon P
; see ?convex_hull
for details on the available algorithms. The vertices of the output polygon are sorted in counter-clockwise fashion.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractPolygon
:
Sorting directions
LazySets.jump2pi
— Function.jump2pi(x::N)::N where {N<:AbstractFloat}
Return $x + 2π$ if $x$ is negative, otherwise return $x$.
Input
x
– real scalar
Output
$x + 2π$ if $x$ is negative, $x$ otherwise.
Examples
julia> using LazySets: jump2pi
julia> jump2pi(0.0)
0.0
julia> jump2pi(-0.5)
5.783185307179586
julia> jump2pi(0.5)
0.5
Base.:<=
— Method.<=(u::AbstractVector{N}, v::AbstractVector{N})::Bool where {N<:AbstractFloat}
Compares two 2D vectors by their direction.
Input
u
– first 2D directionv
– second 2D direction
Output
True iff $\arg(u) [2π] ≤ \arg(v) [2π]$
Notes
The argument is measured in counter-clockwise fashion, with the 0 being the direction (1, 0).
Algorithm
The implementation uses the arctangent function with sign, atan
, which for two arguments implements the atan2
function.
Base.:<=
— Method.<=(u::AbstractVector{N}, v::AbstractVector{N})::Bool where {N<:Real}
Compare two 2D vectors by their direction.
Input
u
– first 2D directionv
– second 2D direction
Output
true
iff $\arg(u) [2π] ≤ \arg(v) [2π]$.
Notes
The argument is measured in counter-clockwise fashion, with the 0 being the direction (1, 0).
Algorithm
The implementation checks the quadrant of each direction, and compares directions using the right-hand rule (see is_right_turn
). In particular, this method does not use the arctangent.
LazySets.quadrant
— Method.quadrant(w::AbstractVector{N})::Int where {N<:Real}
Compute the quadrant where the direction w
belongs.
Input
w
– direction
Output
An integer from 0 to 3, with the following convention:
^
1 | 0
---+-->
2 | 3
Algorithm
The idea is to encode the following logic function: $11 ↦ 0, 01 ↦ 1, 00 ↦ 2, 10 ↦ 3$, according to the convention of above.
This function is inspired from AGPX's answer in: Sort points in clockwise order?
Polyhedra and Polytopes
Constraint representation
Convex polytopes are bounded polyhedra. The types HPolytope
and HPolyhedron
are used to represent polytopes and general polyhedra respectively, the difference being that for HPolytope
there is a running assumption about the boundedness of the set.
LazySets.HPolytope
— Type.HPolytope{N<:Real} <: AbstractPolytope{N}
Type that represents a convex polytope in H-representation.
Fields
constraints
– vector of linear constraintscheck_boundedness
– (optional, default:false
) flag for checking if the constraints make the polytope bounded; (boundedness is a running assumption of this type)
Note
Recall that a polytope is a bounded polyhedron. Boundedness is a running assumption in this type.
LazySets.HPolyhedron
— Type.HPolyhedron{N<:Real} <: AbstractPolyhedron{N}
Type that represents a convex polyhedron in H-representation.
Fields
constraints
– vector of linear constraints
The following methods are shared between HPolytope
and HPolyhedron
.
LazySets.dim
— Method.dim(P::HPoly{N})::Int where {N<:Real}
Return the dimension of a polyhedron in H-representation.
Input
P
– polyhedron in H-representation
Output
The ambient dimension of the polyhedron in H-representation. If it has no constraints, the result is $-1$.
LazySets.ρ
— Method.ρ(d::AbstractVector{N}, P::HPoly{N})::N where {N<:Real}
Evaluate the support function of a polyhedron (in H-representation) in a given direction.
Input
d
– directionP
– polyhedron in H-representation
Output
The support function of the polyhedron. If a polytope is unbounded in the given direction, we throw an error. If a polyhedron is unbounded in the given direction, the result is Inf
.
Algorithm
This implementation uses GLPKSolverLP
as linear programming backend.
LazySets.σ
— Method.σ(d::AbstractVector{N}, P::HPoly{N}) where {N<:Real}
Return the support vector of a polyhedron (in H-representation) in a given direction.
Input
d
– directionP
– polyhedron in H-representation
Output
The support vector in the given direction.
Algorithm
This implementation uses GLPKSolverLP
as linear programming backend.
LazySets.addconstraint!
— Method.addconstraint!(P::HPoly{N},
constraint::LinearConstraint{N})::Nothing where {N<:Real}
Add a linear constraint to a polyhedron in H-representation.
Input
P
– polyhedron in H-representationconstraint
– linear constraint to add
Output
Nothing.
Notes
It is left to the user to guarantee that the dimension of all linear constraints is the same.
LazySets.constraints_list
— Method.constraints_list(P::HPoly{N})::Vector{LinearConstraint{N}} where {N<:Real}
Return the list of constraints defining a polyhedron in H-representation.
Input
P
– polyhedron in H-representation
Output
The list of constraints of the polyhedron.
LazySets.tohrep
— Method.tohrep(P::HPoly{N}) where {N<:Real}
Return a constraint representation of the given polyhedron in constraint representation (no-op).
Input
P
– polyhedron in constraint representation
Output
The same polyhedron instance.
LazySets.tovrep
— Method.tovrep(P::HPoly{N};
[backend]=default_polyhedra_backend(P, N)) where {N<:Real}
Transform a polyhedron in H-representation to a polytope in V-representation.
Input
P
– polyhedron in constraint representationbackend
– (optional, default:default_polyhedra_backend(P, N)
) the backend for polyhedral computations
Output
The VPolytope
which is the vertex representation of the given polyhedron in constraint representation.
Notes
The conversion may not preserve the numeric type (e.g., with N == Float32
) depending on the backend. For further information on the supported backends see Polyhedra's documentation.
Base.isempty
— Method.isempty(P::HPoly{N}, witness::Bool=false; [usepolyhedrainterface]::Bool=false, [solver]=GLPKSolverLP(), [backend]=nothing )::Union{Bool, Tuple{Bool, Vector{N}}} where {N<:Real}
Determine whether a polyhedron is empty.
Input
P
– polyhedronwitness
– (optional, default:false
) compute a witness if activateduse_polyhedra_interface
– (optional, default:false
) iftrue
, we use thePolyhedra
interface for the emptiness testsolver
– (optional, default:GLPKSolverLP()
) LP-solver backendbackend
– (optional, default:nothing
) backend for polyhedral computations inPolyhedra
; its value is set internally (see the Notes below for details)
Output
- If
witness
option is deactivated:true
iff $P = ∅$ - If
witness
option is activated:(true, [])
iff $P = ∅$(false, v)
iff $P ≠ ∅$ and $v ∈ P$
Notes
The default value of the backend
is set internally and depends on whether the use_polyhedra_interface
option is set or not. If the option is set, we use default_polyhedra_backend(P, N)
.
Witness production is not supported if use_polyhedra_interface
is true
.
Algorithm
The algorithm sets up a feasibility LP for the constraints of P
. If use_polyhedra_interface
is true
, we call Polyhedra.isempty
. Otherwise, we set up the LP internally.
LazySets.translate
— Method.translate(P::PT, v::AbstractVector{N}; share::Bool=false
) where {N<:Real, PT<:HPoly{N}}
Translate (i.e., shift) a polyhedron in constraint representation by a given vector.
Input
P
– polyhedron in constraint representationv
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated polyhedron in constraint representation.
Notes
The normal vectors of the constraints (vector a
in a⋅x ≤ b
) are shared with the original constraints if share == true
.
Algorithm
We translate every constraint.
LazySets.cartesian_product
— Method.cartesian_product(P1::HPoly{N}, P2::HPoly{N};
[backend]=default_polyhedra_backend(P1, N)
) where {N<:Real}
Compute the Cartesian product of two polyhedra in H-representaion.
Input
P1
– polyhedronP2
– another polyhedronbackend
– (optional, default:default_polyhedra_backend(P1, N)
) the polyhedral computations backend
Output
The polyhedron obtained by the concrete cartesian product of P1
and P2
.
Notes
For further information on the supported backends see Polyhedra's documentation.
Polyhedra.polyhedron
— Method.polyhedron(P::HPoly{N};
[backend]=default_polyhedra_backend(P, N)) where {N<:Real}
Return an HRep
polyhedron from Polyhedra.jl
given a polytope in H-representation.
Input
P
– polytopebackend
– (optional, default: calldefault_polyhedra_backend(P, N)
) the polyhedral computations backend
Output
An HRep
polyhedron.
Notes
For further information on the supported backends see Polyhedra's documentation.
LazySets.remove_redundant_constraints
— Method.remove_redundant_constraints(P::PT;
backend=GLPKSolverLP()
)::Union{PT, EmptySet{N}} where {N<:Real,
PT<:HPoly{N}}
Remove the redundant constraints in a polyhedron in H-representation.
Input
P
– polyhedronbackend
– (optional, default:GLPKSolverLP
) the numeric LP solver backend
Output
A polyhedron equivalent to P
but with no redundant constraints, or an empty set if P
is detected to be empty, which may happen if the constraints are infeasible.
Algorithm
See remove_redundant_constraints!(::Vector{LinearConstraint{<:Real}})
for details.
LazySets.remove_redundant_constraints!
— Method.remove_redundant_constraints!(P::HPoly{N};
backend=GLPKSolverLP())::Bool where {N<:Real}
Remove the redundant constraints in a polyhedron in H-representation; the polyhedron is updated in-place.
Input
P
– polyhedronbackend
– (optional, default:GLPKSolverLP
) the numeric LP solver backend
Output
true
if the method was successful and the polyhedron P
is modified by removing its redundant constraints, and false
if P
is detected to be empty, which may happen if the constraints are infeasible.
Algorithm
See remove_redundant_constraints!(::Vector{LinearConstraint{<:Real}})
for details.
Inherited from LazySet
:
Inherited from AbstractPolyhedron
:
∈
- [
constrained_dimensions
](@ref constrained_dimensions(::AbstractPolyhedron) linear_map
Polytopes
The following methods are specific for HPolytope
.
Base.rand
— Method.rand(::Type{HPolytope}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::HPolytope{N}
Create a random polytope in constraint representation.
Input
HPolytope
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseedingnum_vertices
– (optional, default:-1
) upper bound on the number of vertices of the polytope (see comment below)
Output
A random polytope in constraint representation.
Algorithm
We create a random polytope in vertex representation and convert it to constraint representation (hence the argument num_vertices
). See rand(::Type{VPolytope})
.
LazySets.vertices_list
— Method.vertices_list(P::HPolytope{N};
[backend]=default_polyhedra_backend(P, N),
[prunefunc]=removevredundancy!)::Vector{Vector{N}} where
{N<:Real}
Return the list of vertices of a polytope in constraint representation.
Input
P
– polytope in constraint representationbackend
– (optional, default:default_polyhedra_backend(P, N)
) the polyhedral computations backendprunefunc
– (optional, default:removevredundancy!
) function to post-process the output ofvreps
Output
List of vertices.
Notes
For further information on the supported backends see Polyhedra's documentation.
LazySets.isbounded
— Function.isbounded(P::HPolytope, [use_type_assumption]::Bool=true)::Bool
Determine whether a polytope in constraint representation is bounded.
Input
P
– polytope in constraint representationuse_type_assumption
– (optional, default:true
) flag for ignoring the type assumption that polytopes are bounded
Output
true
if use_type_assumption
is activated. Otherwise, true
iff P
is bounded.
Algorithm
If !use_type_assumption
, we convert P
to an HPolyhedron
P2
and then use isbounded(P2)
.
Inherited from AbstractPolytope
:
Polyhedra
The following methods are specific for HPolyhedron
.
Base.rand
— Method.rand(::Type{HPolyhedron}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::HPolyhedron{N}
Create a polyhedron.
Input
HPolyhedron
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimension (is ignored)rng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A polyhedron.
Algorithm
We first create a random polytope and then randomly remove some of the constraints.
LazySets.isbounded
— Method.isbounded(P::HPolyhedron)::Bool
Determine whether a polyhedron in constraint representation is bounded.
Input
P
– polyhedron in constraint representation
Output
true
iff the polyhedron is bounded.
Algorithm
We first check if the polyhedron has more than max(dim(P), 1)
constraints, which is a necessary condition for boundedness. If so, we check boundedness via isbounded_unit_dimensions
.
LazySets.vertices_list
— Method.vertices_list(P::HPolyhedron{N}) where {N<:Real}
Return the list of vertices of a polyhedron in constraint representation.
Input
P
– polyhedron in constraint representation
Output
This function returns an error because the polyhedron is possibly unbounded. If P
is known to be bounded, try converting to HPolytope
first:
julia> P = HPolyhedron([HalfSpace([1.0, 0.0], 1.0),
HalfSpace([0.0, 1.0], 1.0),
HalfSpace([-1.0, 0.0], 1.0),
HalfSpace([0.0, -1.0], 1.0)]);
julia> P_as_polytope = convert(HPolytope, P);
LazySets.singleton_list
— Method.singleton_list(P::HPolyhedron{N})::Vector{Singleton{N}} where {N<:Real}
Return the vertices of a polyhedron in H-representation as a list of singletons.
Input
P
– polytope in constraint representation
Output
This function returns an error because the polyhedron is possibly unbounded. If P
is known to be bounded, try converting to HPolytope
first.
Vertex representation
LazySets.VPolytope
— Type.VPolytope{N<:Real} <: AbstractPolytope{N}
Type that represents a convex polytope in V-representation.
Fields
vertices
– the list of vertices
LazySets.dim
— Method.dim(P::VPolytope)::Int
Return the dimension of a polytope in V-representation.
Input
P
– polytope in V-representation
Output
The ambient dimension of the polytope in V-representation. If it is empty, the result is $-1$.
Examples
julia> v = VPolytope();
julia> dim(v) > 0
false
julia> v = VPolytope([ones(3)])
VPolytope{Float64}(Array{Float64,1}[[1.0, 1.0, 1.0]])
julia> dim(v) == 3
true
LazySets.σ
— Method.σ(d::AbstractVector{N}, P::VPolytope{N}; algorithm="hrep") where {N<:Real}
Return the support vector of a polyhedron (in V-representation) in a given direction.
Input
d
– directionP
– polyhedron in V-representationalgorithm
– (optional, default:'hrep'
) method to compute the support vector
Output
The support vector in the given direction.
Base.:∈
— Method.∈(x::AbstractVector{N}, P::VPolytope{N};
solver=GLPKSolverLP(method=:Simplex))::Bool where {N<:Real}
Check whether a given point is contained in a polytope in vertex representation.
Input
x
– point/vectorP
– polytope in vertex representation
Output
true
iff $x ∈ P$.
Algorithm
We check, using linear programming, the definition of a convex polytope that a point is in the set if and only if it is a convex combination of the vertices.
Let $\{v_j\}$ be the $m$ vertices of P
. Then we solve the following $m$-dimensional linear program.
Base.rand
— Method.rand(::Type{VPolytope}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::VPolytope{N}
Create a random polytope in vertex representation.
Input
VPolytope
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseedingnum_vertices
– (optional, default:-1
) upper bound on the number of vertices of the polytope (see comment below)
Output
A random polytope in vertex representation.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1.
The number of vertices can be controlled with the argument num_vertices
. For a negative value we choose a random number in the range dim:5*dim
(except if dim == 1
, in which case we choose in the range 1:2
). Note that we do not guarantee that the vertices are not redundant.
LazySets.translate
— Method.translate(P::VPolytope{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a polytope in vertex representation by a given vector.
Input
P
– polytope in vertex representationv
– translation vector
Output
A translated polytope in vertex representation.
Algorithm
We add the vector to each vertex of the polytope.
LazySets.vertices_list
— Method.vertices_list(P::VPolytope{N})::Vector{Vector{N}} where {N<:Real}
Return the list of vertices of a polytope in V-representation.
Input
P
– polytope in vertex representation
Output
List of vertices.
LazySets.remove_redundant_vertices
— Method.remove_redundant_vertices(P::VPolytope{N};
[backend]=default_polyhedra_backend(P, N))::VPolytope{N} where {N<:Real}
Return the polytope obtained by removing the redundant vertices of the given polytope.
Input
P
– polytope in vertex representationbackend
– (optional, default:default_polyhedra_backend(P1, N)
) the polyhedral computations backend, see Polyhedra's documentation for further information
Output
A new polytope such that its vertices are the convex hull of the given polytope.
LazySets.constraints_list
— Method.constraints_list(P::VPolytope{N})::Vector{LinearConstraint{N}} where {N<:Real}
Return the list of constraints defining a polytope in V-representation.
Input
P
– polytope in V-representation
Output
The list of constraints of the polytope.
Algorithm
First the H-representation of $P$ is computed, then its list of constraints is returned.
LazySets.tohrep
— Method.tohrep(P::VPolytope{N};
[backend]=default_polyhedra_backend(P, N)) where {N<:Real}
Transform a polytope in V-representation to a polytope in H-representation.
Input
P
– polytope in vertex representationbackend
– (optional, default:default_polyhedra_backend(P, N)
) the backend for polyhedral computations
Output
The HPolytope
which is the constraint representation of the given polytope in vertex representation.
Notes
The conversion may not preserve the numeric type (e.g., with N == Float32
) depending on the backend. For further information on the supported backends see Polyhedra's documentation.
LazySets.tovrep
— Method.tovrep(P::VPolytope)
Return a vertex representation of the given polytope in vertex representation (no-op).
Input
P
– polytope in vertex representation
Output
The same polytope instance.
LazySets.cartesian_product
— Method.cartesian_product(P1::VPolytope{N}, P2::VPolytope{N};
[backend]=default_polyhedra_backend(P1, N)) where {N}
Compute the Cartesian product of two polytopes in V-representation.
Input
P1
– polytopeP2
– another polytopebackend
– (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 Cartesian product of P1
and P2
.
Polyhedra.polyhedron
— Method.polyhedron(P::VPolytope{N};
[backend]=default_polyhedra_backend(P, N)) where {N<:Real}
Return an VRep
polyhedron from Polyhedra.jl
given a polytope in V-representation.
Input
P
– polytopebackend
– (optional, default:default_polyhedra_backend(P, N)
) the polyhedral computations backend, see Polyhedra's documentation for further information
Output
A VRep
polyhedron.
LazySets.linear_map
— Method.linear_map(M::AbstractMatrix{N}, P::VPolytope{N}) where {N<:Real}
Concrete linear map of a polytope in vertex representation.
Input
M
– matrixP
– polytope in vertex representation
Output
A polytope in vertex representation.
Algorithm
The linear map $M$ is applied to each vertex of the given set $P$, obtaining a polytope in V-representation. The output type is again a VPolytope
.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Singleton
LazySets.Singleton
— Type.Singleton{N<:Real, VN<:AbstractVector{N}} <: AbstractSingleton{N}
Type that represents a singleton, that is, a set with a unique element.
Fields
element
– the only element of the set
Base.rand
— Method.rand(::Type{Singleton}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Singleton{N}
Create a random singleton.
Input
Singleton
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
A random singleton.
Algorithm
The element is a normally distributed vector with entries of mean 0 and standard deviation 1.
LazySets.element
— Method.element(S::Singleton{N}) where {N<:Real}
Return the element of a singleton.
Input
S
– singleton
Output
The element of the singleton.
LazySets.element
— Method.element(S::Singleton{N}, i::Int)::N where {N<:Real}
Return the i-th entry of the element of a singleton.
Input
S
– singletoni
– dimension
Output
The i-th entry of the element of the singleton.
LazySets.translate
— Method.translate(S::Singleton{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a singleton by a given vector.
Input
S
– singletonv
– translation vector
Output
A translated singleton.
Algorithm
We add the vector to the point in the singleton.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
:
Inherited from AbstractHyperrectangle
:
Inherited from AbstractSingleton
:
Universe
LazySets.Universe
— Type.Universe{N<:Real} <: AbstractPolyhedron{N}
Type that represents the universal set, i.e., the set of all elements.
LazySets.dim
— Method.dim(U::Universe)
Return the dimension of a universe.
Input
U
– universe
Output
The dimension of a universe.
LazySets.ρ
— Method.ρ(d::AbstractVector{N}, U::Universe{N}) where {N<:Real}
Return the support function of a universe.
Input
d
– directionU
– universe
Output
The support function in the given direction.
Algorithm
If the direction is all zero, the result is zero. Otherwise, the result is Inf
.
LazySets.σ
— Method.σ(d::AbstractVector{N}, U::Universe{N}) where {N<:Real}
Return the support vector of a universe.
Input
d
– directionU
– universe
Output
A vector with infinity values, except in dimensions where the direction is zero.
Base.:∈
— Method.∈(x::AbstractVector{N}, U::Universe{N})::Bool where {N<:Real}
Check whether a given point is contained in a universe.
Input
x
– point/vectorU
– universe
Output
The output is always true
.
Examples
julia> ∈([1.0, 0.0], Universe(2))
true
Base.rand
— Method.rand(::Type{Universe}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Universe{N}
Create a universe (note that there is nothing to randomize).
Input
Universe
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
The (only) universe of the given numeric type and dimension.
LazySets.an_element
— Method.an_element(U::Universe{N}) where {N<:Real}
Return some element of a universe.
Input
U
– universe
Output
The origin.
Base.isempty
— Method.isempty(U::Universe)::Bool
Return if a universe is empty or not.
Input
U
– universe
Output
false
.
LazySets.isbounded
— Method.isbounded(U::Universe)::Bool
Determine whether a universe is bounded.
Input
S
– universe
Output
false
as the universe is unbounded.
LinearAlgebra.norm
— Function.norm(U::Universe, [p]::Real=Inf)
Return the norm of a universe. It is the norm of the enclosing ball (of the given $p$-norm) of minimal volume that is centered in the origin.
Input
U
– universep
– (optional, default:Inf
) norm
Output
An error.
LazySets.radius
— Function.radius(U::Universe, [p]::Real=Inf)
Return the radius of a universe. It is the radius of the enclosing ball (of the given $p$-norm) of minimal volume with the same center.
Input
U
– universep
– (optional, default:Inf
) norm
Output
An error.
LazySets.diameter
— Function.diameter(U::Universe, [p]::Real=Inf)
Return the diameter of a universe. It is the maximum distance between any two elements of the set, or, equivalently, the diameter of the enclosing ball (of the given $p$-norm) of minimal volume with the same center.
Input
U
– universep
– (optional, default:Inf
) norm
Output
An error.
LazySets.constraints_list
— Method.constraints_list(U::Universe{N})::Vector{LinearConstraint{N}}
where {N<:Real}
Return the list of constraints defining a universe.
Input
U
– universe
Output
The empty list of constraints, as the universe is unconstrained.
LazySets.constrained_dimensions
— Method.constrained_dimensions(U::Universe)::Vector{Int}
Return the indices in which a universe is constrained.
Input
U
– universe
Output
The empty vector, as the universe is unconstrained in every dimension.
LazySets.translate
— Method.translate(U::Universe{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a universe by a given vector.
Input
U
– universev
– translation vector
Output
The universe.
Zero set
LazySets.ZeroSet
— Type.ZeroSet{N<:Real} <: AbstractSingleton{N}
Type that represents the zero set, i.e., the set that only contains the origin.
Fields
dim
– the ambient dimension of this zero set
LazySets.dim
— Method.dim(Z::ZeroSet)::Int
Return the ambient dimension of this zero set.
Input
Z
– a zero set, i.e., a set that only contains the origin
Output
The ambient dimension of the zero set.
LazySets.σ
— Method.σ(d::AbstractVector{N}, Z::ZeroSet{N}) where {N<:Real}
Return the support vector of a zero set.
Input
Z
– a zero set, i.e., a set that only contains the origin
Output
The returned value is the origin since it is the only point that belongs to this set.
Base.:∈
— Method.∈(x::AbstractVector{N}, Z::ZeroSet{N})::Bool where {N<:Real}
Check whether a given point is contained in a zero set.
Input
x
– point/vectorZ
– zero set
Output
true
iff $x ∈ Z$.
Examples
julia> Z = ZeroSet(2);
julia> ∈([1.0, 0.0], Z)
false
julia> ∈([0.0, 0.0], Z)
true
Base.rand
— Method.rand(::Type{ZeroSet}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::ZeroSet{N}
Create a zero set (note that there is nothing to randomize).
Input
ZeroSet
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseeding
Output
The (only) zero set of the given numeric type and dimension.
LazySets.element
— Method.element(S::ZeroSet{N})::Vector{N} where {N<:Real}
Return the element of a zero set.
Input
S
– zero set
Output
The element of the zero set, i.e., a zero vector.
LazySets.element
— Method.element(S::ZeroSet{N}, ::Int)::N where {N<:Real}
Return the i-th entry of the element of a zero set.
Input
S
– zero seti
– dimension
Output
The i-th entry of the element of the zero set, i.e., 0.
LazySets.linear_map
— Method.linear_map(M::AbstractMatrix{N}, Z::ZeroSet{N}) where {N<:Real}
Concrete linear map of a zero set.
Input
M
– matrixZ
– zero set
Output
The zero set whose dimension matches the output dimension of the given matrix.
LazySets.translate
— Method.translate(Z::ZeroSet{N}, v::AbstractVector{N}) where {N<:Real}
Translate (i.e., shift) a zero set by a given vector.
Input
Z
– zero setv
– translation vector
Output
A singleton containing the vector v
.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
:
Inherited from AbstractHyperrectangle
:
Inherited from AbstractSingleton
:
Zonotope
LazySets.Zonotope
— Type.Zonotope{N<:Real} <: AbstractCentrallySymmetricPolytope{N}
Type that represents a zonotope.
Fields
center
– center of the zonotopegenerators
– matrix; each column is a generator of the zonotope
Notes
Mathematically, a zonotope is defined as the set
where $c \in \mathbb{R}^n$ is its center and $\{g_i\}_{i=1}^p$, $g_i \in \mathbb{R}^n$, is the set of generators. This characterization defines a zonotope as the finite Minkowski sum of line segments. Zonotopes can be equivalently described as the image of a unit infinity-norm ball in $\mathbb{R}^n$ by an affine transformation.
Zonotope(center::AbstractVector{N}, generators::AbstractMatrix{N}) where {N<:Real}
Zonotope(center::AbstractVector{N}, generators_list::AbstractVector{VN} ) where {N<:Real, VN<:AbstractVector{N}}
The optional argument remove_zero_generators
controls whether we remove zero columns from the generators
matrix. This option is active by default.
Examples
A two-dimensional zonotope with given center and set of generators:
julia> Z = Zonotope([1.0, 0.0], [0.1 0.0; 0.0 0.1])
Zonotope{Float64}([1.0, 0.0], [0.1 0.0; 0.0 0.1])
julia> dim(Z)
2
Compute its vertices:
julia> vertices_list(Z)
4-element Array{Array{Float64,1},1}:
[1.1, 0.1]
[0.9, 0.1]
[1.1, -0.1]
[0.9, -0.1]
Evaluate the support vector in a given direction:
julia> σ([1., 1.], Z)
2-element Array{Float64,1}:
1.1
0.1
Alternative constructor: A zonotope in two dimensions with three generators:
julia> Z = Zonotope(ones(2), [[1., 0.], [0., 1.], [1., 1.]])
Zonotope{Float64}([1.0, 1.0], [1.0 0.0 1.0; 0.0 1.0 1.0])
julia> Z.generators
2×3 Array{Float64,2}:
1.0 0.0 1.0
0.0 1.0 1.0
LazySets.ρ
— Method.ρ(d::AbstractVector{N}, Z::Zonotope{N}) where {N<:Real}
Return the support function of a zonotope in a given direction.
Input
d
– directionZ
– zonotope
Output
The support function of the zonotope in the given direction.
Algorithm
The support value is $cᵀ d + ‖Gᵀ d‖₁$ where $c$ is the center and $G$ is the generator matrix of Z
.
LazySets.σ
— Method.σ(d::AbstractVector{N}, Z::Zonotope{N}) where {N<:Real}
Return the support vector of a zonotope in a given direction.
Input
d
– directionZ
– zonotope
Output
Support vector in the given direction. If the direction has norm zero, the vertex with $ξ_i = 1 \ \ ∀ i = 1,…, p$ is returned.
Base.:∈
— Method.∈(x::AbstractVector{N}, Z::Zonotope{N};
solver=GLPKSolverLP(method=:Simplex))::Bool where {N<:Real}
Check whether a given point is contained in a zonotope.
Input
x
– point/vectorZ
– zonotopesolver
– (optional, default:GLPKSolverLP(method=:Simplex)
) the backend used to solve the linear program
Output
true
iff $x ∈ Z$.
Examples
julia> Z = Zonotope([1.0, 0.0], [0.1 0.0; 0.0 0.1]);
julia> ∈([1.0, 0.2], Z)
false
julia> ∈([1.0, 0.1], Z)
true
Algorithm
The membership problem is computed by stating and solving the following linear program with the simplex method. Let $p$ and $n$ be the number of generators and ambient dimension, respectively. We consider the minimization of $x_0$ in the $p+1$-dimensional space of elements $(x_0, ξ_1, …, ξ_p)$ constrained to $0 ≤ x_0 ≤ ∞$, $ξ_i ∈ [-1, 1]$ for all $i = 1, …, p$, and such that $x-c = Gξ$ holds. If a feasible solution exists, the optimal value $x_0 = 0$ is achieved.
Notes
This function is parametric in the number type N
. For exact arithmetic use an appropriate backend, e.g. solver=GLPKSolverLP(method=:Exact)
.
Base.rand
— Method.rand(::Type{Zonotope}; [N]::Type{<:Real}=Float64, [dim]::Int=2,
[rng]::AbstractRNG=GLOBAL_RNG, [seed]::Union{Int, Nothing}=nothing
)::Zonotope{N}
Create a random zonotope.
Input
Zonotope
– type for dispatchN
– (optional, default:Float64
) numeric typedim
– (optional, default: 2) dimensionrng
– (optional, default:GLOBAL_RNG
) random number generatorseed
– (optional, default:nothing
) seed for reseedingnum_generators
– (optional, default:-1
) number of generators of the zonotope (see comment below)
Output
A random zonotope.
Algorithm
All numbers are normally distributed with mean 0 and standard deviation 1.
The number of generators can be controlled with the argument num_generators
. For a negative value we choose a random number in the range dim:2*dim
(except if dim == 1
, in which case we only create a single generator).
LazySets.vertices_list
— Method.vertices_list(Z::Zonotope{N})::Vector{Vector{N}} where {N<:Real}
Return the vertices of a zonotope.
Input
Z
– zonotope
Output
List of vertices as a vector of vectors.
Algorithm
If the zonotope has $p$ generators, each of the $2^p$ vertices is computed by taking the sum of the center and a linear combination of generators, where the combination factors are $ξ_i ∈ \{-1, 1\}$.
Notes
For high dimensions, it would be preferable to develop a vertex_iterator
approach.
LazySets.constraints_list
— Method.constraints_list(P::Zonotope{N}
)::Vector{LinearConstraint{N}} where {N<:Real}
Return the list of constraints defining a zonotope.
Input
Z
– zonotope
Output
The list of constraints of the zonotope.
Algorithm
This is the (inefficient) fallback implementation for rational numbers. It first computes the vertices and then converts the corresponding polytope to constraint representation.
LazySets.constraints_list
— Method.constraints_list(Z::Zonotope{N}
)::Vector{LinearConstraint{N}} where {N<:AbstractFloat}
Return the list of constraints defining a zonotope.
Input
Z
– zonotope
Output
The list of constraints of the zonotope.
Notes
The algorithm assumes that no generator is redundant. The result has $2 \binom{p}{n-1}$ (with $p$ being the number of generators and $n$ being the ambient dimension) constraints, which is optimal under this assumption.
If $p < n$, we fall back to the (slower) computation based on the vertex representation.
Algorithm
We follow the algorithm presented in Althoff, Stursberg, Buss: Computing Reachable Sets of Hybrid Systems Using a Combination of Zonotopes and Polytopes. 2009.
The one-dimensional case is not covered by that algorithm; we manually handle this case, assuming that there is only one generator.
LazySets.center
— Method.center(Z::Zonotope{N})::Vector{N} where {N<:Real}
Return the center of a zonotope.
Input
Z
– zonotope
Output
The center of the zonotope.
LazySets.order
— Method.order(Z::Zonotope)::Rational
Return the order of a zonotope.
Input
Z
– zonotope
Output
A rational number representing the order of the zonotope.
Notes
The order of a zonotope is defined as the quotient of its number of generators and its dimension.
LazySets.minkowski_sum
— Method.minkowski_sum(Z1::Zonotope{N}, Z2::Zonotope{N}) where {N<:Real}
Concrete Minkowski sum of a pair of zonotopes.
Input
Z1
– one zonotopeZ2
– another zonotope
Output
The zonotope obtained by summing the centers and concatenating the generators of $Z_1$ and $Z_2$.
LazySets.linear_map
— Method.linear_map(M::AbstractMatrix{N}, Z::Zonotope{N}) where {N<:Real}
Concrete linear map of a zonotope.
Input
M
– matrixZ
– zonotope
Output
The zonotope obtained by applying the linear map to the center and generators of $Z$.
LazySets.translate
— Method.translate(Z::Zonotope{N}, v::AbstractVector{N}; share::Bool=false
) where {N<:Real}
Translate (i.e., shift) a zonotope by a given vector.
Input
Z
– zonotopev
– translation vectorshare
– (optional, default:false
) flag for sharing unmodified parts of the original set representation
Output
A translated zonotope.
Notes
The generator matrix is shared with the original zonotope if share == true
.
Algorithm
We add the vector to the center of the zonotope.
LazySets.scale
— Method.scale(α::Real, Z::Zonotope)
Concrete scaling of a zonotope.
Input
α
– scalarZ
– zonotope
Output
The zonotope obtained by applying the numerical scale to the center and generators of $Z$.
LazySets.ngens
— Method.ngens(Z::Zonotope)::Int
Return the number of generators of a zonotope.
Input
Z
– zonotope
Output
Integer representing the number of generators.
LazySets.reduce_order
— Method.reduce_order(Z::Zonotope, r)::Zonotope
Reduce the order of a zonotope by overapproximating with a zonotope with less generators.
Input
Z
– zonotoper
– desired order
Output
A new zonotope with less generators, if possible.
Algorithm
This function implements the algorithm described in A. Girard's Reachability of Uncertain Linear Systems Using Zonotopes, HSCC. Vol. 5. 2005.
If the desired order is smaller than one, the zonotope is not reduced.
Base.split
— Method.split(Z::Zonotope, j::Int)
Return two zonotopes obtained by splitting the given zonotope.
Input
Z
– zonotopej
– index of the generator to be split
Output
The zonotope obtained by splitting Z
into two zonotopes such that their union is Z
and their intersection is possibly non-empty.
Algorithm
This function implements [Prop. 3, 1], that we state next. The zonotope $Z = ⟨c, g^{(1, …, p)}⟩$ is split into:
such that $Z₁ ∪ Z₂ = Z$ and $Z₁ ∩ Z₂ = Z^*$, where
[1] Althoff, M., Stursberg, O., & Buss, M. (2008). Reachability analysis of nonlinear systems with uncertain parameters using conservative linearization. In Proc. of the 47th IEEE Conference on Decision and Control.
Inherited from LazySet
:
Inherited from AbstractPolytope
:
Inherited from AbstractCentrallySymmetricPolytope
: