OrdinaryDiffEqRosenbrock

Rosenbrock methods are semi-implicit Runge-Kutta methods designed for small to medium-sized stiff systems of differential equations. These methods combine the stability properties of implicit methods with computational efficiency by using approximate Jacobians. They work particularly well for strict tolerance requirements.

Key Properties

Rosenbrock methods provide:

  • Excellent efficiency for small to medium systems (<1000 ODEs)
  • L-stable and A-stable variants for stiff problems
  • W-method structure making them robust to inaccurate Jacobians
  • Automatic differentiation compatibility for Jacobian computation
  • High-order accuracy with embedded error estimation
  • Excellent performance for strict tolerance requirements
  • Stiffly accurate variants for enhanced stability

When to Use Rosenbrock Methods

Rosenbrock methods are recommended for:

  • Small to medium stiff systems (10 to 1000 equations)
  • Problems where Jacobians are available or can be computed efficiently
  • Medium tolerance requirements (1e-8 to 1e-2)
  • Stiff ODEs arising from reaction kinetics and chemical systems
  • Moderately stiff PDEs after spatial discretization
  • Problems requiring reliable error control for stiff systems

Solver Selection Guide

Low tolerance (>1e-2)

  • Rosenbrock23: Second/third-order method, recommended for low tolerance requirements

Medium tolerance (1e-8 to 1e-2)

  • Rodas5P: Fifth-order method, most efficient for many problems
  • Rodas4P: Fourth-order method, more reliable than Rodas5P
  • Rodas5Pe: Enhanced fifth-order variant with stiffly accurate embedded estimate for better adaptivity on highly stiff equations
  • Rodas5Pr: Fifth-order variant with residual test for robust error estimation that guarantees accuracy on interpolation

Performance Guidelines

  • Rodas5P: Best overall efficiency at medium tolerances
  • Rodas4P: Most reliable when Rodas5P fails
  • Rosenbrock23: Fastest at high tolerances (>1e-2)

When to Choose Alternatives

Consider other methods when:

  • System size > 1000: Use BDF methods (QNDF, FBDF)
  • Matrix-free methods: If using Krylov solvers for the linear solver (matrix-free methods), SDIRK or BDF methods are preferred

Advantages of Rosenbrock Methods

  • Very low tolerances: Rosenbrock23 performs well even at very low tolerances
  • Very high stiffness: Rosenbrock methods are often more stable than SDIRK or BDF methods because other methods can diverge due to bad initial guesses for the Newton method (leading to nonlinear solver divergence)

Installation

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

using Pkg
Pkg.add("OrdinaryDiffEqRosenbrock")

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 OrdinaryDiffEqRosenbrock

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, Rodas5P())

Full list of solvers

