Interface
High-level API
#
LyapunovExponents.ContinuousLEProblem
— Type.
ContinuousLEProblem(phase_dynamics!, u0, tspan [, p [, num_attr]];
<keyword arguments>)
This is a short-hand notation for:
LEProblem(ODEProblem(...) [, num_attr]; ...)
If tspan
is a Real
instead of a Tuple
, then (0, tspan)
is passed as the tspan
argument of ODEProblem
.
For the list of usable keyword arguments, see LEProblem
.
#
LyapunovExponents.DiscreteLEProblem
— Type.
DiscreteLEProblem(phase_dynamics!, u0, tspan [, p [, num_attr]];
<keyword arguments>)
This is a short-hand notation for:
LEProblem(DiscreteProblem(...) [, num_attr]; ...)
If tspan
is a Integer
instead of a Tuple
, then (0, tspan)
is passed as the tspan
argument of DiscreteProblem
.
For the list of usable keyword arguments, see LEProblem
.
#
LyapunovExponents.lyapunov_exponents
— Function.
lyapunov_exponents(solver)
Get the result of Lyapunov exponents calculation stored in solver
.
lyapunov_exponents(phase_dynamics!, u0, tspan; <keyword arguments>)
Calculate Lyapunov exponents of a dynamical system.
Low-level API
#
LyapunovExponents.AbstractLEProblem
— Type.
The Problem type represents the setup of the Lyapunov exponents calculation.
#
LyapunovExponents.AbstractRelaxer
— Type.
The Relaxer type represents the calculation required for throwing away the transient part of the dynamics.
#
LyapunovExponents.AbstractLESolver
— Type.
The Solver type represents the core Lyapunov Exponents (LE) calculation. The LE calculation is done by calling the in-place mutating function solve!
.
The methods connecting the three principal types (Problem, Relaxer and Solver) for the LE calculation are shown in the following diagram:
┌─ Problem (AbstractLEProblem
)
│ │
│ │ get_relaxer
, relaxed
│ ▼
│ Relaxer (AbstractRelaxer
) ┄┄ ⟲ relax!
│ │
│ │ init
│ │
│┄┄┄┄┄┄┄ init
, solve
│ │
│ ▼
└▶ Solver (AbstractLESolver
) ┄┄ ⟲ solve!
, step!
#
LyapunovExponents.LEProblem
— Type.
LEProblem(phase_prob, num_attr; <keyword arguments>)
LEProblem(phase_prob; num_attr, <keyword arguments>)
Arguments
phase_prob
: Phase space dynamics represented in the form ofODEProblem
orDiscreteProblem
from DifferentialEquations.jl.phase_prob.tspan
represents the inter-orthonormalization-interval.num_attr::Integer
: Number of orthonormalizations (or some kind of "renormalisation") for calculating Lyapunov Exponents. In general, this is the number of points (considered to be) on the attractor used for the solver. The simulated time of the system for this calculation is given bynum_attr * (tspan[1] - tspan[0])
. This argument is always required and can be given as positional or keyword argument.num_tran::Integer
: Number of iterations to through away to get rid of the transient dynamics.dim_lyap::Integer
: Number of Lyapunov exponents to be calculated. Default to the full system dimension.Q0::Array
: The initial guess of the Gram-Schmidt "Lyapunov vectors".tangent_dynamics::Function
: A vector field for solving phase space evolution and tangent space evolution together. If this is not provided,tangent_dynamics
is derived fromphase_prob.f
. See alsoPhaseTangentDynamics
.
#
LyapunovExponents.LESolver
— Type.
LESolver(integrator; <keyword arguments>)
A type representing the main calculation of Lyapunov Exponents (LE). This struct holds all temporary state required for LE calculation.
#
LyapunovExponents.get_relaxer
— Function.
get_relaxer(prob::AbstractLEProblem; <keyword arguments>) :: AbstractRelaxer
Get a relaxer for a LE problem.
#
LyapunovExponents.relaxed
— Function.
relaxed(prob::AbstractLEProblem; <keyword arguments>) :: AbstractRelaxer
Throwaway the transient part of the phase space dynamics of the LE problem prob
.
That is to say, convert a LE problem (AbstractLEProblem
) to a relaxer (AbstractRelaxer
) and then call relax!
.
#
LyapunovExponents.relax!
— Function.
relax!(relaxer::AbstractRelaxer; <keyword arguments>)
Throwaway the transient part of the phase space dynamics of the LE problem prob
.
#
DiffEqBase.init
— Function.
init(prob::AbstractLEProblem; <keyword arguments>) :: AbstractLESolver
init(relaxer::AbstractRelaxer; <keyword arguments>) :: AbstractLESolver
Run phase space simulation to throw away the transient and then construct a LE solver.
#
DiffEqBase.step!
— Function.
step!(solver::AbstractLESolver)
Evolve the dynamics and then do an orthonormalization.
#
DiffEqBase.solve!
— Function.
solve!(solver::AbstractLESolver; <keyword arguments>)
Do solver.num_attr
times of orthonormalization step!(solver)
.
solve!(demo::LEDemo; progress=-1, <keyword arguments>)
Initialize demo.solver
from demo.prob
and run solve!(demo.solver)
to calculate the Lyapunov exponents.
#
DiffEqBase.solve
— Function.
solve(prob::AbstractLEProblem; <keyword arguments>)
:: AbstractLESolver
Initialize the solver (init
) and then go through the LE calculation (solve!
).
#
LyapunovExponents.PhaseTangentDynamics
— Type.
Auto-generated dynamics for solving phase and tangent dynamics together.
API for CLV solver
#
LyapunovExponents.CovariantVectors.CLVProblem
— Type.
CLVProblem(phase_prob, num_clv; <keyword arguments>)
CLVProblem(le_prob::LEProblem; <keyword arguments>)
Covariant Lyapunov vector (CLV) problem. This is a struc that holds the dynamical system definition (phase_prob
and tangent_dynamics!
) and the configuration parameters for the algorithm (num_clv
, etc.).
The CLVs are calculated using the 'dynamical' algorithm proposed by Ginelli et al. (2007, 2013).
Arguments
num_clv::Int
: Number of points at which CLV are sampled. It is0.8 * le_prob.num_attr
when constructed fromLEProblem
.num_forward_tran::Int
,num_backward_tran::Int
: Forward and backward transient dynamics. They are0.1 * le_prob.num_attr
when constructed fromLEProblem
.- See
LEProblem
forphase_prob
,tangent_dynamics!
andQ0
.
Examples (in the online documentation)
See: https://tkf.github.io/LyapunovExponents.jl/latest/gallery/
Reference
- Ginelli, F., Poggi, P., Turchi, A., Chaté, H., Livi, R., & Politi, A. (2007). Characterizing Dynamics with Covariant Lyapunov Vectors. Physical Review Letters, 99(13), 1–4. http://doi.org/10.1103/PhysRevLett.99.130601
- Ginelli, F., Chaté, H., Livi, R., & Politi, A. (2013). Covariant Lyapunov vectors. Journal of Physics A: Mathematical and Theoretical, 46(25), 254005. http://doi.org/10.1088/1751-8113/46/25/254005
#
LyapunovExponents.CovariantVectors.CLV
— Module.
Convenience methods for accessing CLV matrices.
Methods CLV.C
, CLV.R
and CLV.G
are the accessor to the matrices $C$, $R$ and $G$ in
(Eq. 32, Ginelli et al., 2013).
#
LyapunovExponents.CovariantVectors.CLV.M
— Function.
$M_n$ (the cocycle)
#
LyapunovExponents.CovariantVectors.CLV.G
— Function.
$G_n$ (Q
from the QR decomposition of $M_{k,n-k} G_{n-k}$)
#
LyapunovExponents.CovariantVectors.CLV.R
— Function.
$R_n$ (R
from the QR decomposition of $M_{k,n-k} G_{n-k}$)
#
LyapunovExponents.CovariantVectors.CLV.R_prev
— Function.
$R_{n-k}$
#
LyapunovExponents.CovariantVectors.CLV.C
— Function.
$C_n$
#
LyapunovExponents.CovariantVectors.CLV.D
— Function.
$D_n$
#
LyapunovExponents.phase_state
— Function.
phase_state(solverish) :: Vector
Get current phase-space state stored in solverish
.
#
LyapunovExponents.CovariantVectors.forward_dynamics!
— Function.
forward_dynamics!(solver::CLVSolver) :: ForwardDynamics
Solve the CLV problem up to the forward dynamics stage and return an iterator to step through the forward dynamics. See also: backward_dynamics!
, goto!
.
#
LyapunovExponents.CovariantVectors.backward_dynamics!
— Function.
backward_dynamics!(solver::CLVSolver) :: BackwardDynamics
Solve the CLV problem up to the (final) backward dynamics stage and return an iterator to step through the backward dynamics. See also forward_dynamics!
, goto!
.
Note that finishing iteration of the returned iterator does not finalize all the solver stages (namely, recording to solver.sol
, if non-default backward_dynamics
is used). In this case, solve!(solver)
has to be called after the iteration.
Example
angles = [acos(abs(dot(C[:, 1], C[:, 2]))) * 2 / π for C
in backward_dynamics!(solver)]
#
LyapunovExponents.CovariantVectors.indexed_forward_dynamics!
— Function.
indexed_forward_dynamics!(solver::CLVSolver)
indexed_forward_dynamics!(stage::ForwardDynamics)
Just a short-hand for enumerate(forward_dynamics!(solver))
. It's for symmetry with indexed_backward_dynamics!
.
#
LyapunovExponents.CovariantVectors.indexed_backward_dynamics!
— Function.
indexed_backward_dynamics!(solver::CLVSolver)
indexed_backward_dynamics!(stage::BackwardDynamics)
It is equivalent to zip(some_counter, backward_dynamics!(solver))
where some_counter
is an iterator over integers such that the same indices returned by indexed_forward_dynamics!
indicate that those events are at the same time point of the backward and forward passes. This is useful when combining matrices CLV.G
and CLV.C
to obtain the CLV in the (original) tangent space.
For such example, see: Covariant Lyapunov vectors on the Lorenz attractor in the online manual.
Note that indexed_backward_dynamics!(solver::CLVSolver)
is equivalent to
backward = backward_dynamics!(solver)
indexed_backward_dynamics!(backward)
Separately calling backward_dynamics!
is useful when the quantities other than CLV.G
(e.g., CLV.R
) are required.
See also indexed_forward_dynamics!
.
Low-level API for CLV solver
#
LyapunovExponents.CovariantVectors.goto!
— Function.
goto!(solver::CLVSolver, stage_type::Type{T}) :: T
Solve the CLV problem up to the stage of type stage_type
and return it. Instead of directly calling goto!(solver, ForwardDynamics)
and goto!(solver, BackwardDynamics)
, use the convenience methods forward_dynamics!
and backward_dynamics!
.