Miscellaneous Methods

This page covers specialized methods for particular types of problems or applications.

Composite Algorithms

StochasticCompositeAlgorithm - Multi-Method Solving

StochasticDiffEq.StochasticCompositeAlgorithmType
StochasticCompositeAlgorithm(algs, choice_function)

StochasticCompositeAlgorithm: Multi-Method Composite Algorithm

Composite algorithm that automatically switches between multiple SDE solvers based on problem characteristics.

Method Properties

  • Approach: Multi-method solving with automatic switching
  • Adaptivity: Changes methods during integration
  • Flexibility: Combines strengths of different algorithms

Parameters

  • algs::Tuple: Tuple of algorithms to switch between
  • choice_function::Function: Function determining which algorithm to use

When to Use

  • Problems with changing characteristics during integration
  • When different regions require different solution approaches
  • Combining methods for different regimes (e.g., stiff/nonstiff)
  • When no single method is optimal for entire domain

Choice Function

The choice_function(integrator) should return an integer indicating which algorithm from algs to use:

function choice_function(integrator)
    if stiff_region(integrator.u, integrator.t)
        return 1  # Use first algorithm (e.g., implicit)
    else
        return 2  # Use second algorithm (e.g., explicit)
    end
end

Algorithm Features

  • Automatic method switching during integration
  • Maintains continuity across method transitions
  • Combines computational efficiency with robustness
  • Can handle complex multi-scale problems

References

  • Composite algorithm methodology for SDEs

RODE Methods (Random ODEs)

RandomEM - Random Euler Method

StochasticDiffEq.RandomEMType
RandomEM()

RandomEM: Random Euler Method (RODE)

Euler method for Random Ordinary Differential Equations (RODEs) with random parameters.

Method Properties

  • Problem type: Random ODEs (RODEs)
  • Strong Order: 1.0 (for deterministic part)
  • Randomness: Handles random parameters, not Brownian motion
  • Time stepping: Fixed step size

When to Use

  • Random ODEs with random parameters but no Brownian motion
  • Uncertainty quantification with parameter randomness
  • Problems with random coefficients or initial conditions
  • Monte Carlo simulation of deterministic systems with random inputs

RODE vs SDE

  • RODE: Random parameters, deterministic evolution
  • SDE: Fixed parameters, stochastic (Brownian) evolution

References

  • Random ordinary differential equation methods

RandomHeun - Random Heun Method

StochasticDiffEq.RandomHeunType
RandomHeun()

RandomHeun: Random Heun Method (RODE)

Heun method for Random Ordinary Differential Equations with improved accuracy.

Method Properties

  • Problem type: Random ODEs (RODEs)
  • Strong Order: 2.0 (for deterministic part)
  • Randomness: Handles random parameters
  • Time stepping: Fixed step size

When to Use

  • RODEs requiring higher accuracy than RandomEM
  • When computational cost per step is acceptable
  • Random parameter problems needing second-order accuracy

References

  • Higher-order methods for random ODEs

RandomTamedEM - Tamed Random Euler

StochasticDiffEq.RandomTamedEMType
RandomTamedEM()

RandomTamedEM: Tamed Random Euler Method (RODE)

Tamed Euler method for RODEs with potentially explosive behavior.

Method Properties

  • Problem type: Random ODEs with potential blow-up
  • Approach: Taming to prevent numerical explosion
  • Stability: Enhanced stability for unstable random systems
  • Time stepping: Fixed step size with taming

When to Use

  • RODEs that may exhibit explosive growth
  • When RandomEM gives unstable or explosive solutions
  • Random systems with strong nonlinearities
  • Problems requiring enhanced numerical stability

Taming Mechanism

Applies taming technique to prevent numerical blow-up while maintaining accuracy for well-behaved solutions.

References

  • Tamed methods for random differential equations

Langevin Dynamics

BAOAB - Langevin Integrator

StochasticDiffEq.BAOABType
BAOAB(;gamma=1.0, scale_noise=true)

BAOAB: Langevin Dynamics Integrator (Specialized)

Specialized integrator for Langevin dynamics in molecular dynamics simulations, particularly effective for configurational sampling.

Method Properties

  • Problem type: Langevin dynamics (second-order SDEs)
  • Structure: Position-velocity formulation
  • Sampling: Designed for equilibrium sampling
  • Time stepping: Fixed step size
  • Conservation: Preserves equilibrium distributions

