Large deviation theory

Action functionals

Freidlin-Wentzell action

CriticalTransitions.fw_actionFunction
fw_action(sys::CoupledSDEs, path, time) -> Any

Calculates the Freidlin-Wentzell action of a given path with time points time in a drift field specified by the deterministic dynamics f = dynamic_rule(sys) and (normalized) noise covariance matrix covariance_matrix(sys).

The path must be a (D x N) matrix, where D is the dimensionality of the system sys and N is the number of path points. The time array must have length N.

Returns a single number, which is the value of the action functional

$S_T[\phi_t] = \frac{1}{2} \int_0^T || \dot \phi_t - f(\phi_t) ||^2_Q \text{d}t$

where $\phi_t$ denotes the path in state space, $b$ the drift field, and $T$ the total time of the path. The subscript $Q$ refers to the generalized norm $||a||_Q^2 := \langle a, Q^{-1} b \rangle$ (see anorm). Here $Q$ is the noise covariance matrix normalized by $D/L_1(Q)$, with $L_1$ being the L1 matrix norm.

source

Geometric Freidlin-Wentzell action

CriticalTransitions.geometric_actionFunction
geometric_action(sys::CoupledSDEs, path) -> Any
geometric_action(sys::CoupledSDEs, path, arclength) -> Any

Calculates the geometric action of a given path with specified arclength for the drift field specified by the deterministic dynamics f = dynamic_rule(sys) and (normalized) noise covariance matrix covariance_matrix(sys).

For a given path $\varphi$, the geometric action $\bar S$ corresponds to the minimum of the Freidlin-Wentzell action $S_T(\varphi)$ over all travel times $T>0$, where $\varphi$ denotes the path's parameterization in physical time (see fw_action). It is given by the integral

$\bar S[\varphi] = \int_0^L \left( ||\varphi'||_Q \, ||f(\varphi)||_Q - \langle \varphi', \, f(\varphi) \rangle_Q \right) \, \text{d}s$

where $s$ is the arclength coordinate, $L$ the arclength, $f$ the drift field, and the subscript $Q$ refers to the generalized dot product $\langle a, b \rangle_Q := a^{\top} \cdot Q^{-1} b$ (see anorm). Here $Q$ is the noise covariance matrix normalized by $D/L_1(Q)$, with $L_1$ being the L1 matrix norm.

Returns the value of the geometric action $\bar S$.

source

Onsager-Machlup action

CriticalTransitions.om_actionFunction
om_action(sys::CoupledSDEs, path, time, noise_strength)

Calculates the Onsager-Machlup action of a given path with time points time for the drift field f = dynamic_rule(sys) at given noise_strength.

The path must be a (D x N) matrix, where D is the dimensionality of the system sys and N is the number of path points. The time array must have length N.

Returns a single number, which is the value of the action functional

$S^{\sigma}_T[\phi_t] = \frac{1}{2} \int_0^T \left( || \dot \phi_t - f(\phi_t) ||^2_Q + \sigma^2 \nabla \cdot f \right) \, \text{d} t$

where $\phi_t$ denotes the path in state space, $b$ the drift field, $T$ the total time of the path, and $\sigma$ the noise strength. The subscript $Q$ refers to the generalized norm $||a||_Q^2 := \langle a, Q^{-1} b \rangle$ (see anorm). Here $Q$ is the noise covariance matrix normalized by $D/L_1(Q)$, with $L_1$ being the L1 matrix norm.

source

For convenience, a general action function is available where the type of functional is set as an argument:

CriticalTransitions.actionFunction
action(
    sys::CoupledSDEs,
    path::Matrix,
    time,
    functional;
    noise_strength
) -> Any

Computes the action functional specified by functional for a given CoupledSDEs sys and path parameterized by time.

  • functional = "FW": Returns the Freidlin-Wentzell action (fw_action)
  • functional = "OM": Returns the Onsager-Machlup action (om_action)
source

Minimum action paths

We provide the following two methods to calculate instantons, or minimum action paths, between two states of a CoupledSDEs system.

Minimum action method (MAM)

Minimization of the Freidlin-Wentzell action using the L-BFGS algorithm of Optim.jl.

CriticalTransitions.min_action_methodFunction
min_action_method(
    sys::CoupledSDEs,
    x_i,
    x_f,
    N::Int64,
    T::Real;
    functional,
    maxiter,
    blocks,
    method,
    save_iterations,
    show_progress,
    verbose
) -> Any

Runs the Minimum Action Method (MAM) to find the minimum action path (instanton) between an initial state x_i and final state x_f.

This algorithm uses the minimizers of the Optim package to minimize the Freidlin-Wentzell action functional (see fw_action) for the given CoupledSDEs sys. The path is initialized as a straight line between x_i and x_f, parameterized in time via N equidistant points and total time T. Thus, the time step between discretized path points is $\Delta t = T/N$. To set an initial path different from a straight line, see the multiple dispatch method

  • min_action_method(sys::CoupledSDEs, init::Matrix, T::Real; kwargs...).

The minimization can be performed in blocks to save intermediate results.

