Library Reference
CUDAExt.GPUSampledBoxMap
GAIO.Box
GAIO.BoxMap
GAIO.BoxMeasure
GAIO.BoxPartition
GAIO.BoxSet
GAIO.IntervalBoxMap
GAIO.Node
GAIO.SampledBoxMap
GAIO.TransferOperator
GAIO.TreePartition
SIMDExt.CPUSampledBoxMap
Arpack.eigs
Arpack.eigs
Arpack.eigs
Arpack.svds
Arpack.svds
Arpack.svds
Base.:∘
Base.sum
GAIO.AdaptiveBoxMap
GAIO.GridBoxMap
GAIO.GridBoxMap
GAIO.GridBoxMap
GAIO.MonteCarloBoxMap
GAIO.MonteCarloBoxMap
GAIO.MonteCarloBoxMap
GAIO.PointDiscretizedBoxMap
GAIO.PointDiscretizedBoxMap
GAIO.PointDiscretizedBoxMap
GAIO.adaptive_newton_step
GAIO.approx_lipschitz
GAIO.armijo_rule
GAIO.bounded_point_to_key
GAIO.box_dimension
GAIO.center
GAIO.center
GAIO.chain_recurrent_set
GAIO.cover
GAIO.cover_manifold
GAIO.cover_roots
GAIO.density
GAIO.depth
GAIO.expon
GAIO.find_at_depth
GAIO.finite_time_lyapunov_exponents
GAIO.finite_time_lyapunov_exponents
GAIO.fixqr!
GAIO.hidden_keys
GAIO.index_pair
GAIO.index_quad
GAIO.index_to_key
GAIO.key_to_box
GAIO.key_to_index
GAIO.leaves
GAIO.linreg
GAIO.marginal
GAIO.marginal
GAIO.marginal
GAIO.maximal_backward_invariant_set
GAIO.maximal_forward_invariant_set
GAIO.maximal_invariant_set
GAIO.morse_adjacencies
GAIO.morse_adjacencies_and_tiles
GAIO.morse_component_map
GAIO.morse_graph
GAIO.morse_map
GAIO.morse_tiles
GAIO.neighborhood
GAIO.nth_iterate_jacobian
GAIO.point_to_box
GAIO.point_to_key
GAIO.preimage
GAIO.preimage
GAIO.radius
GAIO.rescale
GAIO.rescale
GAIO.rescale
GAIO.rk4
GAIO.rk4_flow_map
GAIO.sample_adaptive
GAIO.seba
GAIO.seba
GAIO.subdivide
GAIO.subdivide
GAIO.subdivide!
GAIO.symmetric_image
GAIO.tree_search
GAIO.unstable_set
GAIO.vertices
GAIO.volume
MakieExt.plotboxes
GAIO.@save
Exported functions
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.BoxMap
— MethodBoxMap(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.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).
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.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.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.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.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.
.
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...
.
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.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.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
.
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
.
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.AdaptiveBoxMap
— MethodBoxMap(: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.GridBoxMap
— MethodBoxMap(: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.
GAIO.MonteCarloBoxMap
— MethodBoxMap(:montecarlo, map, domain::Box{N}; n_points=16*N) -> SampledBoxMap
Construct a SampledBoxMap
that uses n_points
Monte-Carlo test points.
GAIO.PointDiscretizedBoxMap
— MethodBoxMap(: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
.
GAIO.adaptive_newton_step
— Functionadaptive_newton_step(g, g_jacobian, x, k=1)
Return one step of the adaptive Newton algorithm for the point x
.
GAIO.approx_lipschitz
— Methodapprox_lipschitz(f, center::SVector, radius::SVector, accel=nothing) -> Matrix
Compute an approximation of the Lipschitz matrix, i.e. the matrix that satisifies
\[| f(x) - f(y) | \leq L | x - y | \quad \forall \, x,y \in \text{Box(center, radius)}\]
componentwise.
GAIO.armijo_rule
— Functionarmijo_rule(g, Dg, x, d, σ=1e-4, ρ=0.8, α₀=0.05, α₁=1.0)
Find a step size multiplier $\alpha \in (\alpha_0, \alpha_1]$ such that
\[g(x + \alpha d) - g(x) \leq \alpha \sigma \, Dg(x) \cdot d\]
This is done by initializing $\alpha = 1$ and testing the above condition. If it is not satisfied, scale $\alpha$ by some constant $\rho < 1$ (i.e. set $\alpha = \rho \cdot \alpha$), and test the condition again.
GAIO.bounded_point_to_key
— Methodbounded_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.box_dimension
— Methodbox_dimension(boxsets) -> D
For an iterator boxsets
of (successively finer) BoxSet
s, compute the box dimension D
.
Example
# F is some BoxMap, S is some BoxSet
box_dimension( relative_attractor(F, S, steps=k) for k in 1:20 )
GAIO.center
— Methodcenter(center, radius)
Return the center of a box as an iterable. Default function for image_points
in SampledBoxMap
s.
GAIO.center
— Methodcenter(b::Box)
Return the center of a box.
GAIO.chain_recurrent_set
— Methodchain_recurrent_set(F::BoxMap, B::BoxSet; steps=12) -> BoxSet
Compute the chain recurrent set over the box set B
. B
should be a (coarse) covering of the relative attractor, e.g. B = cover(P, :)
for a partition P
.
GAIO.cover
— FunctionBoxSet
constructors:- set of all boxes in partition / box set
P
:
julia B = cover(P, :)
- cover the point
x
, or pointsx = [x_1, x_2, x_3] # etc ...
using boxes fromP
julia B = cover(P, x)
- a covering of
S
using boxes fromP
julia S = [Box(center_1, radius_1), Box(center_2, radius_2), Box(center_3, radius_3)] # etc... B = cover(P, S)
- set of all boxes in partition / box set
Return a subset of the partition or box set P
based on the second argument.
GAIO.cover_manifold
— Methodcover_manifold(f, B::BoxSet; steps=12)
Use interval arithmetic to compute a covering of an implicitly defined manifold $M$ of the form
\[f(M) \equiv 0\]
for some function $f : \mathbb{R}^N \to \mathbb{R}$.
The starting BoxSet B
should (coarsely) cover the manifold.
GAIO.cover_roots
— Methodcover_roots(g, Dg, B::BoxSet; steps=12) -> BoxSet
Compute a covering of the roots of g
within the partition P
. Generally, B
should be a box set containing the whole partition P
, ie B = cover(P, :)
, and should contain a root of g
.
GAIO.density
— Methoddensity(μ::BoxMeasure) -> Function
Return the measure μ
as a callable density g
, i.e.
\[\int f(x) \, d\mu (x) = \int f(x)g(x) \, dx . \]
GAIO.depth
— Methoddepth(tree::TreePartition)
Return the depth of the tree structure.
GAIO.find_at_depth
— Methodfind_at_depth(tree, depth)
Return all node indices at a specified depth.
GAIO.finite_time_lyapunov_exponents
— Methodfinite_time_lyapunov_exponents(f, Df, μ::BoxMeasure; n=8) -> σ
Compute the Lyapunov exponents using a spatial integration method [1] based on Birkhoff's ergodic theorem. Computes
\[\sigma_j = \frac{1}{n} \int \log R_{jj}( Df^n (x) ) \, dμ (x), \quad j = 1, \ldots, d\]
with respect to an ergodic invariant measure $\mu$.
[1] Beyn, WJ., Lust, A. A hybrid method for computing Lyapunov exponents. Numer. Math. 113, 357–375 (2009). https://doi.org/10.1007/s00211-009-0236-4
GAIO.finite_time_lyapunov_exponents
— Methodfinite_time_lyapunov_exponents(F::SampledBoxMap, boxset::BoxSet) -> BoxMeasure
Compute the Finite Time Lyapunov Exponent for every box in boxset
, where F
represents a time-T
integration of some continuous dynamical system. It is assumed that all boxes in boxset
have radii of some fixed order ϵ.
hidden_keys(tree)
Return all keys within the tree, including keys not corresponding to leaf nodes.
GAIO.index_pair
— Methodindex_pair(F::BoxMap, N::BoxSet) -> (P₁, P₀)
Compute an index pair of BoxSet
s P₀ ⊆ P₁ ⊆ M where M = N ∪ nbhd(N).
GAIO.index_quad
— Methodindex_quad(F::BoxMap, N::BoxSet) -> (P₁, P₀, P̄₁, P̄₀)
Compute a tuple of index pairs such that F: (P₁, P₀) → (P̄₁, P̄₀)
GAIO.index_to_key
— Methodindex_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
.
GAIO.key_to_box
— Methodkey_to_box(P::BoxPartition, key)
Return the box associated with the index within a BoxPartition
.
GAIO.key_to_index
— Methodkey_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.leaves
— Methodleaves(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.
GAIO.marginal
— Methodmarginal(μ::BoxMeasure{Box{N}}; dim) -> BoxMeasure{Box{N-1}}
Compute the marginal distribution of μ along an axis given by its dimension dim
.
GAIO.marginal
— Methodmarginal(P::BoxPartition{N}; dim) -> BoxPartition{N-1}
Construct the projection of a BoxPartition
along an axis given by its dimension dim
.
GAIO.marginal
— Methodmarginal(B::BoxSet{Box{N}}; dim) -> BoxSet{Box{N-1}}
Construct the projection of the BoxSet
along an axis given by its dimension dim
. This means that all boxes are projected to dimension N-1. Overlapping boxes are counted only once.
GAIO.maximal_backward_invariant_set
— Functionrelative_attractor(F::BoxMap, B::BoxSet; steps=12, subdivision=true) -> BoxSet
maximal_backward_invariant_set(F::BoxMap, B::BoxSet; steps=12, subdivision=true) -> BoxSet
Compute the attractor relative to B
. B
should be a (coarse) covering of the relative attractor, e.g. B = cover(P, :)
for a partition P
.
GAIO.maximal_forward_invariant_set
— Functionmaximal_forward_invariant_set(F::BoxMap, B::BoxSet; steps=12, subdivision=true)
Compute the maximal forward invariant set contained in B
. B
should be a (coarse) covering of a forward invariant set, e.g. B = cover(P, :)
for a partition P
.
GAIO.maximal_invariant_set
— Functionmaximal_invariant_set(F::BoxMap, B::BoxSet; steps=12, subdivision=true)
Compute the maximal invariant set contained in B
. B
should be a (coarse) covering of an invariant set, e.g. B = cover(P, :)
for a partition P
.
GAIO.morse_component_map
— MethodConcatenation of the condensation map and morse map. See morse_map
GAIO.morse_map
— MethodGiven a Strong_components_output
from MatrixNetworks
(in particular the component map), compute a second map on the vertices of the condensation graph to the vertices of the morse graph. Vertices of the condensation graph which do not correspond to morse sets, get sent to the (arbitrary) vertex index 0.
origninal condensation morse
graph graph graph
* ──────┐
│
* ──────┴───────→ * ───────────────→ *
* ──────────────→ * ───┐ ┌─────→ *
│ │
* ──────────────→ * ───┴─────┼────→ 0
│
* ─────┬────────→ * ─────────┘
│
* ─────┤
│
* ─────┘
⋮ ==============⟹ ⋮ ==============⟹ ⋮
condensation morse
map map
.
GAIO.morse_tiles
— MethodGiven a transfer operator and a morse component map (see morse_component_map
), compute the boxes corresponding to the vertices of the morse graph.
GAIO.neighborhood
— Methodneighborhood(B::BoxSet) -> BoxSet
nbhd(B::BoxSet) -> BoxSet
Return a one-box wide neighborhood of a BoxSet B
.
GAIO.nth_iterate_jacobian
— Methodnth_iterate_jacobian(f, Df, x, n; return_QR=false) -> Z[, R]
Compute the Jacobian of the n
-times iterated function f ∘ f ∘ ... ∘ f
at x
using a QR iteration based on [1]. Requires an approximation Df
of the jacobian of f
, e.g. Df(x) = ForwardDiff.jacobian(f, x)
. Optionally, return the QR decomposition.
[1] Dieci, L., Russell, R. D., Van Vleck, E. S.: "On the Computation of Lyapunov Exponents for Continuous Dynamical Systems," submitted to SIAM J. Numer. Ana. (1993).
GAIO.point_to_box
— Methodpoint_to_box(P::AbstractBoxPartition, point)
Find the box within a BoxPartition
containing a point.
GAIO.point_to_key
— Methodpoint_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.preimage
— Methodpreimage(F::BoxMap, B::BoxSet, Q::BoxSet) -> BoxSet
Compute the (restricted to Q
) preimage of B
under F
, i.e.
\[F^{-1} (B) \cap Q . \]
Note that the larger $Q$ is, the more calculation time required.
GAIO.preimage
— Methodpreimage(F::BoxMap, B::BoxSet) -> BoxSet
Efficiently compute
\[F^{-1} (B) \cap B . \]
Significantly faster than calling preimage(F, B, B)
.
preimage(F, B)
computes the RESTRICTED preimage $F^{-1} (B) \cap B$, NOT the full preimage $F^{-1} (B)$.
GAIO.radius
— Methodradius(b::Box)
Return the radius of a box.
GAIO.rescale
— Methodrescale(center, radius, points)
Return an iterable which calls rescale(center, radius, point)
for each point in points
.
GAIO.rescale
— Methodrescale(points)
Return a function
(center, radius) -> rescale(center, radius, points)
Used in domain_points
for BoxMap
, PointDiscretizedMap
.
GAIO.rescale
— Methodrescale(box, point::Union{<:StaticVector{N,T}, <:NTuple{N,T}})
rescale(center, radius, point::Union{<:StaticVector{N,T}, <:NTuple{N,T}})
Scale a point
within the unit box $[-1, 1]^N$ to lie within box = Box(center, radius)
.
GAIO.rk4
— Methodrk4(f, x, τ)
Compute one step with step size τ
of the classic fourth order Runge-Kutta method.
GAIO.rk4_flow_map
— Functionrk4_flow_map(f, x, step_size=0.01, steps=20)
Perform steps
steps of the classic Runge-Kutta fourth order method, with step size step_size
.
GAIO.sample_adaptive
— Methodsample_adaptive(f, center::SVector, radius::SVector)
Create a grid of test points using the adaptive technique 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.seba
— Methodseba(V::Matrix{<:Real}, Rinit=nothing, maxiter=5000) -> S, R
Construct a sparse approximation of the basis V
, as described in [1]. Returns matrices $S$, $R$ such that
\[\frac{1}{2} \| V - SR \|_F^2 + \mu \| S \|_{1,1}\]
is minimized, where $\mu \in \mathbb{R}$, $\| \cdot \|_F$ is the Frobenuius-norm, and $\| \cdot \|_{1,1}$ is the element sum norm, and $R$ is orthogonal. See [1] for further information on the argument Rinit
, as well as a description of the algorithm.
[1] Gary Froyland, Christopher P. Rock, and Konstantinos Sakellariou. Sparse eigenbasis approximation: multiple feature extraction across spatiotemporal scales with application to coherent set identification. Communications in Nonlinear Science and Numerical Simulation, 77:81-107, 2019. https://arxiv.org/abs/1812.02787
GAIO.seba
— Methodseba(V::Vector{<:BoxMeasure}, Rinit=nothing; which=partition_disjoint, maxiter=5000) -> S, A
Construct a sparse eigenbasis approximation of V
, as described in [1]. Returns an Array
of BoxMeasure
s corresponding to the eigenbasis, as well as a maximum-likelihood BoxMeasure
that maps a box to the element of S
which has the largest value over the support.
The keyword which
is used to set the threshholding heuristic, which is used to extract a partition of the supports from the sparse basis. Builtin options are
partition_unity, partition_disjoint, partition_likelihood
which are all exported functions.
GAIO.subdivide!
— Methodsubdivide!(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.subdivide
— Methodsubdivide(B::BoxSet{<:Any,<:Any,<:TreePartition}) -> BoxSet
Bisect every box in boxset
along an axis, giving rise to a new partition of the domain, with double the amount of boxes. Axis along which to bisect depends on the depth of the nodes.
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.symmetric_image
— Methodsymmetric_image(F::BoxMap, B::BoxSet) -> BoxSet
Efficiently compute
\[F (B) \cap B \cap F^{-1} (B) . \]
Internally performs the following computation (though more efficiently)
# create a measure with support over B
μ = BoxMeasure(B)
# compute transfer weights (restricted to B)
T = TransferOperator(F, B, B)
C⁺ = BoxSet(T*μ) # support of pushforward measure
C⁻ = BoxSet(T'μ) # support of pullback measure
C = C⁺ ∩ C⁻
GAIO.tree_search
— Methodtree_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.unstable_set
— Methodunstable_set(F::BoxMap, B::BoxSet) -> BoxSet
Compute the unstable set for a box set B
. Generally, B
should be a small box surrounding a fixed point of F
. The partition must be fine enough, since no subdivision occurs in this algorithm.
GAIO.vertices
— Methodvertices(box)
Return an iterator over the vertices of a box = Box(center, radius)
.
GAIO.volume
— Methodvolume(box::Box)
Compute the volume of a box.
GAIO.@save
— Macro@save boxset prefix="./" suffix=".boxset" -> filename
@save boxset filename -> filename
Save a BoxSet
as a list of keys. The default file name is the variable name.
Note that this does not include information on the partition of the BoxSet
, just the keys.
.
@save boxmap source prefix="./" suffix=".boxmap" -> filename
@save boxmap source filename -> filename
Save a BoxMap
as a list of source-keys and their image-keys in the form
key_1 -> {image_1, image_2, image_3}
key_2 -> {image_2, image_4, image_8, image_6}
⋮
.
@save transfer_operator prefix="./" suffix=".boxmap" -> filename
@save transfer_operator filename -> filename
Save a TransferOperator
as a list of keys and their image-keys in the form
key_1 -> {image_1, image_2, image_3}
key_2 -> {image_2, image_4, image_8, image_6}
⋮
MakieExt.plotboxes
— Methodplot(boxset::BoxSet)
plot(boxmeas::BoxMeasure)
plot!(boxset::BoxSet)
plot!(boxmeas::BoxMeasure)
Plot a BoxSet
or BoxMeasure
.
Special Attributes:
projection = x -> x[1:3]
If the dimension of the system is larger than 3, use this function to project to 3-d space.
color = :red
Color used for the boxes.
colormap = :default
Colormap used for plotting BoxMeasure
s values.
marker = HyperRectangle(GeometryBasics.Vec3f0(0), GeometryBasics.Vec3f0(1))
The marker for an individual box. Only works if using Makie for plotting.
All other attributes are taken from MeshScatter.
Nonexported functions
GAIO.Node
— TypeNode structure used for TreePartition
s
Fields:
left
andright
refer to indices w.r.t.
trp.nodes
for a TreePartition
trp
.
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
.
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)$.
GAIO.expon
— Functionexpon(h, k=1, ϵ=0.2, δ=0.1)
Return a rough estimate of how many Newton steps should be taken, given a step size h.
GAIO.fixqr!
— Methodfixqr!(Q, R)
Adjust a QR-decomposition such that the R-factor has positive diagonal entries.
GAIO.linreg
— Methodlinreg(xs, ys)
Simple one-dimensional lunear regression used to approximate box dimension.
GAIO.morse_adjacencies
— MethodGiven a strong_components_output
from MatrixNetworks
(in particular the component map) as well as the morse map (see morse_map
), compute the adjacency matrix for the morse graph.
GAIO.morse_adjacencies_and_tiles
— MethodGiven a transfer operator (interpreted as a transfer graph), compute the adjacency matrix for the mose graph as well as the boxes representing the vertices for the morse graph.
GAIO.morse_graph
— Methodmorse_graph(F::BoxMap, B::BoxSet) -> MetaGraph
morse_graph(F♯::TransferOperator) -> MetaGraph
Construct the morse graph
SIMDExt.CPUSampledBoxMap
— TypeBoxMap(:cpu, map, domain; n_points) -> CPUSampledBoxMap
Transforms a $map: Q → Q$ defined on points in the domain $Q ⊂ ℝᴺ$ to a CPUSampledBoxMap
defined on Box
es.
Uses the CPU's SIMD acceleration capabilities.
By default uses a grid of sample points.
BoxMap(:sampled, :cpu, boxmap, idx_base)
Type representing a discretization of a map using sample points which are explicitly vectorized. This type performs roughly 2x as many floating point operations per second as standard SampledBoxMap
s.
Fields:
boxmap
:SampledBoxMap
with one restriction:boxmap.domain_points(c, r)
must return an iterable with eltypeSVector{N, SIMD.Vec{S,T}}
whereN
is the dimension,S
is the cpu's SIMD operation capacity, e.g.4
, andT
is the individual element type, e.g.Float64
.idx_base
:SIMD.Vec{S,Int}
which is used to transform aVector{SVector{N, SIMD.Vec{S,T}}}
into aVector{SVector{N,T}}
.
.
GAIO.GridBoxMap
— MethodBoxMap(: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.
GAIO.MonteCarloBoxMap
— MethodBoxMap(: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.
GAIO.PointDiscretizedBoxMap
— MethodBoxMap(: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
.
CUDAExt.GPUSampledBoxMap
— TypeBoxMap(:gpu, map, domain; n_points) -> GPUSampledBoxMap
Transforms a $map: Q → Q$ defined on points in the domain $Q ⊂ ℝᴺ$ to a GPUSampledBoxMap
defined on Box
es.
Uses the GPU's acceleration capabilities.
By default uses a grid of sample points.
BoxMap(:sampled, :gpu, boxmap)
Type representing a dicretization of a map using sample points, which are mapped on the gpu. This type performs orders of magnitude faster than standard SampledBoxMap
s.
GPUSampledBoxMap
makes NO use of the image_points
field in SampledBoxMap
s.
Fields:
boxmap
:SampledBoxMap
with one restriction:boxmap.image_points
will not be used.
Requires a CUDA-capable gpu.
GAIO.GridBoxMap
— MethodBoxMap(: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
— MethodBoxMap(: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.PointDiscretizedBoxMap
— MethodBoxMap(: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.