Common Set Operations

Common Set Operations

This section of the manual describes the basic symbolic types describing operations between sets.

Cartesian Product

Binary Cartesian Product

CartesianProduct{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents a Cartesian product of two convex sets.

Fields

  • X – first convex set

  • Y – second convex set

Notes

The Cartesian product of three elements is obtained recursively. See also CartesianProductArray for an implementation of a Cartesian product of many sets without recursion, instead using an array.

  • CartesianProduct{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}}(X1::S1, X2::S2) – default constructor

  • CartesianProduct(Xarr::Vector{S}) where {S<:LazySet} – constructor from an array of convex sets

source
LazySets.dimMethod.
dim(cp::CartesianProduct)::Int

Return the dimension of a Cartesian product.

Input

  • cp – Cartesian product

Output

The ambient dimension of the Cartesian product.

source
LazySets.σMethod.
σ(d::AbstractVector{<:Real}, cp::CartesianProduct)::AbstractVector{<:Real}

Return the support vector of a Cartesian product.

Input

  • d – direction

  • cp – Cartesian product

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the product sets.

source
Base.:*Method.
    *(X::LazySet, Y::LazySet)::CartesianProduct

Return the Cartesian product of two convex sets.

Input

  • X – convex set

  • Y – convex set

Output

The Cartesian product of the two convex sets.

source
Base.:∈Method.
∈(x::AbstractVector{<:Real}, cp::CartesianProduct)::Bool

Check whether a given point is contained in a Cartesian product set.

Input

  • x – point/vector

  • cp – Cartesian product

Output

true iff $x ∈ cp$.

source

$n$-ary Cartesian Product

CartesianProductArray{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the Cartesian product of a finite number of convex sets.

Fields

  • sfarray – array of sets

Notes

  • CartesianProductArray(sfarray::Vector{<:LazySet}) – default constructor

  • CartesianProductArray() – constructor for an empty Cartesian product

  • CartesianProductArray(n::Int, [N]::Type=Float64) – constructor for an empty Cartesian product with size hint and numeric type

source
LazySets.dimMethod.
dim(cpa::CartesianProductArray)::Int

Return the dimension of a Cartesian product of a finite number of convex sets.

Input

  • cpa – Cartesian product array

Output

The ambient dimension of the Cartesian product of a finite number of convex sets.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, cpa::CartesianProductArray{N, <:LazySet{N}}
 )::AbstractVector{N} where {N<:Real}

Support vector of a Cartesian product.

Input

  • d – direction

  • cpa – Cartesian product array

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the product sets.

source
Base.:*Method.
    *(X::LazySet, Y::LazySet)::CartesianProduct

Return the Cartesian product of two convex sets.

Input

  • X – convex set

  • Y – convex set

Output

The Cartesian product of the two convex sets.

source
    *(cpa::CartesianProductArray, S::LazySet)::CartesianProductArray

Multiply a convex set to a Cartesian product of a finite number of convex sets from the right.

Input

  • cpa – Cartesian product array (is modified)

  • S – convex set

Output

The modified Cartesian product of a finite number of convex sets.

source
Base.:*Method.
    *(X::LazySet, Y::LazySet)::CartesianProduct

Return the Cartesian product of two convex sets.

Input

  • X – convex set

  • Y – convex set

Output

The Cartesian product of the two convex sets.

source
    *(S::LazySet, cpa::CartesianProductArray)::CartesianProductArray

Multiply a convex set to a Cartesian product of a finite number of convex sets from the left.

Input

  • S – convex set

  • cpa – Cartesian product array (is modified)

Output

The modified Cartesian product of a finite number of convex sets.

source
Base.:*Method.
    *(X::LazySet, Y::LazySet)::CartesianProduct

Return the Cartesian product of two convex sets.

Input

  • X – convex set

  • Y – convex set

Output

The Cartesian product of the two convex sets.

source
    *(cpa::CartesianProductArray, S::LazySet)::CartesianProductArray

Multiply a convex set to a Cartesian product of a finite number of convex sets from the right.

Input

  • cpa – Cartesian product array (is modified)

  • S – convex set

Output

The modified Cartesian product of a finite number of convex sets.

source
    *(S::LazySet, cpa::CartesianProductArray)::CartesianProductArray

