Data Structures
GAIO.Box
— TypeBox{N,T}(center, radius)
Box(center, radius)
A generalized box in dimension N
with element type T
. Mathematically, this is a set
\[[center_1 - radius_1,\ center_1 + radius_1) \ \times \ \ldots \ \times \ [center_N - radius_N,\ center_N + radius_N)\]
Fields:
center
: vector where the box's center is locatedradius
: vector of radii, length of the box in each dimension
Methods implemented:
:(==), in #, etc ...
.
GAIO.volume
— Methodvolume(box::Box)
Compute the volume of a box.
GAIO.center
— Functioncenter(b::Box)
Return the center of a box.
center(center, radius)
Return the center of a box as an iterable. Default function for image_points
in SampledBoxMap
s.
GAIO.radius
— Functionradius(b::Box)
Return the radius of a box.
GAIO.BoxPartition
— TypeBoxPartition(domain::Box{N}, dims::NTuple{N,<:Integer} = ntuple(_->1, N))
Data structure to partition a domain into a dims[1] x dims[2] x ... dims[N]
equidistant box grid.
Fields:
domain
: box defining the entire domainleft
: leftmost / bottom edge of the domainscale
: 1 / diameter of each box in the new partition (componentwise)dims
: tuple, number of boxes in each dimension
Methods implemented:
:(==), ndims, size, length, keys, keytype #, etc ...
.
GAIO.point_to_key
— Functionpoint_to_key(P::BoxPartition, point)
Find the index for the box within a BoxPartition
contatining a point, or nothing
if the point does not lie in the domain.
GAIO.bounded_point_to_key
— Functionbounded_point_to_key(P::BoxPartition, point)
Find the cartesian index of the nearest box within a BoxPartition
to a point. Conicides with point_to_key
if the point lies in the partition. Default behavior is to set NaN = Inf
if NaN
s are present in point
.
GAIO.key_to_box
— Functionkey_to_box(P::BoxPartition, key)
Return the box associated with the index within a BoxPartition
.
GAIO.point_to_box
— Functionpoint_to_box(P::AbstractBoxPartition, point)
Find the box within a BoxPartition
containing a point.
GAIO.subdivide
— Methodsubdivide(P::BoxPartition, dim) -> BoxPartition
subdivide(B::BoxSet, dim) -> BoxSet
Bisect every box in the BoxPartition
or BoxSet
along the axis dim
, giving rise to a new partition of the domain, with double the amount of boxes.
GAIO.TreePartition
— TypeTreePartition(domain::Box)
Binary tree structure to partition domain
into (variably sized) boxes.
Fields:
domain
:Box
denoting the full domain.nodes
: vector ofNode
s. Each node holds two indices pointing to other nodes in the vector, or 0 if the node is a leaf.
Methods implemented:
copy, keytype, keys, subdivide #, etc...
.
GAIO.subdivide!
— Functionsubdivide!(tree::TreePartition, key::keytype(tree)) -> TreePartition
subdivide!(tree::TreePartition, depth::Integer) -> TreePartition
subdivide!(boxset::BoxSet{<:Any,<:Any,<:TreePartition}, key) -> BoxSet
subdivide!(boxset::BoxSet{<:Any,<:Any,<:TreePartition}, depth) -> BoxSet
Subdivide a TreePartition
at key
. Dimension along which the node is subdivided depends on the depth of the node.
GAIO.depth
— Functiondepth(tree::TreePartition)
Return the depth of the tree structure.
GAIO.tree_search
— Functiontree_search(tree, point, max_depth=Inf) -> key, node_idx
Find the key and correspoinding node index within the tree data structure containing a point.
GAIO.find_at_depth
— Functionfind_at_depth(tree, depth)
Return all node indices at a specified depth.
GAIO.leaves
— Functionleaves(tree, initial_node_idx=1)
Return the node indices of all leaves. Begins search at initial_node_idx
, i.e. only returns node indices of nodes below initial_node_idx
within the tree.
hidden_keys(tree)
Return all keys within the tree, including keys not corresponding to leaf nodes.
GAIO.BoxSet
— TypeBoxSet(partition, indices::AbstractSet)
Internal data structure to hold boxes within a partition.
Constructors:
- set of all boxes in partition / box set
P
:
B = cover(P, :)
- cover the point
x
, or pointsx = [x_1, x_2, x_3] # etc ...
using boxes fromP
B = cover(P, x)
- a covering of
S
using boxes fromP
S = [Box(center_1, radius_1), Box(center_2, radius_2), Box(center_3, radius_3)] # etc...
B = cover(P, S)
Fields:
partition
: the partition that the set is defined overset
: set of partition-keys corresponding to the boxes in the set
Most set operations such as
union, intersect, setdiff, symdiff, issubset, isdisjoint, issetequal, isempty, length, # etc...
are supported.
GAIO.neighborhood
— Functionneighborhood(B::BoxSet) -> BoxSet
nbhd(B::BoxSet) -> BoxSet
Return a one-box wide neighborhood of a BoxSet B
.
GAIO.BoxMap
— TypeBoxMap(map, domain) -> SampledBoxMap
Transforms a $map: Q → Q$ defined on points in the domain $Q ⊂ ℝᴺ$ to a SampledBoxMap
defined on Box
es.
By default uses adaptive test-point sampling. For SIMD- and GPU-accelerated BoxMap
s, uses a grid of test points by default.
GAIO.SampledBoxMap
— TypeBoxMap(:sampled, map, domain::Box, domain_points, image_points)
Type representing a discretization of a map using sample points.
Fields:
map
: map that defines the dynamical system.domain
: domain of the map,B
.domain_points
: the spread of test points to be mapped forward in intersection algorithms. Must have the signaturedomain_points(center, radius)
and return an iterator of points withinBox(center, radius)
.image_points
: the spread of test points for comparison in intersection algorithms. Must have the signaturedomain_points(center, radius)
and return an iterator of points withinBox(center, radius)
.
.
GAIO.AdaptiveBoxMap
— FunctionBoxMap(:adaptive, f, domain::Box) -> SampledBoxMap
Construct a SampledBoxMap
which uses sample_adaptive
to generate test points as described in
Oliver Junge. “Rigorous discretization of subdivision techniques”. In: International Conference on Differential Equations. Ed. by B. Fiedler, K. Gröger, and J. Sprekels. 1999.
GAIO.PointDiscretizedBoxMap
— FunctionBoxMap(:pointdiscretized, map, domain, points) -> SampledBoxMap
Construct a SampledBoxMap
that uses the iterator points
as test points. points
must be an array or iterator of test points within the unit cube [-1,1]^N
.
BoxMap(:pointdiscretized, :simd, map, domain, points) -> CPUSampledBoxMap
Construct a CPUSampledBoxMap
that uses the iterator points
as test points. points
must have eltype SVector{N, SIMD.Vec{S,T}}
and be within the unit cube [-1,1]^N
.
BoxMap(:pointdiscretized, :gpu, map, domain::Box{N}, points) -> GPUSampledBoxMap
Construct a GPUSampledBoxMap
that uses the Vector points
as test points. points
must be a VECTOR of test points within the unit cube [-1,1]^N
.
Requires a CUDA-capable gpu.
GAIO.GridBoxMap
— FunctionBoxMap(:grid, map, domain::Box{N}; n_points::NTuple{N} = ntuple(_->16, N)) -> SampledBoxMap
Construct a SampledBoxMap
that uses a grid of test points. The size of the grid is defined by n_points
, which is a tuple of length equal to the dimension of the domain.
BoxMap(:grid, :simd, map, domain::Box{N}; n_points::NTuple{N} = ntuple(_->16, N)) -> CPUSampledBoxMap
Construct a CPUSampledBoxMap
that uses a grid of test points. The size of the grid is defined by n_points
, which is a tuple of length equal to the dimension of the domain. The number of points is rounded up to the nearest mutiple of the cpu's SIMD capacity.
BoxMap(:grid, :gpu, map, domain::Box{N}; n_points::NTuple{N} = ntuple(_->16, N)) -> GPUSampledBoxMap
Construct a GPUSampledBoxMap
that uses a grid of test points. The size of the grid is defined by n_points
, which is a tuple of length equal to the dimension of the domain.
Requires a CUDA-capable gpu.
GAIO.MonteCarloBoxMap
— FunctionBoxMap(:montecarlo, map, domain::Box{N}; n_points=16*N) -> SampledBoxMap
Construct a SampledBoxMap
that uses n_points
Monte-Carlo test points.
BoxMap(:montecarlo, :simd, map, domain::Box{N}; n_points=16*N) -> SampledBoxMap
Construct a CPUSampledBoxMap
that uses n_points
Monte-Carlo test points. The number of points is rounded up to the nearest multiple of the cpu's SIMD capacity.
BoxMap(:montecarlo, :gpu, map, domain::Box{N}; n_points=16*N) -> GPUSampledBoxMap
Construct a GPUSampledBoxMap
that uses n_points
Monte-Carlo test points.
Requires a CUDA-capable gpu.
GAIO.IntervalBoxMap
— TypeBoxMap(:interval, map, domain::Box{N}; n_subintervals::NTuple{N} = ntuple(_->4, N)) -> IntervalBoxMap
BoxMap(:interval, map, domain::Box{N}; n_subintervals::Function) -> IntervalBoxMap
Type representing a discretization of a map using interval arithmetic to construct rigorous outer coverings of map images. n_subintervals
describes how many times a given box will be subdivided before mapping. n_subintervals
is a Function which has the signature n_subintervals(center, radius)
and returns a tuple. If a tuple is passed directly for n_subintervals
, then this is converted to a constant Function (_, _) -> n_subintervals
Fields:
map
: Map that defines the dynamical system.domain
: Domain of the map,B
.n_subintervals
: Function with the signaturen_subintervals(center, radius)
which returns a tuple describing how many times a box is subdivided in each dimension before mapping.
.
GAIO.BoxMeasure
— TypeBoxMeasure(partition, vals)
Discrete measure with support on partition.domain
, and a density with respect to the volume measure which is piecewise constant on the boxes of its support.
Implemented as a dictionary mapping partition keys to weights.
Constructors:
- BoxMeasure with constant weight 0 of Type
T
(default Float64)
supported over a BoxSet
B
:
μ = BoxMeasure(B, T)
- BoxMeasure with specified weights per key
P = B.partition
weights = Dict( key => 1 for key in keys(B) )
BoxMeasure(P, weights)
- BoxMeasure with vector of weights supportted over a
BoxSet
B
:
weights = rand(length(B))
μ = BoxMeasure(B, weights)
(Note that since Boxset
s do not have a deterministic iteration order by default, this may have unintented results. This constructor should therefore only be used with BoxSet{<:Any, <:Any, <:OrderedSet}
types)
Fields:
partition
: AnAbstractBoxPartition
whose indices are used
for vals
vals
: A dictionary whose keys are the box indices from
partition
, and whose values are the measure fo the corresponding box.
Methods implemented:
length, sum, iterate, values, isapprox, ∘, LinearAlgebra.norm, LinearAlgebra.normalize!
The p-norm of a BoxMeasure is the L^p function norm of its density (w.r.t. the Lebesgue measure).
Base.sum
— Methodsum(f, μ::BoxMeasure)
sum(f, μ::BoxMeasure, B::BoxSet)
μ(B) = sum(x->1, μ, B)
Approximate the value of
\[\int_Q f \, d\mu .\]
If a BoxSet B
is passed as the third argument, then the integration is restricted to the boxes in B
\[\int_{Q \cap \bigcup_{b \in B} b} f \, d\mu .\]
The notation μ(B)
is offered to compute $\mu (\bigcup_{b \in B} b)$.
Base.:∘
— Method∘(f, boxmeas::BoxMeasure) -> BoxMeasure
∘(boxmeas::BoxMeasure, F::BoxMap) -> BoxMeasure
Postcompose the function f
with the boxmeas
, or precompose a BoxMap F
with the boxmeas
(by applying the Koopman operator). Note that the support of BoxMeasure
must be forward-invariant under F
.
GAIO.TransferOperator
— TypeTransferOperator(map::BoxMap, domain::BoxSet)
TransferOperator(map::BoxMap, domain::BoxSet, codomain::BoxSet)
Discretization of the Perron-Frobenius operator, or transfer operator. Implemented as a sparse matrix with indices referring to two BoxSet
s: domain
and codomain
.
There exists two constructors:
- only provide a
boxmap
and adomain
. In this case, thecodomain
is generated as the image ofdomain
under theboxmap
.julia> P = BoxPartition( Box((0,0), (1,0)), (10,10) ) 10 x 10 - element BoxPartition julia> domain = BoxSet( P, Set((1,2), (2,3), (3,4)) ) 3 - element Boxset over 10 x 10 - element BoxPartition julia> T = TransferOperator(boxmap, domain) TransferOperator over [...]
- provide
domain
andcodomain
. In this case, the size of the transition matrix is given.julia> codomain = domain 3 - element Boxset over 10 x 10 - element BoxPartition julia> T = TransferOperator(boxmap, domain, codomain) TransferOperator over [...]
Fields:
mat
:SparseMatrixCSC
containing transfer weights. The indexT.mat[i,j]
represents the transfer weight FROM thej
'th box incodomain
TO thei
'th box indomain
.boxmap
:SampledBoxMap
map which dictates the transfer weights.domain
:BoxSet
which contains keys for the already calculated transfers. Effectively, these are column pointers, i.e. thej
th column ofT.mat
contains transfer weights FROM box Bj, where Bj is thej
th box ofdomain
.codomain
:BoxSet
which contains keys for the already calculated transfers. Effectively, these are row pointers, i.e. thei
th row ofT.mat
contains transfer weights TO box Bi, where Bi is thei
th box ofcodomain
.
domain -->
codomain . . . . .
| . . . . .
| . . . . .
v . . mat . .
. . . . .
. . . . .
. . . . .
. . . . .
It is important to note that TranferOperator
is only supported over the box set domain
, but if one lets a TranferOperator
act on a BoxMeasure
, e.g. by multiplication, then the domain
is extended "on the fly" to include the support of the BoxMeasure
.
Methods Implemented:
:(==), size, eltype, getindex, setindex!, SparseArrays.sparse, Arpack.eigs, LinearAlgebra.mul! #, etc ...
Implementation detail:
The reader may have noticed that the matrix representation depends on the order of boxes in support
. For this reason an OrderedSet
is used. BoxSet
s using regular Set
s will be copied and converted to OrderedSet
s.
.
Arpack.eigs
— Methodeigs(gstar::TransferOperator [; kwargs...]) -> (d,[v,],nconv,niter,nmult,resid)
Compute a set of eigenvalues d
and eigenmeasures v
of gstar
. Works with the adjoint Koopman operator as well. All keyword arguments from Arpack.eigs
can be passed. See the documentation for Arpack.eigs
.
Arpack.svds
— Methodsvds(gstar::TransferOperator [; kwargs...]) -> ([U,], σ, [V,], nconv, niter, nmult, resid)
Compute a set of
- singular values
σ
- left singular vectors
U
- right singular vectors
V
of gstar
, where U
and V
are Vector
s of BoxMeasure
s. Works with the adjoint Koopman operator as well. All keyword arguments from Arpack.svds
can be passed. See the documentation for Arpack.svds
.
GAIO.key_to_index
— Functionkey_to_index(iterable, key)
Find the index in 1..length(iterable)
which holds key
, or return nothing
. Used to enumerate BoxSet
s as $\left\{ B_1, B_2, \ldots, B_n \right\}$ in TransferOperator
, BoxGraph
.
GAIO.index_to_key
— Functionindex_to_key(iterable, i)
Return the object held in the i
th position of iterable
. Used to enumerate BoxSet
s as $\left\{ B_1, B_2, \ldots, B_n \right\}$ in TransferOperator
, BoxGraph
.