Parameters

  • gamma::Real = 1.0: Friction coefficient
  • scale_noise::Bool = true: Whether to scale noise appropriately

System Structure

Designed for Langevin systems:

du = v dt
dv = f(v,u) dt - γv dt + g(u)√(2γ) dW

where:

  • u: position coordinates
  • v: velocity coordinates
  • γ: friction coefficient
  • f(v,u): force function
  • g(u): noise scaling function

When to Use

  • Molecular dynamics simulations with Langevin thermostat
  • Configurational sampling of molecular systems
  • Equilibrium sampling from canonical ensemble
  • Second-order SDEs with damping and noise

Algorithm Features

  • BAOAB splitting: B(kick) - A(drift) - O(Ornstein-Uhlenbeck) - A(drift) - B(kick)
  • Preserves correct equilibrium distribution
  • Robust and efficient for molecular sampling
  • Well-suited for long-time integration

References

  • Leimkuhler B., Matthews C., "Robust and efficient configurational molecular sampling via Langevin dynamics", J. Chem. Phys. 138, 174102 (2013)

Predictor-Corrector Methods

PCEuler - Predictor-Corrector Euler

StochasticDiffEq.PCEulerType
PCEuler(ggprime; theta=1/2, eta=1/2)

PCEuler: Predictor-Corrector Euler Method (Nonstiff)

A predictor-corrector variant of the Euler-Maruyama method requiring analytic derivatives of the diffusion term, with adjustable implicitness parameters for drift-diffusion coupling.

Method Properties

  • Strong Order: 0.5 (in the Itô sense)
  • Weak Order: 1.0
  • Time stepping: Fixed time step only
  • Noise types: General noise with available derivative information
  • SDE interpretation: Itô only

Parameters

  • ggprime::Function: The required derivative of the diffusion term
    • For scalar problems: ggprime = g * ∂g/∂x
    • For multi-dimensional problems: ggprime_k = Σ_{j=1...M, i=1...D} g^(j)_i * ∂g^(j)_k/∂x_i
    • where g^(j) corresponds to the noise vector due to the j-th noise channel
    • Must match the in-place/out-of-place specification of the problem
  • theta::Real = 0.5: Degree of implicitness in the drift term (default: 0.5)
  • eta::Real = 0.5: Degree of implicitness in the diffusion term (default: 0.5)

When to Use

  • Problems requiring specific drift-diffusion coupling
  • When analytical ggprime function is available
  • Specialized predictor-corrector applications
  • When the derivative g*g' provides stability or accuracy benefits

Algorithm Description

The method uses a predictor-corrector approach with the specific requirement of computing the derivative of the diffusion coefficient. This additional information allows for improved handling of drift-diffusion interactions through the adjustable parameters θ and η.

Limitations

  • Requires analytical computation of ggprime (cannot be approximated)
  • Fixed time step only (no adaptive versions available)
  • Limited to Itô interpretation

References

  • Jentzen, A., Kloeden, P.E., "The numerical approximation of stochastic partial differential equations", Milan J. Math. 77, 205–244 (2009). https://doi.org/10.1007/s00032-009-0100-0
  • Originally introduced in PR #88 (commit 42e2510) by Tatsuhiro Onodera (2018)
Warning

The derivative ggprime must be computed analytically for correctness. The original paper contains a typo in the definition of ggprime - this implementation follows the corrected formulation.

Integro-Integral-Form (IIF) Methods

IIF1M, IIF2M, IIF1Mil - IIF Methods

StochasticDiffEq.IIF1MType
IIF1M(;nlsolve=NLSOLVEJL_SETUP())

IIF1M: Integrating Factor Method 1 (Semi-Linear)

First-order integrating factor method for semi-linear SDEs with stiff linear parts.

Method Properties

  • Strong Order: 1.0
  • Weak Order: 1.0
  • Time stepping: Fixed or adaptive
  • Problem type: Semi-linear SDEs with stiff linear components
  • Treatment: Exponential integrator approach

Parameters

  • nlsolve: Nonlinear solver configuration

When to Use

  • Semi-linear SDEs: du = (L*u + N(u))dt + g(u)dW where L is stiff linear operator
  • Problems amenable to integrating factor techniques
  • When exponential integrators are appropriate
  • Stiff linear parts with nonlinear perturbations