Multiply a convex set to a Cartesian product of a finite number of convex sets from the left.

Input

  • S – convex set

  • cpa – Cartesian product array (is modified)

Output

The modified Cartesian product of a finite number of convex sets.

source
    *(cpa1::CartesianProductArray, cpa2::CartesianProductArray)::CartesianProductArray

Multiply a finite Cartesian product of convex sets to another finite Cartesian product.

Input

  • cpa1 – first Cartesian product array (is modified)

  • cpa2 – second Cartesian product array

Output

The modified first Cartesian product.

source
Base.:∈Method.
∈(x::AbstractVector{N}, cpa::CartesianProductArray{N, <:LazySet{N}}
 )::Bool  where {N<:Real}

Check whether a given point is contained in a Cartesian product of a finite number of sets.

Input

  • x – point/vector

  • cpa – Cartesian product array

Output

true iff $x ∈ \text{cpa}$.

source

Convex Hull

Binary Convex Hull

ConvexHull{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents the convex hull of the union of two convex sets.

Fields

  • X – convex set

  • Y – convex set

Examples

Convex hull of two 100-dimensional Euclidean balls:

julia> b1, b2 = Ball2(zeros(100), 0.1), Ball2(4*ones(100), 0.2);

julia> c = ConvexHull(b1, b2);

julia> typeof(c)
LazySets.ConvexHull{Float64,LazySets.Ball2{Float64},LazySets.Ball2{Float64}}
source
LazySets.CHType.
CH

Alias for ConvexHull.

source
LazySets.dimMethod.
dim(ch::ConvexHull)::Int

Return the dimension of a convex hull of two convex sets.

Input

  • ch – convex hull of two convex sets

Output

The ambient dimension of the convex hull of two convex sets.

source
LazySets.σMethod.
σ(d::AbstractVector{<:Real}, ch::ConvexHull)::AbstractVector{<:Real}

Return the support vector of a convex hull of two convex sets in a given direction.

Input

  • d – direction

  • ch – convex hull of two convex sets

source

$n$-ary Convex Hull

ConvexHullArray{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the symbolic convex hull of a finite number of convex sets.

Fields

  • array – array of sets

Examples

Convex hull of 100 two-dimensional balls whose centers follows a sinusoidal:

julia> b = [Ball2([2*pi*i/100, sin(2*pi*i/100)], 0.05) for i in 1:100];

julia> c = ConvexHullArray(b);
source
CHArray

Alias for ConvexHullArray.

source
LazySets.dimMethod.
dim(cha::ConvexHullArray)::Int

Return the dimension of the convex hull of a finite number of convex sets.

Input

  • cha – convex hull array

Output

The ambient dimension of the convex hull of a finite number of convex sets.

source
LazySets.σMethod.
σ(d::AbstractVector{<:Real}, cha::ConvexHullArray)::AbstractVector{<:Real}

Return the support vector of a convex hull array in a given direction.

Input

  • d – direction

  • cha – convex hull array

source

Convex Hull Algorithms

LazySets.convex_hullFunction.
convex_hull(points::Vector{S}; [algorithm]::String="monotone_chain"
           )::Vector{S} where {S<:AbstractVector{N}} where {N<:Real}

Compute the convex hull of points in the plane.

Input

  • points – list of 2D vectors

  • algorithm – (optional, default: "monotone_chain") the convex hull algorithm, valid options are:

    • "monotone_chain"

    • "monotone_chain_sorted"

Output

The convex hull as a list of 2D vectors with the coordinates of the points.

Examples

Compute the convex hull of a random set of points:

julia> points = [randn(2) for i in 1:30]; # 30 random points in 2D

julia> hull = convex_hull(points);

julia> typeof(hull)
Array{Array{Float64,1},1}

Plot both the random points and the computed convex hull polygon:

julia> using Plots;

julia> plot([Tuple(pi) for pi in points], seriestype=:scatter);

julia> plot!(VPolygon(hull), alpha=0.2);
source
LazySets.convex_hull!Function.
convex_hull!(points::Vector{S}; [algorithm]::String="monotone_chain"
            )::Vector{S} where {S<:AbstractVector{N}} where {N<:Real}

Compute the convex hull of points in the plane, in-place.

Input

  • points – list of 2D vectors (is modified)

  • algorithm – (optional, default: "monotone_chain") the convex hull algorithm; valid options are:

    • "monotone_chain"

    • "monotone_chain_sorted"

Output

The convex hull as a list of 2D vectors with the coordinates of the points.

Notes

See the non-modifying version convex_hull for more details.

source
LazySets.right_turnFunction.
right_turn(O::AbstractVector{N}, A::AbstractVector{N}, B::AbstractVector{N}
          )::N where {N<:Real}

Determine if the acute angle defined by the three points O, A, B in the plane is a right turn (counter-clockwise) with respect to the center O.

Input

  • O – 2D center point

  • A – 2D one point

  • B – 2D another point

Output

Scalar representing the rotation.

Algorithm

The cross product is used to determine the sense of rotation. If the result is 0, the points are collinear; if it is positive, the three points constitute a positive angle of rotation around O from A to B; otherwise they constitute a negative angle.

source
monotone_chain!(points::Vector{S}; sort::Bool=true
               )::Vector{S} where {S<:AbstractVector{N}} where {N<:Real}

Compute the convex hull of points in the plane using Andrew's monotone chain method.

Input

  • points – list of 2D vectors; is sorted in-place inside this function

  • sort – (optional, default: true) flag for sorting the vertices lexicographically; sortedness is required for correctness

Output

List of vectors containing the 2D coordinates of the corner points of the convex hull.

Notes

For large sets of points, it is convenient to use static vectors to get maximum performance. For information on how to convert usual vectors into static vectors, see the type SVector provided by the StaticArrays package.

Algorithm

This function implements Andrew's monotone chain convex hull algorithm to construct the convex hull of a set of $n$ points in the plane in $O(n \log n)$ time. For further details see Monotone chain

source

Intersection

Intersection{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents the intersection of two convex sets.

Fields

  • X – convex set

  • Y – convex set

source
LazySets.dimMethod.
dim(cap::Intersection)::Int

Return the dimension of an intersection of two convex sets.

Input

  • cap – intersection of two convex sets

Output

The ambient dimension of the intersection of two convex sets.

source
LazySets.σMethod.
σ(d::AbstractVector{<:Real}, cap::Intersection)::AbstractVector{<:Real}

Return the support vector of an intersection of two convex sets in a given direction.

Input

  • d – direction

  • cap – intersection of two convex sets

Output

The support vector in the given direction.

source
Base.:∈Method.
∈(x::AbstractVector{N}, cap::Intersection{N})::Bool where {N<:Real}

Check whether a given point is contained in a intersection of two convex sets.

Input

  • x – point/vector

  • cap – intersection of two convex sets

Output

true iff $x ∈ cap$.

source
Base.isemptyMethod.
isempty(cap::Intersection)::Bool

Return if the intersection is empty or not.

Input

  • cap – intersection of two convex sets

Output

true iff the intersection is empty.

source

Minkowski Sum

Binary Minkowski Sum

MinkowskiSum{N<:Real, S1<:LazySet{N}, S2<:LazySet{N}} <: LazySet{N}

Type that represents the Minkowski sum of two convex sets.

Fields

  • X – first convex set

  • Y – second convex set

source
LazySets.dimMethod.
dim(ms::MinkowskiSum)::Int

Return the dimension of a Minkowski sum.

Input

  • ms – Minkowski sum

Output

The ambient dimension of the Minkowski sum.

source
LazySets.σMethod.
σ(d::AbstractVector{<:Real}, ms::MinkowskiSum)::AbstractVector{<:Real}

Return the support vector of a Minkowski sum.

Input

  • d – direction

  • ms – Minkowski sum

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the summand sets.

source
Base.:+Method.
X + Y

Convenience constructor for Minkowski sum.

Input

  • X – a convex set

  • Y – another convex set

Output

The symbolic Minkowski sum of $X$ and $Y$.

source
LazySets.:⊕Function.
⊕(X::LazySet, Y::LazySet)

Unicode alias constructor oplus for the Minkowski sum operator +(X, Y).

source

$n$-ary Minkowski Sum

MinkowskiSumArray{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the Minkowski sum of a finite number of convex sets.

Fields

  • sfarray – array of convex sets

Notes

This type assumes that the dimensions of all elements match.

  • MinkowskiSumArray(sfarray::Vector{<:LazySet}) – default constructor

  • MinkowskiSumArray() – constructor for an empty sum

  • MinkowskiSumArray(n::Int, [N]::Type=Float64) – constructor for an empty sum with size hint and numeric type

source
LazySets.dimMethod.
dim(msa::MinkowskiSumArray)::Int

Return the dimension of a Minkowski sum of a finite number of sets.

Input

  • msa – Minkowski sum array

Output

The ambient dimension of the Minkowski sum of a finite number of sets.

source
LazySets.σMethod.
σ(d::AbstractVector{<:Real}, msa::MinkowskiSumArray)::Vector{<:Real}

Return the support vector of a Minkowski sum of a finite number of sets in a given direction.

Input

  • d – direction

  • msa – Minkowski sum array

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the summand sets.

source
Base.:+Method.
X + Y

Convenience constructor for Minkowski sum.

Input

  • X – a convex set

  • Y – another convex set

Output

The symbolic Minkowski sum of $X$ and $Y$.

source
+(msa::MinkowskiSumArray, S::LazySet)::MinkowskiSumArray

Add a convex set to a Minkowski sum of a finite number of convex sets from the right.

Input

  • msa – Minkowski sum array (is modified)

  • S – convex set

Output

The modified Minkowski sum of a finite number of convex sets.

source
Base.:+Method.
X + Y

Convenience constructor for Minkowski sum.

Input

  • X – a convex set

  • Y – another convex set

Output

The symbolic Minkowski sum of $X$ and $Y$.

source
+(S::LazySet, msa::MinkowskiSumArray)::MinkowskiSumArray

Add a convex set to a Minkowski sum of a finite number of convex sets from the left.

Input

  • S – convex set

  • msa – Minkowski sum array (is modified)

Output

The modified Minkowski sum of a finite number of convex sets.

source
Base.:+Method.
X + Y

Convenience constructor for Minkowski sum.

Input

  • X – a convex set

  • Y – another convex set

Output

The symbolic Minkowski sum of $X$ and $Y$.

source
+(msa::MinkowskiSumArray, S::LazySet)::MinkowskiSumArray

Add a convex set to a Minkowski sum of a finite number of convex sets from the right.

Input

  • msa – Minkowski sum array (is modified)

  • S – convex set

Output

The modified Minkowski sum of a finite number of convex sets.

source
+(S::LazySet, msa::MinkowskiSumArray)::MinkowskiSumArray

Add a convex set to a Minkowski sum of a finite number of convex sets from the left.

Input

  • S – convex set

  • msa – Minkowski sum array (is modified)

Output

The modified Minkowski sum of a finite number of convex sets.

source
+(msa1::MinkowskiSumArray, msa2::MinkowskiSumArray)::MinkowskiSumArray

Add the elements of a finite Minkowski sum of convex sets to another finite Minkowski sum.

Input

  • msa1 – first Minkowski sum array (is modified)

  • msa2 – second Minkowski sum array

Output

The modified first Minkowski sum of a finite number of convex sets.

source
Base.:+Method.
X + Y

Convenience constructor for Minkowski sum.

Input

  • X – a convex set

  • Y – another convex set

Output

The symbolic Minkowski sum of $X$ and $Y$.

source
+(msa::MinkowskiSumArray, S::LazySet)::MinkowskiSumArray

Add a convex set to a Minkowski sum of a finite number of convex sets from the right.

Input

  • msa – Minkowski sum array (is modified)

  • S – convex set

Output

The modified Minkowski sum of a finite number of convex sets.

source
+(msa::MinkowskiSumArray, Z::ZeroSet)::MinkowskiSumArray

Returns the original array because addition with an empty set is a no-op.

Input

  • msa – Minkowski sum array

  • Z – a Zero set

source

Maps

Linear Map

LinearMap{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents a linear transformation $M⋅S$ of a convex set $S$.

Fields

  • M – matrix/linear map

  • sf – convex set

source
LazySets.dimMethod.
dim(lm::LinearMap)::Int

Return the dimension of a linear map.

Input

  • lm – linear map

Output

The ambient dimension of the linear map.

source
LazySets.σMethod.
σ(d::AbstractVector{<:Real}, lm::LinearMap)::AbstractVector{<:Real}

Return the support vector of the linear map.

Input

  • d – direction

  • lm – linear map

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the wrapped set.

Notes

If $L = M⋅S$, where $M$ is a matrix and $S$ is a convex set, it follows that $σ(d, L) = M⋅σ(M^T d, S)$ for any direction $d$.

source
Base.:*Method.
    *(M::AbstractMatrix{<:Real}, S::LazySet)

Return the linear map of a convex set.

Input

  • M – matrix/linear map

  • S – convex set

Output

If the matrix is null, a ZeroSet is returned; otherwise a lazy linear map.

source
Base.:*Method.
    *(a::N, X::S)::LinearMap{N, S} where {S<:LazySet{N}} where {N<:Real}

Return a linear map of a convex set by a scalar value.

Input

  • a – real scalar

  • S – convex set

Output

The linear map of the convex set.

source
Base.:∈Method.
∈(x::AbstractVector{N}, lm::LinearMap{N, <:LazySet})::Bool where {N<:Real}

Check whether a given point is contained in a linear map of a convex set.

Input

  • x – point/vector

  • lm – linear map of a convex set

Output

true iff $x ∈ lm$.

Algorithm

Note that $x ∈ M⋅S$ iff $M^{-1}⋅x ∈ S$. This implementation does not explicitly invert the matrix, which is why it also works for non-square matrices.

Examples

julia> lm = LinearMap([2.0 0.0; 0.0 1.0], BallInf([1., 1.], 1.));

julia> ∈([5.0, 1.0], lm)
false
julia> ∈([3.0, 1.0], lm)
true

An example with non-square matrix:

julia> B = BallInf(zeros(4), 1.);

julia> M = [1. 0 0 0; 0 1 0 0]/2;

julia> ∈([0.5, 0.5], M*B)
true
source

Exponential Map

ExponentialMap{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the action of an exponential map on a convex set.

Fields

  • spmexp – sparse matrix exponential

  • X – convex set

source
LazySets.dimMethod.
dim(em::ExponentialMap)::Int

Return the dimension of an exponential map.

Input

  • em – an ExponentialMap

Output

The ambient dimension of the exponential map.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, em::ExponentialMap)::AbstractVector{N} where {N<:Real}

Return the support vector of the exponential map.

Input

  • d – direction

  • em – exponential map

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the wrapped set.

Notes

If $E = \exp(M)⋅S$, where $M$ is a matrix and $S$ is a convex set, it follows that $σ(d, E) = \exp(M)⋅σ(\exp(M)^T d, S)$ for any direction $d$.

source
Base.:∈Method.
∈(x::AbstractVector{N}, em::ExponentialMap{<:LazySet{N}})::Bool where {N<:Real}

Check whether a given point is contained in an exponential map of a convex set.

Input

  • x – point/vector

  • em – linear map of a convex set

Output

true iff $x ∈ em$.

Algorithm

This implementation exploits that $x ∈ \exp(M)⋅S$ iff $\exp(-M)⋅x ∈ S$. This follows from $\exp(-M)⋅\exp(M) = I$ for any $M$.

Examples

julia> em = ExponentialMap(SparseMatrixExp(SparseMatrixCSC([2.0 0.0; 0.0 1.0])),
                           BallInf([1., 1.], 1.));

julia> ∈([5.0, 1.0], em)
false
julia> ∈([1.0, 1.0], em)
true
source
ExponentialProjectionMap{N<:Real, S<:LazySet{N}} <: LazySet{N}

Type that represents the application of a projection of a sparse matrix exponential to a convex set.

Fields

  • spmexp – projection of a sparse matrix exponential

  • X – convex set

source
LazySets.dimMethod.
dim(eprojmap::ExponentialProjectionMap)::Int

Return the dimension of a projection of an exponential map.

Input

  • eprojmap – projection of an exponential map

Output

The ambient dimension of the projection of an exponential map.

source
LazySets.σMethod.
σ(d::AbstractVector{N},
  eprojmap::ExponentialProjectionMap)::AbstractVector{N} where {N<:Real}

Return the support vector of a projection of an exponential map.

Input

  • d – direction

  • eprojmap – projection of an exponential map

Output

The support vector in the given direction. If the direction has norm zero, the result depends on the wrapped set.

Notes

If $S = (L⋅M⋅R)⋅X$, where $L$ and $R$ are matrices, $M$ is a matrix exponential, and $X$ is a set, it follows that $σ(d, S) = L⋅M⋅R⋅σ(R^T⋅M^T⋅L^T⋅d, X)$ for any direction $d$.

source
SparseMatrixExp{N}

Type that represents the matrix exponential, $\exp(M)$, of a sparse matrix.

Fields

  • M – sparse matrix

Notes

This type is provided for use with very large and very sparse matrices. The evaluation of the exponential matrix action over vectors relies on the Expokit package.

source
Base.:*Method.
    *(spmexp::SparseMatrixExp, X::LazySet)::ExponentialMap

Return the exponential map of a convex set from a sparse matrix exponential.

Input

  • spmexp – sparse matrix exponential

  • X – convex set

Output

The exponential map of the convex set.

source
ProjectionSparseMatrixExp{N<:Real}

Type that represents the projection of a sparse matrix exponential, i.e., $L⋅\exp(M)⋅R$ for a given sparse matrix $M$.

Fields

  • L – left multiplication matrix

  • E – sparse matrix exponential

  • R – right multiplication matrix

source
Base.:*Method.
    *(projspmexp::ProjectionSparseMatrixExp,
      X::LazySet)::ExponentialProjectionMap

Return the application of a projection of a sparse matrix exponential to a convex set.

Input

  • projspmexp – projection of a sparse matrix exponential

  • X – convex set

Output

The application of the projection of a sparse matrix exponential to the convex set.

source

Symmetric Interval Hull

SymmetricIntervalHull{N<:Real, S<:LazySet{N}} <: AbstractHyperrectangle{N}

Type that represents the symmetric interval hull of a convex set.

Fields

  • X – convex set

  • cache – partial storage of already computed bounds, organized as mapping from dimension to tuples (bound, valid), where valid is a flag indicating if the bound entry has been computed

Notes

The symmetric interval hull can be computed with $2n$ support vector queries of unit vectors, where $n$ is the dimension of the wrapped set (i.e., two queries per dimension). When asking for the support vector for a direction $d$, one needs $2k$ such queries, where $k$ is the number of non-zero entries in $d$.

However, if one asks for many support vectors in a loop, the number of computations may exceed $2n$. To be most efficient in such cases, this type stores the intermediately computed bounds in the cache field.

source
LazySets.dimMethod.
dim(P::AbstractPointSymmetricPolytope)::Int

Return the ambient dimension of a point symmetric set.

Input

  • P – set

Output

The ambient dimension of the set.

source
dim(sih::SymmetricIntervalHull)::Int

Return the dimension of a symmetric interval hull of a convex set.

Input

  • sih – symmetric interval hull of a convex set

Output

The ambient dimension of the symmetric interval hull of a convex set.

source
LazySets.σMethod.
σ(d::AbstractVector{N}, H::AbstractHyperrectangle{N}
 )::AbstractVector{N} where {N<:Real}

Return the support vector of a hyperrectangular set in a given direction.

Input

  • d – direction

  • H – hyperrectangular set

Output

The support vector in the given direction. If the direction has norm zero, the vertex with biggest values is returned.

source
σ(d::AbstractVector{N}, sih::SymmetricIntervalHull
 )::AbstractVector{<:Real} where {N<:Real}

Return the support vector of a symmetric interval hull of a convex set in a given direction.

Input

  • d – direction

  • sih – symmetric interval hull of a convex set

Output

The support vector of the symmetric interval hull of a convex set in the given direction. If the direction has norm zero, the origin is returned.

Algorithm

For each non-zero entry in d we need to either look up the bound (if it has been computed before) or compute it, in which case we store it for future queries. One such computation just asks for the support vector of the underlying set for both the positive and negative unit vector in the respective dimension.

source
an_element(P::AbstractPointSymmetricPolytope{N})::Vector{N} where {N<:Real}

Return some element of a point symmetric polytope.

Input

  • P – point symmetric polytope

Output

The center of the point symmetric polytope.

source