Keyword arguments

  • functional = "FW": type of action functional to minimize. Defaults to fw_action, alternative: om_action.
  • maxiter = 100: maximum number of iterations before the algorithm stops.
  • action_tol=1e-5: relative tolerance of action value to determine convergence
  • abstol=1e-8: absolute tolerance of action gradient to determine convergence
  • reltol=1e-8: relative tolerance of action gradient to determine convergence
  • blocks = 1: number of iterative optimization blocks
  • method = LBFGS(): minimization algorithm (see Optim)
  • save_iterations = true: whether to save Optim information
  • verbose = true: whether to print Optim information during the run
  • show_progress = false: whether to print a progress bar

Output

If save_iterations, returns Optim.OptimizationResults. Else, returns only the optimizer (path). If blocks > 1, a vector of results/optimizers is returned.

source
min_action_method(
    sys::CoupledSDEs,
    init::Matrix,
    T::Real;
    functional,
    maxiter,
    blocks,
    method,
    action_tol,
    abstol,
    reltol,
    save_iterations,
    verbose,
    show_progress
) -> Any

Runs the Minimum Action Method (MAM) to find the minimum action path (instanton) from an initial condition init, given a system sys and total path time T.

The initial path init must be a matrix of size (D, N), where D is the dimension of the system and N is the number of path points. The physical time of the path is specified by T, such that the time step between consecutive path points is $\Delta t = T/N$.

For more information see the main method, min_action_method(sys::CoupledSDEs, x_i, x_f, N::Int, T::Real; kwargs...).

source

Geometric minimum action method (gMAM)

Minimization of the geometric action following Heymann and Vanden-Eijnden, PRL (2008). The gMAM reformulates MAM to avoid numerical stiffness by reparametrizing the path in terms of arc length. This ensures an even distribution of points along the path, enhancing numerical stability and accuracy. By solving the reparametrized integral, gMAM accurately captures the geometry of the action functional and is well-suited for systems with complex energy landscapes or intricate dynamics.

CriticalTransitions.geometric_min_action_methodFunction
geometric_min_action_method(
    sys::CoupledSDEs,
    x_i,
    x_f;
    N,
    kwargs...
) -> Tuple{Vector{Matrix}, Vector{Float64}}

Computes the minimizer of the geometric Freidlin-Wentzell action based on the geometric minimum action method (gMAM), using optimizers of Optim.jl or the original formulation by Heymann and Vanden-Eijnden[1].

To set an initial path different from a straight line, see the multiple dispatch method

  • geometric_min_action_method(sys::CoupledSDEs, init::Matrix, arclength::Float64; kwargs...).

Keyword arguments

  • maxiter::Int=100: maximum number of optimization iterations before the alogrithm stops
  • action_tol=1e-5: relative tolerance of action value to determine convergence
  • abstol=1e-8: absolute tolerance of action gradient to determine convergence
  • reltol=1e-8: relative tolerance of action gradient to determine convergence
  • method=LBFGS(): optimizer method (see Optim.jl)
  • iter_per_batch=1: number of iterations per optimization batch
  • tau=0.1: parameter in HeymannVandenEijnden method
  • verbose=false: if true, print additional output
  • show_progress=true: if true, display a progress bar

Optimization algorithms

The method keyword argument takes solver methods of the Optim.jl package; alternatively, the option solver = "HeymannVandenEijnden" uses the original gMAM algorithm[1].

source
geometric_min_action_method(
    sys::CoupledSDEs,
    init::Matrix;
    maxiter,
    abstol,
    reltol,
    action_tol,
    method,
    tau,
    iter_per_batch,
    verbose,
    show_progress
) -> Tuple{Vector{Matrix}, Vector{Float64}}

Runs the geometric Minimum Action Method (gMAM) to find the minimum action path (instanton) from an initial condition init, given a system sys and total arc length arclength.

The initial path init must be a matrix of size (D, N), where D is the dimension of the system and N is the number of path points.

For more information see the main method, geometric_min_action_method(sys::CoupledSDEs, x_i, x_f, arclength::Float64; kwargs...).

source

Simple Geometric minimum action method (sgMAM)

Simplified minimization of the geometric action following Heymann and Vanden-Eijnden, PRL (2008). The sgMAM is a streamlined version of the gMAM that simplifies the computation by avoiding explicit reparametrization of the path. Instead, it introduces an implicit reparametrization by rewriting the problem as a constrained optimization.

The implementation below perform a constrained gradient descent where it assumes an autonomous system with additive Gaussian noise.

CriticalTransitions.Sgmam.SgmamSystemType

A structure representing a system with Hamiltonian functions Hx and Hp.

This system operates in an extended phase space where the Hamiltonian is assumed to be quadratic in the extended momentum. The phase space coordinates x are doubled to form a 2n-dimensional extended phase space.

source
CriticalTransitions.Sgmam.sgmamFunction
sgmam(
    sys::SgmamSystem,
    x_initial;
    ϵ,
    iterations,
    show_progress,
    reltol
) -> Tuple{Any, Float64, Any, Any, Any}

Performs the simplified geometric Minimal Action Method (sgMAM) on the given system sys.

This method computes the optimal path in the phase space of a Hamiltonian system that minimizes the action. The Hamiltonian functions H_x and H_p define the system's dynamics in a doubled phase. The initial state x_initial is evolved iteratively using constrained gradient descent with step size parameter ϵ over a specified number of iterations. The method can display a progress meter and will stop early if the relative tolerance reltol is achieved.

The function returns a tuple containing the final state, the action value, the Lagrange multipliers, the momentum, and the state derivatives.

source