OrdinaryDiffEqAdamsBashforthMoulton
Adams-Bashforth and Adams-Moulton multistep methods for non-stiff differential equations. Note that Runge-Kutta methods generally come out as more efficient in benchmarks, except when the ODE function f
is expensive to evaluate or the problem is very smooth. These methods can achieve high accuracy with fewer function evaluations per step than Runge-Kutta methods in those specific cases.
Key Properties
Adams-Bashforth-Moulton methods provide:
- Reduced function evaluations compared to Runge-Kutta methods
- High efficiency for expensive-to-evaluate functions
- Multistep structure using information from previous timesteps
- Variable step and order capabilities for adaptive integration
- Predictor-corrector variants for enhanced accuracy and stability
- Good stability properties for non-stiff problems
When to Use Adams-Bashforth-Moulton Methods
These methods are recommended for:
- Expensive function evaluations where minimizing calls to
f
is critical - Non-stiff smooth problems with regular solution behavior
- Long-time integration where efficiency over many steps matters
- Problems with expensive Jacobian computations that cannot use implicit methods efficiently
- Scientific computing applications with computationally intensive right-hand sides
- Systems where startup cost of multistep methods is amortized over long integration
Method Types
Explicit Adams-Bashforth (AB)
Pure explicit multistep methods using only past information:
- Lower computational cost per step
- Less stability than predictor-corrector variants
- Good for mildly stiff problems
Predictor-Corrector Adams-Bashforth-Moulton (ABM)
Implicit corrector step for enhanced accuracy:
- Better accuracy than pure explicit methods
- Improved stability properties
- Slightly higher cost but often worth it
Solver Selection Guide
Primary recommendation
VCABM
: Main recommendation - adaptive order variable-step Adams-Bashforth-Moulton, best overall choice for Adams methods
Variable-step predictor-corrector methods
VCABM3
: Third-order variable-step Adams-Bashforth-MoultonVCABM4
: Fourth-order variable-step Adams-Bashforth-MoultonVCABM5
: Fifth-order variable-step Adams-Bashforth-Moulton
Variable-step Adams-Bashforth methods
VCAB3
: Third-order variable-step Adams-BashforthVCAB4
: Fourth-order variable-step Adams-BashforthVCAB5
: Fifth-order variable-step Adams-Bashforth
Fixed-step predictor-corrector methods
ABM32
: Third-order Adams-Bashforth-MoultonABM43
: Fourth-order Adams-Bashforth-MoultonABM54
: Fifth-order Adams-Bashforth-Moulton
Fixed-step explicit methods
AB3
: Third-order Adams-BashforthAB4
: Fourth-order Adams-BashforthAB5
: Fifth-order Adams-Bashforth
Performance Considerations
- Most efficient when function evaluation dominates computational cost
- Startup phase requires initial steps from single-step method
- Memory efficient compared to high-order Runge-Kutta methods
- Best for smooth problems - avoid for problems with discontinuities
Installation
To be able to access the solvers in OrdinaryDiffEqAdamsBashforthMoulton
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqAdamsBashforthMoulton")
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 OrdinaryDiffEqAdamsBashforthMoulton
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, VCABM())
Full list of solvers
Explicit Multistep Methods
OrdinaryDiffEqAdamsBashforthMoulton.AB3
— TypeAB3(; thread = OrdinaryDiffEq.False())
Adams-Bashforth Explicit Method The 3-step third order multistep method. Ralston's Second Order Method is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.AB4
— TypeAB4(; thread = OrdinaryDiffEq.False())
Adams-Bashforth Explicit Method The 4-step fourth order multistep method. Runge-Kutta method of order 4 is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.AB5
— TypeAB5(; thread = OrdinaryDiffEq.False())
Adams-Bashforth Explicit Method The 5-step fifth order multistep method. Ralston's 3rd order Runge-Kutta method is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
Predictor-Corrector Methods
OrdinaryDiffEqAdamsBashforthMoulton.ABM32
— TypeABM32(; thread = OrdinaryDiffEq.False())
Adams-Bashforth Explicit Method It is third order method. In ABM32, AB3 works as predictor and Adams Moulton 2-steps method works as Corrector. Ralston's Second Order Method is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.ABM43
— TypeABM43(; thread = OrdinaryDiffEq.False())
Adams-Bashforth Explicit Method It is fourth order method. In ABM43, AB4 works as predictor and Adams Moulton 3-steps method works as Corrector. Runge-Kutta method of order 4 is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.ABM54
— TypeABM54(; thread = OrdinaryDiffEq.False())
Adams-Bashforth Explicit Method It is fifth order method. In ABM54, AB5 works as predictor and Adams Moulton 4-steps method works as Corrector. Runge-Kutta method of order 4 is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.VCAB3
— TypeVCAB3(; thread = OrdinaryDiffEq.False())
Adams explicit Method The 3rd order Adams method. Bogacki-Shampine 3/2 method is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.VCAB4
— TypeVCAB4(; thread = OrdinaryDiffEq.False())
Adams explicit Method The 4th order Adams method. Runge-Kutta 4 is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.VCAB5
— TypeVCAB5(; thread = OrdinaryDiffEq.False())
Adams explicit Method The 5th order Adams method. Runge-Kutta 4 is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.VCABM3
— TypeVCABM3(; thread = OrdinaryDiffEq.False())
Adams explicit Method The 3rd order Adams-Moulton method. Bogacki-Shampine 3/2 method is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.VCABM4
— TypeVCABM4(; thread = OrdinaryDiffEq.False())
Adams explicit Method The 4th order Adams-Moulton method. Runge-Kutta 4 is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.VCABM5
— TypeVCABM5(; thread = OrdinaryDiffEq.False())
Adams explicit Method The 5th order Adams-Moulton method. Runge-Kutta 4 is used to calculate starting values.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1
OrdinaryDiffEqAdamsBashforthMoulton.VCABM
— TypeVCABM(; thread = OrdinaryDiffEq.False())
adaptive order Adams explicit Method An adaptive order adaptive time Adams Moulton method. It uses an order adaptivity algorithm is derived from Shampine's DDEABM.
Keyword Arguments
thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
E. Hairer, S. P. Norsett, G. Wanner, Solving Ordinary Differential Equations I, Nonstiff Problems. Computational Mathematics (2nd revised ed.), Springer (1996) doi: https://doi.org/10.1007/978-3-540-78862-1