API
Problems
LagrangianDescriptors.augmentprob
— Functionaugmentprob(prob::ODEProblem, M; direction::Symbol=:both)
Create an augmented ODEProblem
by aggregating the forward and backward Lagrangian descriptors to the given ODEProblem
.
More precisely, with given f = prob.f
, u0 = prob.u0
, and tspan = prob.tspan
, an augmented ODEProblem
is created that solves the forward and backward ODE and the forward and backward Lagrangian descriptors simultaneously.
The (global) Lagrangian descriptors are based on the provided local descriptor M=M(du, u, p, t)
. The (forward/backward) Lagrangian descriptor is the time-integration of the local descriptor along the (forward/backward) trajectory.
If direction=:forward
or direction=:backward
are given, then only the forward or backard problem is constructed in the augmented system. The default is direction=:both
, where both directions are considered.
augmentprob(prob::RODEProblem, M; direction::Symbol=:both)
Create an augmented RODEProblem
by aggregating the forward and backward Lagrangian descriptors to the given RODEProblem
.
More precisely, with given f = prob.f
, u0 = prob.u0
, W = prob.noise
, and tspan = prob.tspan
, an augmented RODEProblem
is created that solves the forward and backward RODE and the forward and backward Lagrangian descriptors simultaneously.
The (global) Lagrangian descriptors are based on the provided local descriptor M=M(du, u, p, t, W)
. The (forward/backward) Lagrangian descriptor is the time-integration of the local descriptor along the (forward/backward) trajectory.
If direction=:forward
or direction=:backward
are given, then only the forward or backard problem is constructed in the augmented system. The default is direction=:both
, where both directions are considered.
LagrangianDescriptors.LagrangianDescriptorProblem
— TypeLagrangianDescriptorProblem
Defines a Lagrangian descriptor problem associated with a SciML differential equations problem.
Constructor
LagrangianDescriptorProblem(prob, M, uu0; direction=:both)
LagrangianDescriptorProblem
can be constructed by passing a differential equation problem (currently only ODEProblem
, but more problem types will be added), an infinitesimal Lagrangian descriptor with arguments compatible with the differential equation problem, an array of initial conditions, and, optionally, the direction of the flow.
Arguments
prob
: the differential equation problem (e.g.ODEProblem
,SDEProblem
,RODEProblem
, etc.).M
: infinitesimal Lagrangian descriptor (e.g.M=M(du, u, t, p)
for an ODEProblem).uu0
: collection of initial conditions.direction
: the direction of the flow, with default:both
, but also accepting:forward
and:backward
.
Fields
With the given arguments, the constructor for LagrangianDescriptorProblem
returns a type with the following arguments:
ensprob::T1
: a suitable ensemble problem to be solved with the given collection of initial conditionsuu0
for each solve, with a suitableprob_func
to iterate through the collection and a suitableoutput_func
to only collect the Lagrangian descriptors at the end of the time interval.uu0::T2
: the given collection of initial conditions.direction::T3
: the given or the default direction of the flow.
Example Problem
Here we apply the Lagrangian descriptor method to a periodically-forced Duffing equation.
using OrdinaryDiffEq
using LinearAlgebra: norm
using LagrangianDescriptors
function f!(du, u, p, t)
x, y = u
A, ω = p
du[1] = y
du[2] = x - x^3 + A * cos(ω * t)
end
u0 = [0.5, 2.2]
tspan = (0.0, 13.0)
A = 0.3; ω = π; p = (A, ω)
prob = ODEProblem(f!, u0, tspan, p)
M(du, u, p, t) = norm(du)
uu0 = [[x, y] for y in range(-1.0, 1.0, length=301), x in range(-1.8, 1.8, length=301)]
lagprob = LagrangianDescriptorProblem(prob, M, uu0)
lagsol = solve(lagprob, Tsit5())
plot(lagsol)
Solutions
LagrangianDescriptors.LagrangianDescriptorSolution
— Typestruct LagrangianDescriptorSolution{T1, T2, T3} enssol::T1 uu0::T2 direction::T3 end
Representation of the solution to a LagrangianDescriptorProblem
.
Fields
enssol
: theEnsembleSolution
of the associated ensemble problem inLagrangianDescriptorProblem
.uu0
: the collection of initial conditions given in
the LagrangianDescriptorProblem
.
direction:
the direction given in theLagrangianDescriptorProblem
.
CommonSolve.solve
— Functionsolve(prob::LagrangianDescriptorProblem, alg, args...; kwargs...)
Solve a LagrangianDescriptorProblem
, which amounts to solving the associated EnsembleProblem
in prob.ensprob
and returning a LagrangianDescriptorSolution
.
You should provide the necessary args
and the desired kwargs
for solving the associated ensemble problem for the underlying Differential Equation problem.
Post-processing
LagrangianDescriptors.lagrangian_descriptor
— Functionlagrangian_descriptor(sol::ODESolution, M)
Computes the Lagrangian descriptor of a single trajectory of an ODEProblem
.
It takes the provided infinitesimal Lagrangian descriptor of the form M=M(du, u, p, t)
, along with the time interval given in sol.prob.tspan
, to compute
\[L = \int_{t_0}^{t_f} M(du(t), u(t), p, t) \;\mathrm{d}t,\]
where $t_0$ and $t_f$ are the initial and final times in sol.prob.tspan
, u(t) = sol(t)
, p = sol.prob.p
, and du(t)
is computed from sol
in a way dependent on whether sol.prob
is in place or out of place.
The time integration is done via QuadGK.quadgk
.
lagrangian_descriptor(sol::RODESolution, M)
Computes the Lagrangian descriptor of a single trajectory of a RODEProblem
.
It takes the provided infinitesimal Lagrangian descriptor in the form M=M(du, u, p, t, W)
, along with the time interval given in sol.prob.tspan
, to compute
\[L = \int_{t_0}^{t_f} M(du(t), u(t), p, t, W(t)) \;\mathrm{d}t,\]
where $t_0$ and $t_f$ are the initial and final times in sol.prob.tspan
, u(t) = sol(t)
, p = sol.prob.p
, W(t) = first(sol.W(t))
, and du(t)
is either du = sol.prob.f(sol(t), sol.prob.p, t, first(sol.W(t)))
or given by sol.prob.f(du, sol(t), sol.prob.p, t, first(sol.W(t)))
, depending on whether sol.prob
is in place or out of place.
QuadGK.quadgk
is an adaptive integration method and does not perform well on solutions of random or stochastic equations, so the integration, in these case, is performed via a simple trapezoidal rule on the mesh given by sol.t
, which may be far com accurate depending on the mesh.