Approximations
This section of the manual describes the Cartesian decomposition algorithms and the approximation of high-dimensional convex sets using projections.
LazySets.Approximations
— Module.Module Approximations.jl
– polygonal approximation of convex sets through support vectors.
Cartesian Decomposition
LazySets.Approximations.decompose
— Function.decompose(S::LazySet{N},
partition::AbstractVector{<:AbstractVector{Int}},
block_options
)::CartesianProductArray{N} where {N<:Real}
Decompose a high-dimensional set into a Cartesian product of overapproximations of the projections over the specified subspaces.
Input
S
– setpartition
– vector of blocks (i.e., of vectors of integers) (see the Notes below)block_options
– mapping from block indices inpartition
to a corresponding overapproximation option; we only require access via[⋅]
(but see also the Notes below)
Output
A CartesianProductArray
containing the low-dimensional approximated projections.
Algorithm
For each block a specific project
method is called, dispatching on the corresponding overapproximation option.
Notes
The argument partition
requires some discussion. Typically, the list of blocks should form a partition of the set $\{1, \dots, n\}$ represented as a list of consecutive blocks, where $n$ is the ambient dimension of set S
.
However, technically there is no problem if the blocks are not consecutive, blocks are missing, blocks occur more than once, or blocks are overlapping. This function will, however, stick to the order of blocks, so the resulting set must be interpreted with care in such cases. One use case is the need of a projection consisting of several blocks.
For convenience, the argument block_options
can also be given as a single option instead of a mapping, which is then interpreted as the option for all blocks.
Examples
This function supports different options: one can specify the target set, the degree of accuracy, and template directions. These options are exemplified below, where we use the following example.
julia> using LazySets.Approximations: decompose
julia> S = Ball2(zeros(4), 1.); # set to be decomposed (4D 2-norm unit ball)
julia> P2d = [1:2, 3:4]; # a partition with two blocks of size two
julia> P1d = [[1], [2], [3], [4]]; # a partition with four blocks of size one
Different set types
We can decompose using polygons in constraint representation:
julia> all([ai isa HPolygon for ai in array(decompose(S, P2d, HPolygon))])
true
For decomposition into 1D subspaces, we can use Interval
:
julia> all([ai isa Interval for ai in array(decompose(S, P1d, Interval))])
true
However, if you need to specify different set types for different blocks, the interface presented so far does not apply. See the paragraph Advanced input for different block approximations below for how to do that.
Refining the decomposition I: $ε$-close approximation
The $ε$ option can be used to refine a decomposition, i.e., obtain a more accurate result. We use the Iterative refinement algorithm from the Approximations
module.
To illustrate this, consider again the set S
from above. We decompose into two 2D polygons. Using smaller $ε$ implies a better precision, thus more constraints in each 2D decomposition. In the following example, we look at the number of constraints in the first block.
julia> d(ε, bi) = array(decompose(S, P2d, (HPolygon => ε)))[bi]
d (generic function with 1 method)
julia> [length(constraints_list(d(ε, 1))) for ε in [Inf, 0.1, 0.01]]
3-element Array{Int64,1}:
4
8
32
Refining the decomposition II: template polyhedra
Another way to refine a decomposition is by using template polyhedra. The idea is to specify a set of template directions and then to compute on each block the polytopic overapproximation obtained by evaluating the support function of the given input set over the template directions.
For example, octagonal 2D approximations of the set S
are obtained with:
julia> B = decompose(S, P2d, OctDirections);
julia> length(B.array) == 2 && all(dim(bi) == 2 for bi in B.array)
true
See Template directions for the available template directions. Note that, in contrast to the polygonal $ε$-close approximation from above, this method can be applied to blocks of any size.
julia> B = decompose(S, [1:4], OctDirections);
julia> length(B.array) == 1 && dim(B.array[1]) == 4
true
Advanced input for different block approximations
Instead of defining the approximation option uniformly for each block, we can define different approximations for different blocks. The third argument has to be a mapping from block index (in the partition) to the corresponding approximation option.
For example:
julia> res = array(decompose(S, P2d, Dict(1 => Hyperrectangle, 2 => 0.1)));
julia> typeof(res[1]), typeof(res[2])
(Hyperrectangle{Float64}, HPolygon{Float64})
LazySets.Approximations.project
— Function.project(S::LazySet{N},
block::AbstractVector{Int},
set_type::Type{<:LinearMap},
[n]::Int=dim(S)
)::LinearMap{N} where {N<:Real}
Project a high-dimensional set to a given block by using a lazy linear map.
Input
S
– setblock
– block structure - a vector with the dimensions of interestLinearMap
– used for dispatchn
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A lazy LinearMap
representing a projection of the set S
to block block
.
project(S::LazySet{N},
block::AbstractVector{Int},
set_type::Type{<:LazySet},
[n]::Int=dim(S)
) where {N<:Real}
Project a high-dimensional set to a given block and set type, possibly involving an overapproximation.
Input
S
– setblock
– block structure - a vector with the dimensions of interestset_type
– target set typen
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A set of type set_type
representing an overapproximation of the projection of S
.
Algorithm
- Project the set
S
withM⋅S
, whereM
is the identity matrix in the block
coordinates and zero otherwise.
- Overapproximate the projected lazy set using
overapproximate
and
set_type
.
project(S::LazySet{N},
block::AbstractVector{Int},
set_type_and_precision::Pair{<:UnionAll, <:Real},
[n]::Int=dim(S)
) where {N<:Real}
Project a high-dimensional set to a given block and set type with a certified error bound.
Input
S
– setblock
– block structure - a vector with the dimensions of interestset_type_and_precision
– pair(T, ε)
of a target set typeT
and an error boundε
for approximationn
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A set representing the epsilon-close approximation of the projection of S
.
Notes
Currently we only support HPolygon
as set type, which implies that the set must be two-dimensional.
Algorithm
- Project the set
S
withM⋅S
, whereM
is the identity matrix in the block
coordinates and zero otherwise.
- Overapproximate the projected lazy set with the given error bound
ε
.
project(S::LazySet{N},
block::AbstractVector{Int},
ε::Real,
[n]::Int=dim(S)
) where {N<:Real}
Project a high-dimensional set to a given block and set type with a certified error bound.
Input
S
– setblock
– block structure - a vector with the dimensions of interestε
– error bound for approximationn
– (optional, default:dim(S)
) ambient dimension of the setS
Output
A set representing the epsilon-close approximation of the projection of S
.
Algorithm
- Project the set
S
withM⋅S
, whereM
is the identity matrix in the block
coordinates and zero otherwise.
- Overapproximate the projected lazy set with the given error bound
ε
.
The target set type is chosen automatically.
project(S::LazySet{N},
block::AbstractVector{Int},
directions::Type{<:AbstractDirections},
[n]::Int
) where {N<:Real}
Project a high-dimensional set to a given block using template directions.
Input
S
– setblock
– block structure - a vector with the dimensions of interestdirections
– template directionsn
– (optional, default:dim(S)
) ambient dimension of the setS
Output
The template direction approximation of the projection of S
.
project(H::HalfSpace{N}, block::AbstractVector{Int})
Concrete projection of a half-space.
Input
H
– setblock
– block structure, a vector with the dimensions of interest
Output
A set representing the projection of the half-space H
on the dimensions specified by block
.
Notes
Currently only the case where the unconstrained dimensions of H
are a subset of the block
variables is implemented.
Examples
Consider the half-space $x + y + 0⋅z ≤ 1$, whose ambient dimension is 3
. The (trivial) projection in the three dimensions is achieved letting the block of variables to be [1, 2, 3]
:
julia> H = HalfSpace([1.0, 1.0, 0.0], 1.0)
HalfSpace{Float64}([1.0, 1.0, 0.0], 1.0)
julia> using LazySets.Approximations: project
julia> project(H, [1, 2, 3])
HalfSpace{Float64}([1.0, 1.0, 0.0], 1.0)
Projecting along dimensions 1
and 2
only:
julia> project(H, [1, 2])
HalfSpace{Float64}([1.0, 1.0], 1.0)
In general, use the call syntax project(H, constrained_dimensions(H))
to return the half-space projected on the dimensions where it is constrained only:
julia> project(H, constrained_dimensions(H))
HalfSpace{Float64}([1.0, 1.0], 1.0)
project(P::HPolyhedron{N}, block::AbstractVector{Int}) where {N}
Concrete projection of a polyhedron in half-space representation.
Input
P
– setblock
– block structure, a vector with the dimensions of interest
Output
A set representing the projection of P
on the dimensions specified by block
.
Notes
Currently only the case where the unconstrained dimensions of P
are a subset of the block
variables is implemented.
Examples
Consider the four-dimensional cross-polytope (unit ball in the 1-norm):
julia> using LazySets.Approximations: project
julia> P = convert(HPolyhedron, Ball1(zeros(4), 1.0));
All dimensions are constrained, and computing the (trivial) projection on the whole space behaves as expected:
julia> constrained_dimensions(P)
4-element Array{Int64,1}:
1
2
3
4
julia> P_1234 = project(P, [1, 2, 3, 4]);
julia> P_1234 == P
true
Each constraint of the cross polytope is constrained in all dimensions.
Now let's take a ball in the infinity norm and remove some constraints:
julia> B = BallInf(zeros(4), 1.0);
julia> c = constraints_list(B)[1:2]
2-element Array{HalfSpace{Float64},1}:
HalfSpace{Float64}([1.0, 0.0, 0.0, 0.0], 1.0)
HalfSpace{Float64}([0.0, 1.0, 0.0, 0.0], 1.0)
julia> P = HPolyhedron(c);
julia> constrained_dimensions(P)
2-element Array{Int64,1}:
1
2
Finally we take the concrete projection onto variables 1
and 2
:
julia> project(P, [1, 2]) |> constraints_list
2-element Array{HalfSpace{Float64},1}:
HalfSpace{Float64}([1.0, 0.0], 1.0)
HalfSpace{Float64}([0.0, 1.0], 1.0)
Convenience functions
LazySets.Approximations.uniform_partition
— Function. uniform_partition(n::Int, block_size::Int)
Compute a uniform block partition of the given size.
Input
n
– number of dimensions of the partitionblock_size
– size of each block
Output
A vector of ranges, Vector{UnitRange{Int}}
, such that the size of each block is the same, if possible.
Examples
If the number of dimensions n
is 2, we have two options: either two blocks of size 1
or one block of size 2
:
julia> LazySets.Approximations.uniform_partition(2, 1)
2-element Array{UnitRange{Int64},1}:
1:1
2:2
julia> LazySets.Approximations.uniform_partition(2, 2)
1-element Array{UnitRange{Int64},1}:
1:2
If the block size argument is not compatible with (i.e. does not divide) n
, the output is filled with one block of the size needed to reach n
:
julia> LazySets.Approximations.uniform_partition(3, 1)
3-element Array{UnitRange{Int64},1}:
1:1
2:2
3:3
julia> LazySets.Approximations.uniform_partition(3, 2)
2-element Array{UnitRange{Int64},1}:
1:2
3:3
julia> LazySets.Approximations.uniform_partition(10, 6)
2-element Array{UnitRange{Int64},1}:
1:6
7:10
Overapproximations
LazySets.Approximations.overapproximate
— Function.overapproximate(X::S, ::Type{S}) where {S<:LazySet}
Overapproximating a set of type S
with type S
is a no-op.
Input
X
– setType{S}
– set type
Output
The input set.
overapproximate(S::LazySet{N},
::Type{<:HPolygon},
[ε]::Real=Inf)::HPolygon where {N<:Real}
Return an approximation of a given 2D convex set. If no error tolerance is given, or is Inf
, the result is a box-shaped polygon. Otherwise the result is an ε-close approximation as a polygon.
Input
S
– convex set, assumed to be two-dimensionalHPolygon
– type for dispatchε
– (optional, default:Inf
) error bound
Output
A polygon in constraint representation.
overapproximate(S::LazySet, ε::Real)::HPolygon
Alias for overapproximate(S, HPolygon, ε)
.
overapproximate(S::LazySet,
Type{<:Hyperrectangle})::Union{Hyperrectangle, EmptySet}
Return an approximation of a given set as a hyperrectangle.
Input
S
– setHyperrectangle
– type for dispatch
Output
A hyperrectangle.
overapproximate(S::CartesianProductArray{N, <:AbstractHyperrectangle{N}},
::Type{<:Hyperrectangle}) where {N<:Real}
Return a tight overapproximation of the cartesian product array of a finite number of convex sets with and hyperrectangle.
Input
S
– cartesian product array of a finite number of convex setHyperrectangle
– type for dispatch
Output
A hyperrectangle.
Algorithm
This method falls back to the corresponding convert
method. Since the sets wrapped by the cartesian product array are hyperrectangles, it can be done efficiently without overapproximation.
overapproximate(S::CartesianProduct{N, <:AbstractHyperrectangle{N}, <:AbstractHyperrectangle{N}},
::Type{<:Hyperrectangle}) where {N<:Real}
Return a tight overapproximation of the cartesian product of two hyperrectangles by a new hyperrectangle.
Input
S
– cartesian product of two hyperrectangular setsHyperrectangle
– type for dispatch
Output
A hyperrectangle.
Algorithm
This method falls back to the corresponding convert
method. Since the sets wrapped by the cartesian product are hyperrectangles, it can be done efficiently without overapproximation.
overapproximate(lm::LinearMap{N, <:AbstractHyperrectangle{N}},
::Type{Hyperrectangle}) where {N}
Return a tight overapproximation of the linear map of a hyperrectangular set using a hyperrectangle.
Input
S
– linear map of a hyperrectangular setHyperrectangle
– type for dispatch
Output
A hyperrectangle.
Algorithm
If c
and r
denote the center and vector radius of a hyperrectangle H
, a tight hyperrectangular overapproximation of M * H
is obtained by transforming c ↦ M*c
and r ↦ abs.(M) * c
, where abs.(⋅)
denotes the element-wise absolute value operator.
overapproximate(S::LazySet)::Union{Hyperrectangle, EmptySet}
Alias for overapproximate(S, Hyperrectangle)
.
overapproximate(S::ConvexHull{N, Zonotope{N}, Zonotope{N}},
::Type{<:Zonotope})::Zonotope where {N<:Real}
Overapproximate the convex hull of two zonotopes.
Input
S
– convex hull of two zonotopesZonotope
– type for dispatch
Algorithm
This function implements the method proposed in [1]. The convex hull of two zonotopes $Z₁$ and $Z₂$ of the same order, that we write
for $j = 1, 2$, can be overapproximated as follows:
If the zonotope order is not the same, this algorithm calls reduce_order
to reduce the order to the minimum of the arguments.
It should be noted that the output zonotope is not necessarily the minimal enclosing zonotope, which is in general expensive in high dimensions. This is further investigated in [2].
[1] Reachability of Uncertain Linear Systems Using Zonotopes, A. Girard. HSCC 2005.
[2] Zonotopes as bounding volumes, L. J. Guibas et al, Proc. of Symposium on Discrete Algorithms, pp. 803-812.
overapproximate(Z::Zonotope, ::Type{<:Hyperrectangle})::Hyperrectangle
Return a tight overapproximation of a zonotope with an axis-aligned box.
Input
Z
– zonotopeHyperrectangle
– type for dispatch
Output
A hyperrectangle.
Algorithm
This function implements the method in [Section 5.1.2, 1]. A zonotope $Z = ⟨c, G⟩$ can be overapproximated tightly by an axis-aligned box (i.e. a Hyperrectangle
) such that its center is $c$ and the radius along dimension $i$ is the column-sum of the absolute values of the $i$-th row of $G$ for $i = 1,…, p$, where $p$ is the number of generators of $Z$.
[1] Althoff, M., Stursberg, O., & Buss, M. (2010). Computing reachable sets of hybrid systems using a combination of zonotopes and polytopes. Nonlinear analysis: hybrid systems, 4(2), 233-249.
overapproximate(X::LazySet{N}, dir::AbstractDirections{N})::HPolytope{N}
where {N}
Overapproximating a set with template directions.
Input
X
– setdir
– (concrete) direction representation
Output
An HPolytope
overapproximating the set X
with the directions from dir
.
overapproximate(X::LazySet{N},
dir::Type{<:AbstractDirections})::HPolytope{N} where {N}
Overapproximating a set with template directions.
Input
X
– setdir
– type of direction representation
Output
A HPolytope
overapproximating the set X
with the directions from dir
.
overapproximate(S::LazySet{N}, ::Type{Interval}) where {N<:Real}
Return the overapproximation of a real unidimensional set with an interval.
Input
S
– one-dimensional setInterval
– type for dispatch
Output
An interval.
Algorithm
The method relies on the exact conversion to Interval
. Two support function evaluations are needed in general.
overapproximate(cap::Intersection{N, <:LazySet, <:AbstractPolyhedron{N}},
dir::AbstractDirections{N};
kwargs...
) where {N<:Real}
Return the overapproximation of the intersection between a compact set and a polytope given a set of template directions.
Input
cap
– intersection of a compact set and a polytopedir
– template directionskwargs
– additional arguments that are passed to the support function algorithm
Output
A polytope in H-representation such that the normal direction of each half-space is given by an element of dir
.
Algorithm
Let di
be a direction drawn from the set of template directions dir
. Let X
be the compact set and let P
be the polytope. We overapproximate the set X ∩ H
with a polytope in constraint representation using a given set of template directions dir
.
The idea is to solve the univariate optimization problem ρ(di, X ∩ Hi)
for each half-space in the set P
and then take the minimum. This gives an overapproximation of the exact support function.
This algorithm is inspired from G. Frehse, R. Ray. Flowpipe-Guard Intersection for Reachability Computations with Support Functions.
Notes
This method relies on having available the constraints_list
of the polytope P
.
This method of overapproximations can return a non-empty set even if the original intersection is empty.
overapproximate(cap::Intersection{N, <:HalfSpace{N}, <:AbstractPolytope{N}},
dir::AbstractDirections{N};
[kwargs]...
) where {N<:Real}
Return the overapproximation of the intersection between a half-space and a polytope given a set of template directions.
Input
cap
– intersection of a half-space and a polytopedir
– template directionskwargs
– additional arguments that are passed to the support function algorithm
Output
A polytope in H-representation such that the normal direction of each half-space is given by an element of dir
.
Box Approximations
LazySets.Approximations.ballinf_approximation
— Function.ballinf_approximation(S::LazySet{N};
)::BallInf{N} where {N<:Real}
Overapproximate a convex set by a tight ball in the infinity norm.
Input
S
– convex set
Output
A tight ball in the infinity norm.
Algorithm
The center and radius of the box are obtained by evaluating the support function of the given convex set along the canonical directions.
LazySets.Approximations.box_approximation
— Function.box_approximation(S::LazySet{N})::Union{Hyperrectangle{N}, EmptySet{N}}
where {N<:Real}
Overapproximate a convex set by a tight hyperrectangle.
Input
S
– convex set
Output
A tight hyperrectangle.
Algorithm
The center of the hyperrectangle is obtained by averaging the support function of the given set in the canonical directions, and the lengths of the sides can be recovered from the distance among support functions in the same directions.
LazySets.Approximations.interval_hull
— Function.interval_hull
Alias for box_approximation
.
box_approximation_symmetric(S::LazySet{N}
)::Union{Hyperrectangle{N}, EmptySet{N}}
where {N<:Real}
Overapproximate a convex set by a tight hyperrectangle centered in the origin.
Input
S
– convex set
Output
A tight hyperrectangle centered in the origin.
Algorithm
The center of the box is the origin, and the radius is obtained by computing the maximum value of the support function evaluated at the canonical directions.
LazySets.Approximations.symmetric_interval_hull
— Function.symmetric_interval_hull
Alias for box_approximation_symmetric
.
LazySets.Approximations.box_approximation_helper
— Function.box_approximation_helper(S::LazySet{N};
) where {N<:Real}
Common code of box_approximation
and box_approximation_symmetric
.
Input
S
– convex set
Output
A tuple containing the data that is needed to construct a tightly overapproximating hyperrectangle.
c
– centerr
– radius
Algorithm
The center of the hyperrectangle is obtained by averaging the support function of the given convex set in the canonical directions. The lengths of the sides can be recovered from the distance among support functions in the same directions.
Iterative refinement
LocalApproximation{N<:Real}
Type that represents a local approximation in 2D.
Fields
p1
– first inner pointd1
– first directionp2
– second inner pointd2
– second directionq
– intersection of the lines l1 ⟂ d1 at p1 and l2 ⟂ d2 at p2refinable
– states if this approximation is refinableerr
– error upper bound
Notes
The criteria for being refinable are determined in the method new_approx
.
PolygonalOverapproximation{N<:Real}
Type that represents the polygonal approximation of a convex set.
Fields
S
– convex setapprox_stack
– stack of local approximations that still need to be examinedconstraints
– vector of linear constraints that are already finalized (i.e., they satisfy the given error bound)
LazySets.Approximations.new_approx
— Method.new_approx(S::LazySet, p1::Vector{N}, d1::Vector{N}, p2::Vector{N},
d2::Vector{N}) where {N<:Real}
Create a LocalApproximation
instance for the given excerpt of a polygonal approximation.
Input
S
– convex setp1
– first inner pointd1
– first directionp2
– second inner pointd2
– second direction
Output
A local approximation of S
in the given directions.
addapproximation!(Ω::PolygonalOverapproximation, p1::Vector{N},
d1::Vector{N}, p2::Vector{N}, d2::Vector{N}) where {N<:Real}
Input
Ω
– polygonal overapproximation of a convex setp1
– first inner pointd1
– first directionp2
– second inner pointd2
– second direction
Output
The list of local approximations in Ω
of the set Ω.S
is updated in-place and the new approximation is returned by this function.
LazySets.Approximations.refine
— Method.refine(approx::LocalApproximation, S::LazySet
)::Tuple{LocalApproximation, LocalApproximation}
Refine a given local approximation of the polygonal approximation of a convex set by splitting along the normal direction of the approximation.
Input
approx
– local approximation to be refinedS
– 2D convex set
Output
The tuple consisting of the refined right and left local approximations.
LazySets.Approximations.tohrep
— Method.tohrep(Ω::PolygonalOverapproximation{N})::AbstractHPolygon{N}
where {N<:Real}
Convert a polygonal overapproximation into a concrete polygon.
Input
Ω
– polygonal overapproximation of a convex set
Output
A polygon in constraint representation.
Algorithm
Internally we keep the constraints sorted. Hence we do not need to use addconstraint!
when creating the HPolygon
.
LazySets.Approximations.approximate
— Method.approximate(S::LazySet{N},
ε::N)::PolygonalOverapproximation{N} where {N<:Real}
Return an ε-close approximation of the given 2D convex set (in terms of Hausdorff distance) as an inner and an outer approximation composed by sorted local Approximation2D
.
Input
S
– 2D convex setε
– error bound
Output
An ε-close approximation of the given 2D convex set.
LazySets.Approximations.constraint
— Method.constraint(approx::LocalApproximation)
Convert a local approximation to a linear constraint.
Input
approx
– local approximation
Output
A linear constraint.
Template directions
AbstractDirections{N}
Abstract type for template direction representations.
Notes
All subtypes should implement the standard iterator methods from Base
and the function dim(d<:AbstractDirections)::Int
.
BoxDirections{N} <: AbstractDirections{N}
Box direction representation.
Fields
n
– dimension
OctDirections{N} <: AbstractDirections{N}
Octagon direction representation.
Fields
n
– dimension
Notes
Octagon directions consist of all vectors that are zero almost everywhere except in two dimensions $i$, $j$ (possibly $i = j$) where it is $±1$.
BoxDiagDirections{N} <: AbstractDirections{N}
Box-diagonal direction representation.
Fields
n
– dimension
Notes
Box-diagonal directions can be seen as the union of diagonal directions (all entries are ±1) and box directions (one entry is ±1, all other entries are 0). The iterator first enumerates all diagonal directions, and then all box directions.
PolarDirections{N<:AbstractFloat} <: AbstractDirections{N}
Polar directions representation.
Fields
Nφ
– length of the partition of the polar angle
Notes
The PolarDirections
constructor provides a sample of the unit sphere in $\mathbb{R}^2$, which is parameterized by the polar angles $φ ∈ Dφ := [0, 2π]$ respectively; see the wikipedia entry Polar coordinate system. The domain $Dφ$ is discretized in $Nφ$ pieces. Then the Cartesian components of each direction are obtained with
Examples
The integer passed as an argument is used to discretize $φ$:
julia> using LazySets.Approximations: PolarDirections
julia> pd = PolarDirections(2)
PolarDirections{Float64}(2, Array{Float64,1}[[1.0, 0.0], [-1.0, 1.22465e-16]])
julia> pd.Nφ
2
SphericalDirections{N<:AbstractFloat} <: AbstractDirections{N}
Spherical directions representation.
Fields
Nθ
– length of the partition of the azimuthal angleNφ
– length of the partition of the polar anglestack
– list of computed directions
Notes
The SphericalDirections
constructor provides a sample of the unit sphere in $\mathbb{R}^3$, which is parameterized by the azimuthal and polar angles $θ ∈ Dθ := [0, π]$ and $φ ∈ Dφ := [0, 2π]$ respectively, see the wikipedia entry Spherical coordinate system. The domains $Dθ$ and $Dφ$ are discretized in $Nθ$ and $Nφ$ respectively. Then the Cartesian componentes of each direction are obtained with
The north and south poles are treated separately so that those points are not considered more than once.
Examples
A SphericalDirections
can be built in different ways. If you pass only one integer, it is used to discretize both $θ$ and $φ$:
julia> using LazySets.Approximations: SphericalDirections
julia> sd = SphericalDirections(3)
SphericalDirections{Float64}(3, 3, Array{Float64,1}[[0.0, 0.0, 1.0], [0.0, 0.0, -1.0], [1.0, 0.0, 6.12323e-17], [-1.0, 1.22465e-16, 6.12323e-17]])
julia> sd.Nθ, sd.Nφ
(3, 3)
Pass two integers to control the discretization in $θ$ and in $φ$ separately:
julia> sd_4_5 = SphericalDirections(4, 5);
julia> length(sd_4_5)
10
julia> sd_4_8 = SphericalDirections(4, 8);
julia> length(sd_4_8)
16
See also overapproximate(X::LazySet, dir::AbstractDirections)::HPolytope
.