OrdinaryDiffEqBDF

IMEX BDF (Implicit-Explicit Backward Differentiation Formula) methods for stiff differential equations that can be split into stiff and non-stiff components. These methods apply implicit BDF schemes to the stiff part while treating the non-stiff part explicitly, providing efficient handling of problems with mixed stiffness characteristics.

Key Properties

IMEX BDF methods provide:

  • Implicit-explicit splitting for mixed stiffness problems
  • BDF stability for the stiff component with A-stable and L-stable behavior
  • Explicit treatment of non-stiff terms avoiding unnecessary computational cost
  • High-order accuracy up to 4th order for both components
  • Efficient for large systems where full implicit treatment is expensive
  • Natural for operator splitting problems

When to Use IMEX BDF Methods

These methods are recommended for:

  • Reaction-diffusion systems where reaction terms are stiff and diffusion is moderate
  • Convection-diffusion problems with stiff source terms and explicit convection
  • Parabolic PDEs where diffusion operators are naturally split from other terms
  • Problems with natural stiffness separation where some terms require implicit treatment
  • Large-scale systems where full implicit methods are computationally prohibitive
  • Applications requiring operator splitting methodology

Mathematical Background

IMEX BDF methods split the ODE system du/dt = f(u,t) into: du/dt = f₁(u,t) + f₂(u,t)

where:

  • f₁(u,t) contains stiff terms (treated implicitly with BDF)
  • f₂(u,t) contains non-stiff terms (treated explicitly)

This splitting must be chosen carefully to ensure both stability and efficiency.

Problem Splitting Requirements

These methods require a SplitODEProblem formulation where:

  • First functionf₁ should contain stiff, implicit terms
  • Second functionf₂ should contain non-stiff, explicit terms
  • Splitting strategy significantly affects method performance
  • Stiffness characteristics should align with implicit/explicit treatment

Solver Selection Guide

IMEX Multistep Methods

  • SBDF2: Recommended - Second-order IMEX BDF method, good balance of accuracy and stability
  • SBDF3: Third-order method for higher accuracy requirements
  • SBDF4: Fourth-order method for maximum accuracy in IMEX BDF family
  • SBDF: Adaptive order method (experimental)

IMEX SDIRK Methods

  • IMEXEuler: First-order method for simple problems or debugging
  • IMEXEulerARK: Alternative first-order formulation

Performance Guidelines

When IMEX BDF methods excel

  • Natural stiffness separation where splitting is obvious
  • Large systems where full implicit treatment is expensive
  • Parabolic PDEs with natural operator splitting
  • Reaction-diffusion problems with well-separated timescales
  • Problems where implicit component has efficient linear algebra

Splitting strategy considerations

  • Identify stiff vs non-stiff terms based on eigenvalue analysis
  • Linear stiff terms work well in implicit component
  • Nonlinear non-stiff terms are suitable for explicit treatment
  • Test different splittings to optimize performance

Alternative Approaches

Consider these alternatives:

  • Full implicit methods (BDF, SDIRK) if splitting is unclear or ineffective
  • Standard IMEX Runge-Kutta methods for different accuracy/efficiency trade-offs
  • Exponential integrators for linear stiff problems with nonlinear non-stiff terms
  • Rosenbrock methods for moderately stiff problems without natural splitting

Usage Considerations

  • Careful splitting design is crucial for method effectiveness
  • Stability analysis should verify that explicit treatment doesn't introduce instabilities
  • Timestep restrictions may apply to the explicit component
  • Linear algebra efficiency in the implicit component affects overall performance

Installation

To be able to access the solvers in OrdinaryDiffEqBDF, you must first install them use the Julia package manager:

using Pkg
Pkg.add("OrdinaryDiffEqBDF")

This will only install the solvers listed at the bottom of this page. If you want to explore other solvers for your problem, you will need to install some of the other libraries listed in the navigation bar on the left.

Example usage

using OrdinaryDiffEqBDF

function lorenz!(du, u, p, t)
    du[1] = 10.0 * (u[2] - u[1])
    du[2] = u[1] * (28.0 - u[3]) - u[2]
    du[3] = u[1] * u[2] - (8 / 3) * u[3]
end
u0 = [1.0; 0.0; 0.0]
tspan = (0.0, 100.0)
prob = ODEProblem(lorenz!, u0, tspan)
sol = solve(prob, SBDF2())

Full list of solvers

IMEX Multistep

Missing docstring.

Missing docstring for SBDF. Check Documenter's build log for details.

OrdinaryDiffEqBDF.SBDF2Function
SBDF2(;kwargs...)

The two-step version of the IMEX multistep methods of

  • Uri M. Ascher, Steven J. Ruuth, Brian T. R. Wetton. Implicit-Explicit Methods for Time-Dependent Partial Differential Equations. Society for Industrial and Applied Mathematics. Journal on Numerical Analysis, 32(3), pp 797-823, 1995. doi: https://doi.org/10.1137/0732037

See also SBDF.

source
OrdinaryDiffEqBDF.SBDF3Function
SBDF3(;kwargs...)

The three-step version of the IMEX multistep methods of

  • Uri M. Ascher, Steven J. Ruuth, Brian T. R. Wetton. Implicit-Explicit Methods for Time-Dependent Partial Differential Equations. Society for Industrial and Applied Mathematics. Journal on Numerical Analysis, 32(3), pp 797-823, 1995. doi: https://doi.org/10.1137/0732037

See also SBDF.

source
OrdinaryDiffEqBDF.SBDF4Function
SBDF4(;kwargs...)

The four-step version of the IMEX multistep methods of

  • Uri M. Ascher, Steven J. Ruuth, Brian T. R. Wetton. Implicit-Explicit Methods for Time-Dependent Partial Differential Equations. Society for Industrial and Applied Mathematics. Journal on Numerical Analysis, 32(3), pp 797-823, 1995. doi: https://doi.org/10.1137/0732037

See also SBDF.

source

IMEX SDIRK

Note that Implicit Euler is the 1st order BDF method, and is thus implemented here using the same machinery.

OrdinaryDiffEqBDF.IMEXEulerFunction
IMEXEuler(;kwargs...)

The one-step version of the IMEX multistep methods of

  • Uri M. Ascher, Steven J. Ruuth, Brian T. R. Wetton. Implicit-Explicit Methods for Time-Dependent Partial Differential Equations. Society for Industrial and Applied Mathematics. Journal on Numerical Analysis, 32(3), pp 797-823, 1995. doi: https://doi.org/10.1137/0732037

When applied to a SplitODEProblem of the form

u'(t) = f1(u) + f2(u)

The default IMEXEuler() method uses an update of the form

unew = uold + dt * (f1(unew) + f2(uold))

See also SBDF, IMEXEulerARK.

source
OrdinaryDiffEqBDF.IMEXEulerARKFunction
IMEXEulerARK(;kwargs...)

The one-step version of the IMEX multistep methods of

  • Uri M. Ascher, Steven J. Ruuth, Brian T. R. Wetton. Implicit-Explicit Methods for Time-Dependent Partial Differential Equations. Society for Industrial and Applied Mathematics. Journal on Numerical Analysis, 32(3), pp 797-823, 1995. doi: https://doi.org/10.1137/0732037

When applied to a SplitODEProblem of the form

u'(t) = f1(u) + f2(u)

A classical additive Runge-Kutta method in the sense of Araújo, Murua, Sanz-Serna (1997) consisting of the implicit and the explicit Euler method given by

y1   = uold + dt * f1(y1)
unew = uold + dt * (f1(unew) + f2(y1))

See also SBDF, IMEXEuler.

source