Algorithm Description

Applies integrating factor exp(L*t) to handle stiff linear part exactly while treating nonlinear parts numerically.

References

  • Integrating factor methods for stiff SDEs
StochasticDiffEq.IIF2MType
IIF2M(;nlsolve=NLSOLVEJL_SETUP())

IIF2M: Integrating Factor Method 2 (Semi-Linear)

Second-order integrating factor method for semi-linear SDEs.

Method Properties

  • Strong Order: 2.0
  • Weak Order: 2.0
  • Time stepping: Fixed or adaptive
  • Problem type: Semi-linear SDEs with stiff linear components
  • Treatment: Higher-order exponential integrator

Parameters

  • nlsolve: Nonlinear solver configuration

When to Use

  • When higher accuracy than IIF1M is needed
  • Semi-linear problems requiring second-order convergence
  • More expensive but more accurate than IIF1M

References

  • Higher-order integrating factor methods for SDEs
StochasticDiffEq.IIF1MilType
IIF1Mil(;nlsolve=NLSOLVEJL_SETUP())

IIF1Mil: Integrating Factor Milstein Method (Semi-Linear)

Integrating factor method combined with Milstein correction for semi-linear SDEs.

Method Properties

  • Strong Order: 1.0 (with Milstein correction)
  • Weak Order: 1.0
  • Time stepping: Fixed or adaptive
  • Problem type: Semi-linear SDEs with stiff linear components
  • Treatment: Exponential integrator with Milstein correction

Parameters

  • nlsolve: Nonlinear solver configuration

When to Use

  • Semi-linear SDEs requiring Milstein-type accuracy
  • When both stiff linear treatment and higher-order stochastic accuracy are needed
  • Alternative to IIF1M with enhanced stochastic treatment

References

  • Integrating factor methods with Milstein correction

Simplified Methods

SimplifiedEM - Simplified Euler-Maruyama

StochasticDiffEq.SimplifiedEMType

Kloeden, P.E., Platen, E., Numerical Solution of Stochastic Differential Equations. Springer. Berlin Heidelberg (2011)

SimplifiedEM()

SimplifiedEM: Simplified Euler-Maruyama Method (High Weak Order)

Simplified version of the Euler-Maruyama method optimized for weak convergence with reduced computational cost.

Method Properties

  • Strong Order: Not optimized for strong convergence
  • Weak Order: 1.0
  • Time stepping: Fixed step size
  • Noise types: All forms (diagonal, non-diagonal, scalar, and colored noise)
  • SDE interpretation: Itô
  • Computational cost: Reduced compared to standard EM

When to Use

  • Monte Carlo simulations where weak convergence is sufficient
  • Problems where computational efficiency is more important than strong accuracy
  • Large ensemble simulations
  • When only statistical properties (expectations, moments) are needed

Algorithm Features

  • Simplified implementation reducing computational overhead
  • Maintains weak order 1.0 convergence
  • More efficient than standard EM for weak convergence applications
  • Handles general noise structures

Weak vs Strong Convergence

  • Optimized for E[f(XT)] convergence, not pathwise |XT - X_h| convergence
  • Ideal for computing expectations and statistical properties
  • Less suitable when individual trajectory accuracy is important

References

  • Kloeden, P.E., Platen, E., "Numerical Solution of Stochastic Differential Equations", Springer. Berlin Heidelberg (2011)

When to Use Miscellaneous Methods

StochasticCompositeAlgorithm:

  • When problem characteristics change during integration
  • Combining methods for different regimes
  • Automatic method switching based on conditions

RODE Methods:

  • Random ordinary differential equations
  • Problems with random parameters but no Brownian motion
  • Uncertainty quantification applications

BAOAB:

  • Molecular dynamics simulations
  • Langevin equations with specific structure
  • When preserving equilibrium distributions is important

IIF Methods:

  • Semi-linear problems with stiff linear parts
  • Problems amenable to integrating factor techniques
  • When exponential integrators are appropriate

PCEuler:

  • Problems requiring specific drift-diffusion coupling
  • When analytical ggprime function is available
  • Specialized predictor-corrector applications

These methods serve specific niches in stochastic computation and may be optimal for particular problem structures.