Interface
High-level API
#
LyapunovExponents.ContinuousLEProblem
— Type.
ContinuousLEProblem(phase_dynamics, u0 [, p];
t_attr=<number>, <keyword arguments>)
This is a short-hand notation for:
LEProblem(ODEProblem(phase_dynamics, u0 [, p]), t_attr)
For the list of usable keyword arguments, see LEProblem
.
#
LyapunovExponents.DiscreteLEProblem
— Type.
DiscreteLEProblem(phase_dynamics, u0 [, p];
t_attr=<number>, <keyword arguments>)
This is a short-hand notation for:
LEProblem(DiscreteProblem(phase_dynamics, u0 [, p]), t_attr)
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.LEProblem
— Type.
LEProblem(phase_prob, t_attr; <keyword arguments>)
Arguments
phase_prob
: Phase space dynamics represented in the form ofODEProblem
orDiscreteProblem
from DifferentialEquations.jl.phase_prob.tspan
is ignored.t_attr::Real
: Simulated time on the (presumed) attractor used for the solver. It is used for computation of Lyapunov Exponents. Roughlyt_attr / t_renorm
instantaneous exponents are sampled.t_tran::Real
: Number of iterations to throw away to get rid of the effect from the transient dynamics.t_renorm::Real
: Interval between orthonormalizations (renormalizations).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(prob::LEProblem; record::Bool = false)
Create a solver object for a LEProblem
. Record all finite-time (instantaneous) Lyapunov exponents when record = true
is passed.
#
DiffEqBase.init
— Function.
init(prob::LEProblem; <keyword arguments>) :: LESolver
init(prob::CLVProblem; <keyword arguments>) :: CLVSolver
These are simply the aliases of LESolver(prob; <keyword arguments>)
and CLVSolver(prob; <keyword arguments>)
. See LESolver
and CLVSolver
for supported keyword arguments.
#
DiffEqBase.step!
— Function.
step!(integ::DEIntegrator [, dt [, stop_at_tdt]])
Perform one (successful) step on the integrator.
Alternative, if a dt
is given, then step!
the integrator until there is a temporal difference ≥ dt
in integ.t
. When true
is passed to the optional third argument, the integrator advances exactly dt
.
step!(stage::AbstractRenormalizer)
Evolve the dynamics and then do an orthonormalization.
#
DiffEqBase.solve!
— Function.
solve!(solver::LESolver; <keyword arguments>) :: LESolver
solve!(solver::CLVSolver; <keyword arguments>) :: CLVSolver
Solve pre-initialized problem.
#
DiffEqBase.solve
— Function.
solve(prob::LEProblem; <keyword arguments>) :: LESolution
solve(prob::CLVProblem; <keyword arguments>) :: CLVSolution
Equivalent to:
solver = init(prob; <keyword arguments except progress>)
solve!(solver; progress=progress)
solver.sol
#
LyapunovExponents.PhaseTangentDynamics
— Type.
Auto-generated dynamics for solving phase and tangent dynamics together.
API for CLV solver
#
LyapunovExponents.CovariantVectors.CLVProblem
— Type.
CLVProblem(phase_prob, t_clv; <keyword arguments>)
CLVProblem(le_prob::LEProblem; <keyword arguments>)
Covariant Lyapunov vector (CLV) problem. This is a struct that holds the dynamical system definition (phase_prob
and tangent_dynamics
) and the configuration parameters for the algorithm (t_clv
, etc.).
The CLVs are calculated using the 'dynamical' algorithm proposed by Ginelli et al. (2007, 2013).
Arguments
t_clv::Real
: Time spam for which CLV are sampled. It is0.8 * le_prob.t_attr
when constructed fromLEProblem
.t_forward_tran::Real
,t_backward_tran::Real
: Forward and backward transient dynamics. They are0.1 * le_prob.t_attr
when constructed fromLEProblem
.- See
LEProblem
forphase_prob
,t_renorm
,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.CLVSolver
— Type.
CLVSolver(prob::CLVProblem; <keyword arguments>)
CLVSolver(prob::LEProblem; <keyword arguments>)
The preferred and equivalent method to get a solver for a CLVProblem
is init(prob::CLVProblem)
. Note that CLVSolver(prob::LEProblem)
is equivalent to init(CLVProblem(prob))
.
Arguments
record::Vector{Symbol}
: Variables to be saved. A subset of[:G, :C, :D, :x]
.
#
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(stage) :: Vector
Get current phase-space state stored in stage
.
#
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.Stages.goto!
— Function.
goto!(solver::StagedSolver, stage_type::Type{T}) :: T
Advance the solver
up to the stage of type stage_type
and return it.