OrdinaryDiffEqRosenbrock.Rosenbrock23Type
Rosenbrock23(; autodiff = AutoForwardDiff(),
               concrete_jac = nothing,
               linsolve = nothing,
               step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner-W(olfbrandt) Method. An Order 2/3 L-Stable Rosenbrock-W method which is good for very stiff equations with oscillations at low tolerances. 2nd order stiff-aware interpolation.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Rosenbrock23(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Shampine L.F. and Reichelt M., (1997) The MATLAB ODE Suite, SIAM Journal of Scientific Computing, 18 (1), pp. 1-22.
source
OrdinaryDiffEqRosenbrock.Rosenbrock32Type
Rosenbrock32(; autodiff = AutoForwardDiff(),
               concrete_jac = nothing,
               linsolve = nothing,
               step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner-W(olfbrandt) Method. An Order 3/2 A-Stable Rosenbrock-W method which is good for mildly stiff equations without oscillations at low tolerances. Note that this method is prone to instability in the presence of oscillations, so use with caution. 2nd order stiff-aware interpolation.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Rosenbrock32(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Shampine L.F. and Reichelt M., (1997) The MATLAB ODE Suite, SIAM Journal of Scientific Computing, 18 (1), pp. 1-22.
source
OrdinaryDiffEqRosenbrock.ROS3PType
ROS3P(; autodiff = AutoForwardDiff(),
        concrete_jac = nothing,
        linsolve = nothing,
        step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner Method. 3rd order A-stable and stiffly stable Rosenbrock method. Keeps high accuracy on discretizations of nonlinear parabolic PDEs.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Lang, J. & Verwer, ROS3P—An Accurate Third-Order Rosenbrock Solver Designed for Parabolic Problems J. BIT Numerical Mathematics (2001) 41: 731. doi:10.1023/A:1021900219772
source
OrdinaryDiffEqRosenbrock.Rodas3Type
Rodas3(; autodiff = AutoForwardDiff(),
         concrete_jac = nothing,
         linsolve = nothing,
         step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner Method. 3rd order A-stable and stiffly stable Rosenbrock method.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Sandu, Verwer, Van Loon, Carmichael, Potra, Dabdub, Seinfeld, Benchmarking stiff ode solvers for atmospheric chemistry problems-I. implicit vs explicit, Atmospheric Environment, 31(19), 3151-3166, 1997.
source
OrdinaryDiffEqRosenbrock.Rodas23WType
Rodas23W(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner-W(olfbrandt) Method. An Order 2/3 L-Stable Rosenbrock-W method for stiff ODEs and DAEs in mass matrix form. 2nd order stiff-aware interpolation and additional error test for interpolation.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Rodas23W(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Steinebach G., Rosenbrock methods within OrdinaryDiffEq.jl - Overview, recent developments and applications - Preprint 2024. Proceedings of the JuliaCon Conferences. https://proceedings.juliacon.org/papers/eb04326e1de8fa819a3595b376508a40
source
OrdinaryDiffEqRosenbrock.Rodas3PType
Rodas3P(; autodiff = AutoForwardDiff(),
          concrete_jac = nothing,
          linsolve = nothing,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner Method. 3rd order A-stable and stiffly stable Rosenbrock method with a stiff-aware 3rd order interpolant and additional error test for interpolation. Keeps accuracy on discretizations of linear parabolic PDEs.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Steinebach G., Rosenbrock methods within OrdinaryDiffEq.jl - Overview, recent developments and applications - Preprint 2024. Proceedings of the JuliaCon Conferences. https://proceedings.juliacon.org/papers/eb04326e1de8fa819a3595b376508a40
source
OrdinaryDiffEqRosenbrock.Rodas4Type
Rodas4(; autodiff = AutoForwardDiff(),
         concrete_jac = nothing,
         linsolve = nothing,
         step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner Method. A 4th order A-stable stiffly stable Rosenbrock method with a stiff-aware 3rd order interpolant

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • E. Hairer, G. Wanner, Solving ordinary differential equations II, stiff and differential-algebraic problems. Computational mathematics (2nd revised ed.), Springer (1996)
source
OrdinaryDiffEqRosenbrock.Rodas42Type
Rodas42(; autodiff = AutoForwardDiff(),
          concrete_jac = nothing,
          linsolve = nothing,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner Method. A 4th order A-stable stiffly stable Rosenbrock method with a stiff-aware 3rd order interpolant

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • E. Hairer, G. Wanner, Solving ordinary differential equations II, stiff and differential-algebraic problems. Computational mathematics (2nd revised ed.), Springer (1996)
source
OrdinaryDiffEqRosenbrock.Rodas4PType
Rodas4P(; autodiff = AutoForwardDiff(),
          concrete_jac = nothing,
          linsolve = nothing,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner Method. 4th order A-stable stiffly stable Rosenbrock method with a stiff-aware 3rd order interpolant. 4th order on linear parabolic problems and 3rd order accurate on nonlinear parabolic problems (as opposed to lower if not corrected).

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Steinebach, G., Rentrop, P., An adaptive method of lines approach for modelling flow and transport in rivers. Adaptive method of lines , Wouver, A. Vande, Sauces, Ph., Schiesser, W.E. (ed.),S. 181-205,Chapman & Hall/CRC, 2001,
  • Steinebach, G., Order-reduction of ROW-methods for DAEs and method of lines applications. Preprint-Nr. 1741, FB Mathematik, TH Darmstadt, 1995.
source
OrdinaryDiffEqRosenbrock.Rodas4P2Type
Rodas4P2(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner-W(olfbrandt) Method. A 4th order L-stable stiffly stable Rosenbrock method with a stiff-aware 3rd order interpolant. 4th order on linear parabolic problems and 3rd order accurate on nonlinear parabolic problems. It is an improvement of Rodas4P and in case of inexact Jacobians a second order W method.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Rodas4P2(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Steinebach G., Improvement of Rosenbrock-Wanner Method RODASP, In: Reis T., Grundel S., Schöps S. (eds) Progress in Differential-Algebraic Equations II. Differential-Algebraic Equations Forum. Springer, Cham., 165-184, 2020.
source
OrdinaryDiffEqRosenbrock.Rodas5Type
Rodas5(; autodiff = AutoForwardDiff(),
         concrete_jac = nothing,
         linsolve = nothing,
         step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner Method. A 5th order A-stable stiffly stable Rosenbrock method with a stiff-aware 4th order interpolant.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Di Marzo G. RODAS5(4) – Méthodes de Rosenbrock d'ordre 5(4) adaptées aux problèmes différentiels-algébriques. MSc mathematics thesis, Faculty of Science, University of Geneva, Switzerland.
source
OrdinaryDiffEqRosenbrock.Rodas5PType
Rodas5P(; autodiff = AutoForwardDiff(),
          concrete_jac = nothing,
          linsolve = nothing,
          step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner-W(olfbrandt) Method. A 5th order A-stable stiffly stable Rosenbrock method with a stiff-aware 4th order interpolant. Has improved stability in the adaptive time stepping embedding.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Rodas5P(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Steinebach G. Construction of Rosenbrock–Wanner method Rodas5P and numerical benchmarks within the Julia Differential Equations package. In: BIT Numerical Mathematics, 63(2), 2023. doi:10.1007/s10543-023-00967-x
source
OrdinaryDiffEqRosenbrock.Rodas5PeType
Rodas5Pe(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner-W(olfbrandt) Method. Variant of Rodas5P with modified embedded scheme.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Rodas5Pe(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Steinebach G. Rosenbrock methods within OrdinaryDiffEq.jl - Overview, recent developments and applications - Preprint 2024. Proceedings of the JuliaCon Conferences. https://proceedings.juliacon.org/papers/eb04326e1de8fa819a3595b376508a40
source
OrdinaryDiffEqRosenbrock.Rodas5PrType
Rodas5Pr(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing,
           step_limiter! = OrdinaryDiffEq.trivial_limiter!)

Rosenbrock-Wanner-W(olfbrandt) Method. Variant of Rodas5P with additional residual control.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Rodas5Pr(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.
  • step_limiter!: function of the form limiter!(u, integrator, p, t)

References

  • Steinebach G. Rosenbrock methods within OrdinaryDiffEq.jl - Overview, recent developments and applications - Preprint 2024. Proceedings of the JuliaCon Conferences. https://proceedings.juliacon.org/papers/eb04326e1de8fa819a3595b376508a40
source
OrdinaryDiffEqRosenbrock.RosenbrockW6S4OSType
RosenbrockW6S4OS(; autodiff = AutoForwardDiff(),
                   concrete_jac = nothing,
                   linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. A 4th order L-stable Rosenbrock-W method (fixed step only).

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify RosenbrockW6S4OS(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

https://doi.org/10.1016/j.cam.2009.09.017

source
OrdinaryDiffEqRosenbrock.ROS2Type
ROS2(; autodiff = AutoForwardDiff(),
       concrete_jac = nothing,
       linsolve = nothing)

Rosenbrock-Wanner Method. A 2nd order L-stable Rosenbrock method with 2 internal stages.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • J. G. Verwer et al. (1999): A second-order Rosenbrock method applied to photochemical dispersion problems https://doi.org/10.1137/S1064827597326651
source
OrdinaryDiffEqRosenbrock.ROS2PRType
ROS2PR(; autodiff = AutoForwardDiff(),
         concrete_jac = nothing,
         linsolve = nothing)

Rosenbrock-Wanner Method. 2nd order stiffly accurate Rosenbrock method with 3 internal stages with (Rinf=0). For problems with medium stiffness the convergence behaviour is very poor and it is recommended to use ROS2S instead.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • Rang, Joachim (2014): The Prothero and Robinson example: Convergence studies for Runge-Kutta and Rosenbrock-Wanner methods. https://doi.org/10.24355/dbbs.084-201408121139-0
source
OrdinaryDiffEqRosenbrock.ROS2SType
ROS2S(; autodiff = AutoForwardDiff(),
        concrete_jac = nothing,
        linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. 2nd order stiffly accurate Rosenbrock-Wanner W-method with 3 internal stages with B_PR consistent of order 2 with (Rinf=0).

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify ROS2S(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • Rang, Joachim (2014): The Prothero and Robinson example: Convergence studies for Runge-Kutta and Rosenbrock-Wanner methods. https://doi.org/10.24355/dbbs.084-201408121139-0
source
OrdinaryDiffEqRosenbrock.ROS3Type
ROS3(; autodiff = AutoForwardDiff(),
       concrete_jac = nothing,
       linsolve = nothing)

Rosenbrock-Wanner Method. 3rd order L-stable Rosenbrock method with 3 internal stages with an embedded strongly A-stable 2nd order method.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • E. Hairer, G. Wanner, Solving ordinary differential equations II, stiff and differential-algebraic problems. Computational mathematics (2nd revised ed.), Springer (1996)
source
OrdinaryDiffEqRosenbrock.ROS3PRType
ROS3PR(; autodiff = AutoForwardDiff(),
         concrete_jac = nothing,
         linsolve = nothing)

Rosenbrock-Wanner Method. 3nd order stiffly accurate Rosenbrock method with 3 internal stages with B_PR consistent of order 3, which is strongly A-stable with Rinf~=-0.73.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • Rang, Joachim (2014): The Prothero and Robinson example: Convergence studies for Runge-Kutta and Rosenbrock-Wanner methods. https://doi.org/10.24355/dbbs.084-201408121139-0
source
OrdinaryDiffEqRosenbrock.Scholz4_7Type
Scholz4_7(; autodiff = AutoForwardDiff(),
            concrete_jac = nothing,
            linsolve = nothing)

Rosenbrock-Wanner Method. 3nd order stiffly accurate Rosenbrock method with 3 internal stages with B_PR consistent of order 3, which is strongly A-stable with Rinf~=-0.73. Convergence with order 4 for the stiff case, but has a poor accuracy.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • Rang, Joachim (2014): The Prothero and Robinson example: Convergence studies for Runge-Kutta and Rosenbrock-Wanner methods. https://doi.org/10.24355/dbbs.084-201408121139-0
source
OrdinaryDiffEqRosenbrock.ROS34PW1aType
ROS34PW1a(; autodiff = AutoForwardDiff(),
            concrete_jac = nothing,
            linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. A 4th order L-stable Rosenbrock-W method.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify ROS34PW1a(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • Rang, Joachim and Angermann, L (2005): New Rosenbrock W-methods of order 3 for partial differential algebraic equations of index 1. BIT Numerical Mathematics, 45, 761–787.
source
OrdinaryDiffEqRosenbrock.ROS34PW1bType
ROS34PW1b(; autodiff = AutoForwardDiff(),
            concrete_jac = nothing,
            linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. A 4th order L-stable Rosenbrock-W method.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify ROS34PW1b(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • Rang, Joachim and Angermann, L (2005): New Rosenbrock W-methods of order 3 for partial differential algebraic equations of index 1. BIT Numerical Mathematics, 45, 761–787.
source
OrdinaryDiffEqRosenbrock.ROS34PW2Type
ROS34PW2(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. A 4th order stiffy accurate Rosenbrock-W method for PDAEs.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify ROS34PW2(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • Rang, Joachim and Angermann, L (2005): New Rosenbrock W-methods of order 3 for partial differential algebraic equations of index 1. BIT Numerical Mathematics, 45, 761–787.
source
OrdinaryDiffEqRosenbrock.ROS34PW3Type
ROS34PW3(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. A 4th order strongly A-stable (Rinf~0.63) Rosenbrock-W method.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify ROS34PW3(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • Rang, Joachim and Angermann, L (2005): New Rosenbrock W-methods of order 3 for partial differential algebraic equations of index 1. BIT Numerical Mathematics, 45, 761–787.
source
OrdinaryDiffEqRosenbrock.ROS34PRwType
ROS34PRw(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. 3rd order stiffly accurate Rosenbrock-Wanner W-method with 4 internal stages, B_PR consistent of order 2. The order of convergence decreases if medium stiff problems are considered.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify ROS34PRw(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • Joachim Rang, Improved traditional Rosenbrock–Wanner methods for stiff ODEs and DAEs, Journal of Computational and Applied Mathematics, https://doi.org/10.1016/j.cam.2015.03.010
source
OrdinaryDiffEqRosenbrock.ROS3PRLType
ROS3PRL(; autodiff = AutoForwardDiff(),
          concrete_jac = nothing,
          linsolve = nothing)

Rosenbrock-Wanner Method. 3rd order stiffly accurate Rosenbrock method with 4 internal stages, B_PR consistent of order 2 with Rinf=0. The order of convergence decreases if medium stiff problems are considered, but it has good results for very stiff cases.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • Rang, Joachim (2014): The Prothero and Robinson example: Convergence studies for Runge-Kutta and Rosenbrock-Wanner methods. https://doi.org/10.24355/dbbs.084-201408121139-0
source
OrdinaryDiffEqRosenbrock.ROS3PRL2Type
ROS3PRL2(; autodiff = AutoForwardDiff(),
           concrete_jac = nothing,
           linsolve = nothing)

Rosenbrock-Wanner Method. 3rd order stiffly accurate Rosenbrock method with 4 internal stages, B_PR consistent of order 3. The order of convergence does NOT decreases if medium stiff problems are considered as it does for ROS3PRL.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • Rang, Joachim (2014): The Prothero and Robinson example: Convergence studies for Runge-Kutta and Rosenbrock-Wanner methods. https://doi.org/10.24355/dbbs.084-201408121139-0
source
OrdinaryDiffEqRosenbrock.ROK4aType
ROK4a(; autodiff = AutoForwardDiff(),
        concrete_jac = nothing,
        linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. 4rd order L-stable Rosenbrock-Krylov method with 4 internal stages, with a 3rd order embedded method which is strongly A-stable with Rinf~=0.55. (when using exact Jacobians)

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify ROK4a(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • Tranquilli, Paul and Sandu, Adrian (2014): Rosenbrock–Krylov Methods for Large Systems of Differential Equations https://doi.org/10.1137/130923336
source
OrdinaryDiffEqRosenbrock.RosShamp4Type
RosShamp4(; autodiff = AutoForwardDiff(),
            concrete_jac = nothing,
            linsolve = nothing)

Rosenbrock-Wanner Method. An A-stable 4th order Rosenbrock method.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • L. F. Shampine, Implementation of Rosenbrock Methods, ACM Transactions on Mathematical Software (TOMS), 8: 2, 93-113. doi:10.1145/355993.355994
source
OrdinaryDiffEqRosenbrock.Veldd4Type
Veldd4(; autodiff = AutoForwardDiff(),
         concrete_jac = nothing,
         linsolve = nothing)

Rosenbrock-Wanner Method. A 4th order D-stable Rosenbrock method.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • van Veldhuizen, D-stability and Kaps-Rentrop-methods, M. Computing (1984) 32: 229. doi:10.1007/BF02243574
source
OrdinaryDiffEqRosenbrock.Velds4Type
Velds4(; autodiff = AutoForwardDiff(),
         concrete_jac = nothing,
         linsolve = nothing)

Rosenbrock-Wanner-W(olfbrandt) Method. A 4th order A-stable Rosenbrock method.

Keyword Arguments

  • autodiff: Uses ADTypes.jl to specify whether to use automatic differentiation via ForwardDiff.jl or finite differencing via FiniteDiff.jl. Defaults to AutoForwardDiff() for automatic differentiation, which by default uses chunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To use FiniteDiff.jl, the AutoFiniteDiff() ADType can be used, which has a keyword argument fdtype with default value Val{:forward}(), and alternatives Val{:central}() and Val{:complex}().
  • concrete_jac: Specifies whether a Jacobian should be constructed. Defaults to nothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used for linsolve.
  • linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specify Velds4(linsolve = KLUFactorization()). When nothing is passed, uses DefaultLinearSolver.

References

  • van Veldhuizen, D-stability and Kaps-Rentrop-methods, M. Computing (1984) 32: 229. doi:10.1007/BF02243574
source
OrdinaryDiffEqRosenbrock.GRK4TType
GRK4T(; autodiff = AutoForwardDiff(),
        concrete_jac = nothing,
        linsolve = nothing)

Rosenbrock-Wanner Method. An efficient 4th order Rosenbrock method.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • Kaps, P. & Rentrop, Generalized Runge-Kutta methods of order four with stepsize control for stiff ordinary differential equations. P. Numer. Math. (1979) 33: 55. doi:10.1007/BF01396495
source
OrdinaryDiffEqRosenbrock.GRK4AType
GRK4A(; autodiff = AutoForwardDiff(),
        concrete_jac = nothing,
        linsolve = nothing)

Rosenbrock-Wanner Method. An A-stable 4th order Rosenbrock method. Essentially "anti-L-stable" but efficient.

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • Kaps, P. & Rentrop, Generalized Runge-Kutta methods of order four with stepsize control for stiff ordinary differential equations. P. Numer. Math. (1979) 33: 55. doi:10.1007/BF01396495
source
OrdinaryDiffEqRosenbrock.Ros4LStabType
Ros4LStab(; autodiff = AutoForwardDiff(),
            concrete_jac = nothing,
            linsolve = nothing)

Rosenbrock-Wanner Method. A 4th order A-stable stiffly stable Rosenbrock method with a stiff-aware 3rd order interpolant

Keyword Arguments

  • autodiff: boolean to control if the Jacobian should be computed via AD or not
  • concrete_jac: function of the form jac!(J, u, p, t)
  • linsolve: custom solver for the inner linear systems

References

  • E. Hairer, G. Wanner, Solving ordinary differential equations II, stiff and differential-algebraic problems. Computational mathematics (2nd revised ed.), Springer (1996)
source