Utility functions
Helpers for internal use only
Functions and Macros
LazySets.an_element_helper — Function.an_element_helper(hp::Hyperplane{N},
[nonzero_entry_a]::Int)::Vector{N} where {N<:Real}Helper function that computes an element on a hyperplane's hyperplane.
Input
hp– hyperplanenonzero_entry_a– (optional, default: computes the first index) indexisuch thathp.a[i]is different from 0
Output
An element on a hyperplane.
Algorithm
We compute the point on the hyperplane as follows:
- We already found a nonzero entry of $a$ in dimension, say, $i$.
- We set $x[i] = b / a[i]$.
- We set $x[j] = 0$ for all $j ≠ i$.
LazySets.binary_search_constraints — Function.binary_search_constraints(d::AbstractVector{N},
constraints::Vector{LinearConstraint{N}},
n::Int,
k::Int;
[choose_lower]::Bool=false
)::Int where {N<:Real}Performs a binary search in the constraints.
Input
d– directionconstraints– constraintsn– number of constraintsk– start indexchoose_lower– (optional, default:false) flag for choosing the lower index (see the 'Output' section)
Output
In the default setting, the result is the smallest index k such that d <= constraints[k], or n+1 if no such k exists. If the choose_lower flag is set, the result is the largest index k such that constraints[k] < d, which is equivalent to being k-1 in the normal setting.
LazySets.cross_product — Method.cross_product(M::AbstractMatrix{N}) where {N<:Real}Compute the high-dimensional cross product of $n-1$ $n$-dimensional vectors.
Input
M– $n × n - 1$-dimensional matrix
Output
A vector.
Algorithm
The cross product is defined as follows:
where $M^{[i]}$ is defined as $M$ with the $i$-th row removed. See Althoff, Stursberg, Buss: Computing Reachable Sets of Hybrid Systems Using a Combination of Zonotopes and Polytopes. 2009.
LazySets.delete_zero_columns — Function.delete_zero_columns(A::AbstractMatrix)Remove all columns that only contain zeros from a given matrix.
Input
A– matrixcopy– (optional, default:false) flag to copy the matrix
Output
A matrix.
If the input matrix A does not contain any zero column, we return A unless the option copy is set. If the input matrix contains zero columns, we always return a copy if the option copy is set and otherwise a SubArray via @view.
LazySets.dot_zero — Function.dot_zero(x::AbstractVector{N}, y::AbstractVector{N}) where{N<:Real}Dot product with preference for zero value in the presence of infinity values.
Input
x– first vectory– second vector
Output
The dot product of x and y, but with the rule that 0 * Inf == 0.
LazySets.get_radius! — Function.get_radius!(sih::SymmetricIntervalHull{N},
i::Int,
n::Int=dim(sih))::N where {N<:Real}Compute the radius of a symmetric interval hull of a convex set in a given dimension.
Input
sih– symmetric interval hull of a convex seti– dimension in which the radius should be computedn– (optional, default:dim(sih)) set dimension
Output
The radius of a symmetric interval hull of a convex set in a given dimension.
Algorithm
We ask for the support vector of the underlying set for both the positive and negative unit vector in the dimension i.
LazySets.inner — Function.inner(x::AbstractVector{N}, A::AbstractMatrix{N}, y::AbstractVector{N}
) where {N}Compute the inner product $xᵀ A y$.
Input
x– vector on the leftA– matrixy– vector on the right
Output
The (scalar) result of the multiplication.
LazySets.isinvertible — Function.isinvertible(M::Matrix; [cond_tol]::Number=DEFAULT_COND_TOL)A sufficient check of a matrix being invertible (or nonsingular).
Input
M– matrixcond_tol– (optional, default:DEFAULT_COND_TOL) tolerance of matrix condition
Output
If the result is true, M is invertible. If the result is false, the matrix is non-square or this function could not conclude.
Algorithm
We check whether the matrix is square and whether the matrix condition number cond(M) is below some prescribed tolerance.
LazySets.ispermutation — Function.ispermutation(u::AbstractVector{T}, v::AbstractVector{T})::Bool where TCheck that two vectors contain the same elements up to reordering.
Input
u– first vectorv– second vector
Output
true iff the vectors are identical up to reordering.
Examples
julia> LazySets.ispermutation([1, 2, 2], [2, 2, 1])
true
julia> LazySets.ispermutation([1, 2, 2], [1, 1, 2])
false
LazySets.issquare — Function.issquare(M::AbstractMatrix)::BoolCheck whether a matrix is square.
Input
M– matrix
Output
true iff the matrix is square.
LazySets.is_right_turn — Function.is_right_turn(u::AbstractVector{N},
v::AbstractVector{N})::Bool where {N<:Real}Determine if the acute angle defined by two 2D vectors is a right turn (< 180° counter-clockwise).
Input
u– first 2D directionv– second 2D direction
Output
true iff the two vectors constitute a right turn.
LazySets.is_tighter_same_dir_2D — Function.is_tighter_same_dir_2D(c1::LinearConstraint{N},
c2::LinearConstraint{N}) where {N<:Real}Check if the first of two two-dimensional constraints with equivalent normal direction is tighter.
Input
c1– first linear constraintc2– second linear constraintstrict– (optional; default:false) check for strictly tighter constraints?
Output
true iff the first constraint is tighter.
LazySets.nonzero_indices — Function.nonzero_indices(v::AbstractVector{N})::Vector{Int} where {N<:Real}Return the indices in which a vector is non-zero.
Input
v– vector
Output
A vector of ascending indices i such that the vector is non-zero in dimension i.
LazySets.samedir — Function.samedir(u::AbstractVector{N},
v::AbstractVector{N})::Tuple{Bool, Real} where {N<:Real}Check whether two vectors point in the same direction.
Input
u– first vectorv– second vector
Output
(true, k) iff the vectors are identical up to a positive scaling factor k, and (false, 0) otherwise.
Examples
julia> using LazySets: samedir
julia> samedir([1, 2, 3], [2, 4, 6])
(true, 0.5)
julia> samedir([1, 2, 3], [3, 2, 1])
(false, 0)
julia> samedir([1, 2, 3], [-1, -2, -3])
(false, 0)
LazySets.sign_cadlag — Function.sign_cadlag(x::N)::N where {N<:Real}This function works like the sign function but is $1$ for input $0$.
Input
x– real scalar
Output
$1$ if $x ≥ 0$, $-1$ otherwise.
Notes
This is the sign function right-continuous at zero (see càdlàg function). It can be used with vector-valued arguments via the dot operator.
Examples
julia> LazySets.sign_cadlag.([-0.6, 1.3, 0.0])
3-element Array{Float64,1}:
-1.0
1.0
1.0LazySets._random_zero_sum_vector — Function._random_zero_sum_vector(rng::AbstractRNG, N::Type{<:Real}, n::Int)Create a random vector with entries whose sum is zero.
Input
rng– random number generatorN– numeric typen– length of vector
Output
A random vector of random numbers such that all positive entries come first and all negative entries come last, and such that the total sum is zero.
Algorithm
This is a preprocessing step of the algorithm here based on P. Valtr. Probability that n random points are in convex position.
LazySets.remove_duplicates_sorted! — Function.remove_duplicates_sorted!(v::AbstractVector)Remove duplicate entries in a sorted vector.
Input
v– sorted vector
Output
The input vector without duplicates.
LazySets.reseed — Function.reseed(rng::AbstractRNG, seed::Union{Int, Nothing})::AbstractRNGReset the RNG seed if the seed argument is a number.
Input
rng– random number generatorseed– seed for reseeding
Output
The input RNG if the seed is nothing, and a reseeded RNG otherwise.
LazySets.same_block_structure — Function.same_block_structure(x::AbstractVector{S1}, y::AbstractVector{S2}
)::Bool where {S1<:LazySet, S2<:LazySet}Check whether two vectors of sets have the same block structure, i.e., the $i$-th entry in the vectors have the same dimension.
Input
x– first vectory– second vector
Output
true iff the vectors have the same block structure.
LazySets.substitute — Function.substitute(substitution::Dict{Int, T}, x::AbstractVector{T}) where {T}Apply a substitution to a given vector.
Input
substitution– substitution (a mapping from an index to a new value)x– vector
Output
A fresh vector corresponding to x after substitution was applied.
LazySets.substitute! — Function.substitute!(substitution::Dict{Int, T}, x::AbstractVector{T}) where {T}Apply a substitution to a given vector.
Input
substitution– substitution (a mapping from an index to a new value)x– vector (modified in this function)
Output
The same (but see the Notes below) vector x but after substitution was applied.
Notes
The vector x is modified in-place if it has type Vector or SparseVector. Otherwise, we first create a new Vector from it.
LazySets.σ_helper — Function. σ_helper(d::AbstractVector{N},
hp::Hyperplane{N};
error_unbounded::Bool=true,
[halfspace]::Bool=false) where {N<:Real}Return the support vector of a hyperplane.
Input
d– directionhp– hyperplaneerror_unbounded– (optional, default:true)trueif an error should be thrown whenever the set is unbounded in the given directionhalfspace– (optional, default:false)trueif the support vector should be computed for a half-space
Output
A pair (v, b) where v is a vector and b is a Boolean flag.
The flag b is false in one of the following cases:
- The direction has norm zero.
- The direction is the hyperplane's normal direction.
- The direction is the opposite of the hyperplane's normal direction and
halfspace is false. In all these cases, v is any point on the hyperplane.
Otherwise, the flag b is true, the set is unbounded in the given direction, and v is any vector.
If error_unbounded is true and the set is unbounded in the given direction, this function throws an error instead of returning.
Notes
For correctness, consider the weak duality of LPs: If the primal is unbounded, then the dual is infeasible. Since there is only a single constraint, the feasible set of the dual problem is hp.a ⋅ y == d, y >= 0 (with objective function hp.b ⋅ y). It is easy to see that this problem is infeasible whenever a is not parallel to d.
LazySets.@neutral — Macro.@neutral(SET, NEUT)Create functions to make a lazy set operation commutative with a given neutral element set type.
Input
SET– lazy set operation typeNEUT– set type for neutral element
Output
Nothing.
Notes
This macro generates four functions (possibly two more if @absorbing has been used in advance) (possibly two or four more if @declare_array_version has been used in advance).
Examples
@neutral(MinkowskiSum, N) creates at least the following functions:
neutral(::MinkowskiSum) = NMinkowskiSum(X, N) = XMinkowskiSum(N, X) = XMinkowskiSum(N, N) = N
LazySets.@absorbing — Macro.@absorbing(SET, ABS)Create functions to make a lazy set operation commutative with a given absorbing element set type.
Input
SET– lazy set operation typeABS– set type for absorbing element
Output
Nothing.
Notes
This macro generates four functions (possibly two more if @neutral has been used in advance) (possibly two or four more if @declare_array_version has been used in advance).
Examples
@absorbing(MinkowskiSum, A) creates at least the following functions:
absorbing(::MinkowskiSum) = AMinkowskiSum(X, A) = AMinkowskiSum(A, X) = AMinkowskiSum(A, A) = A
LazySets.@neutral_absorbing — Macro.@neutral_absorbing(SET, NEUT, ABS)Create two functions to avoid method ambiguties for a lazy set operation with respect to neutral and absorbing element set types.
Input
SET– lazy set operation typeNEUT– set type for neutral elementABS– set type for absorbing element
Output
A quoted expression containing the function definitions.
Examples
@neutral_absorbing(MinkowskiSum, N, A) creates the following functions as quoted expressions:
MinkowskiSum(N, A) = AMinkowskiSum(A, N) = A
LazySets.@declare_array_version — Macro.@declare_array_version(SET, SETARR)Create functions to connect a lazy set operation with its array set type.
Input
SET– lazy set operation typeSETARR– array set type
Output
Nothing.
Notes
This macro generates six functions (and possibly up to eight more if @neutral/@absorbing has been used in advance for the base and/or array set type).
Examples
@declare_array_version(MinkowskiSum, MinkowskiSumArray) creates at least the following functions:
array_constructor(::MinkowskiSum) = MinkowskiSumArrayis_array_constructor(::MinkowskiSumArray) = trueMinkowskiSum!(X, Y)MinkowskiSum!(X, arr)MinkowskiSum!(arr, X)MinkowskiSum!(arr1, arr2)
LazySets.@array_neutral — Macro.@array_neutral(FUN, NEUT, SETARR)Create two functions to avoid method ambiguities for a lazy set operation with respect to the neutral element set type and the array set type.
Input
FUN– function nameNEUT– set type for neutral elementSETARR– array set type
Output
A quoted expression containing the function definitions.
Examples
@array_neutral(MinkowskiSum, N, ARR) creates the following functions as quoted expressions:
MinkowskiSum(N, ARR) = ARRMinkowskiSum(ARR, N) = ARR
LazySets.@array_absorbing — Macro.@array_absorbing(FUN, ABS, SETARR)Create two functions to avoid method ambiguities for a lazy set operation with respect to the absorbing element set type and the array set type.
Input
FUN– function nameABS– set type for absorbing elementSETARR– array set type
Output
A quoted expression containing the function definitions.
Examples
@array_absorbing(MinkowskiSum, ABS, ARR) creates the following functions as quoted expressions:
MinkowskiSum(ABS, ARR) = ABSMinkowskiSum(ARR, ABS) = ABS
Types
LazySets.CachedPair — Type.CachedPair{N}A mutable pair of an index and a vector.
Fields
idx– indexvec– vector
UnitVector{T} <: AbstractVector{T}A lazy unit vector with arbitrary one-element.
Fields
i– index of non-zero entryn– vector lengthv– non-zero entry
StrictlyIncreasingIndicesIterator over the vectors of m strictly increasing indices from 1 to n.
Fields
n– size of the index domainm– number of indices to choose (resp. length of the vectors)
Notes
The vectors are modified in-place.
The iterator ranges over $\binom{n}{m}$ (n choose m) possible vectors.
This implementation results in a lexicographic order with the last index growing first.
Examples
julia> for v in LazySets.StrictlyIncreasingIndices(4, 2)
println(v)
end
[1, 2]
[1, 3]
[1, 4]
[2, 3]
[2, 4]
[3, 4]Inspection of set interfaces
InteractiveUtils.subtypes — Method.subtypes(interface, concrete::Bool)Return the concrete subtypes of a given interface.
Input
interface– an abstract type, usually a set interfaceconcrete– iftrue, seek further the inner abstract subtypes of the given interface, otherwise return only the direct subtypes ofinterface
Output
A list with the subtypes of the abstract type interface, sorted alphabetically.
Examples
Consider the AbstractPolytope interface. If we include the abstract subtypes of this interface,
julia> using LazySets: subtypes
julia> subtypes(AbstractPolytope, false)
4-element Array{Any,1}:
AbstractCentrallySymmetricPolytope
AbstractPolygon
HPolytope
VPolytopeWe can use this function to obtain the concrete subtypes of AbstractCentrallySymmetricPolytope and AbstractPolygon (further until all concrete types are obtained), using the concrete flag:
julia> subtypes(AbstractPolytope, true)
14-element Array{Type,1}:
Ball1
BallInf
HPolygon
HPolygonOpt
HPolytope
Hyperrectangle
Interval
LineSegment
Singleton
SymmetricIntervalHull
VPolygon
VPolytope
ZeroSet
Zonotope