OrdinaryDiffEqRosenbrock
Rosenbrock methods for mass matrix differential-algebraic equations (DAEs) and stiff ODEs with singular mass matrices. These methods provide efficient integration for moderately stiff systems with algebraic constraints, offering excellent performance for small to medium-sized DAE problems.
Key Properties
Mass matrix Rosenbrock methods provide:
- DAE capability for index-1 differential-algebraic equations
- W-method efficiency using approximate Jacobians for computational savings
- Mass matrix support for singular and non-diagonal mass matrices
- Moderate to high order accuracy (2nd to 6th order available)
- Good stability properties with stiffly accurate behavior
- Embedded error estimation for adaptive timestepping
When to Use Mass Matrix Rosenbrock Methods
These methods are recommended for:
- Index-1 DAE systems with moderate stiffness
- Small to medium constrained systems (< 1000 equations)
- Semi-explicit DAEs arising from discretized PDEs
- Problems requiring good accuracy with moderate computational cost
- DAEs with moderate nonlinearity where W-methods are efficient
- Electrical circuits and mechanical systems with constraints
In order to use OrdinaryDiffEqRosenbrock with DAEs that require a non-trivial consistent initialization, a nonlinear solver is required and thus using OrdinaryDiffEqNonlinearSolve is required or you must pass an initializealg with a valid nlsolve choice.
Mathematical Background
Mass matrix DAEs have the form: M du/dt = f(u,t)
Rosenbrock methods linearize around the current solution and solve linear systems of the form: (M/γh - J) k_i = ...
where J is the Jacobian of f and γ is a method parameter.
Solver Selection Guide
Recommended Methods by Tolerance
- High tolerances (>1e-2):
Rosenbrock23- efficient low-order method - Medium tolerances (1e-8 to 1e-2):
Rodas5P- most efficient choice, orRodas4Pfor higher reliability - Low tolerances (<1e-8):
Rodas5Peor higher-order alternatives
Method families
Rodas5P: Recommended - Most efficient 5th-order method for general useRodas4P: More reliable 4th-order alternativeRosenbrock23: Good for high tolerance problemsRodas5: Standard 5th-order method without embedded pair optimization
Performance Guidelines
When mass matrix Rosenbrock methods excel
- Small to medium DAE systems (< 1000 equations)
- Moderately stiff problems where full BDF methods are overkill
- Problems with efficient Jacobian computation or finite difference approximation
- Index-1 DAEs with well-conditioned mass matrices
- Semi-explicit index-1 problems from spatial discretizations
System size considerations
- Small systems (< 100): Rosenbrock methods often outperform multistep methods
- Medium systems (100-1000): Good performance with proper linear algebra
- Large systems (> 1000): Consider BDF methods instead
Important DAE Considerations
Initial conditions
- Must be consistent with algebraic constraints
- Consistent initialization may require nonlinear solver
- Index-1 assumption for reliable performance
Mass matrix requirements
- Index-1 DAE structure for optimal performance
- Non-singular leading submatrix for differential variables
- Well-conditioned constraint equations
Alternative Approaches
Consider these alternatives:
- Mass matrix BDF methods for larger or highly stiff DAE systems
- Implicit Runge-Kutta methods for higher accuracy requirements
- Standard Rosenbrock methods for regular ODEs without constraints
- IMEX methods if natural explicit/implicit splitting exists
Example Usage
using LinearAlgebra: Diagonal
function rober(du, u, p, t)
y₁, y₂, y₃ = u
k₁, k₂, k₃ = p
du[1] = -k₁ * y₁ + k₃ * y₂ * y₃
du[2] = k₁ * y₁ - k₃ * y₂ * y₃ - k₂ * y₂^2
du[3] = y₁ + y₂ + y₃ - 1
nothing
end
M = Diagonal([1.0, 1.0, 0]) # Singular mass matrix
f = ODEFunction(rober, mass_matrix = M)
prob_mm = ODEProblem(f, [1.0, 0.0, 0.0], (0.0, 1e5), (0.04, 3e7, 1e4))
sol = solve(prob_mm, Rodas5(), reltol = 1e-8, abstol = 1e-8)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.Rosenbrock23 — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRosenbrock23(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.step_limiter!: function of the formlimiter!(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.
OrdinaryDiffEqRosenbrock.Rosenbrock32 — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRosenbrock32(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.step_limiter!: function of the formlimiter!(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.
OrdinaryDiffEqRosenbrock.ROS3P — Type
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 notconcrete_jac: function of the formjac!(J, u, p, t)linsolve: custom solver for the inner linear systemsstep_limiter!: function of the formlimiter!(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
OrdinaryDiffEqRosenbrock.Rodas3 — Type
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 notconcrete_jac: function of the formjac!(J, u, p, t)linsolve: custom solver for the inner linear systemsstep_limiter!: function of the formlimiter!(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.
OrdinaryDiffEqRosenbrock.Rodas23W — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRodas23W(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.step_limiter!: function of the formlimiter!(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
OrdinaryDiffEqRosenbrock.Rodas3P — Type
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 notconcrete_jac: function of the formjac!(J, u, p, t)linsolve: custom solver for the inner linear systemsstep_limiter!: function of the formlimiter!(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
OrdinaryDiffEqRosenbrock.Rodas4 — Type
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 notconcrete_jac: function of the formjac!(J, u, p, t)linsolve: custom solver for the inner linear systemsstep_limiter!: function of the formlimiter!(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)
OrdinaryDiffEqRosenbrock.Rodas42 — Type
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 notconcrete_jac: function of the formjac!(J, u, p, t)linsolve: custom solver for the inner linear systemsstep_limiter!: function of the formlimiter!(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)
OrdinaryDiffEqRosenbrock.Rodas4P — Type
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 notconcrete_jac: function of the formjac!(J, u, p, t)linsolve: custom solver for the inner linear systemsstep_limiter!: function of the formlimiter!(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.
OrdinaryDiffEqRosenbrock.Rodas4P2 — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRodas4P2(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.step_limiter!: function of the formlimiter!(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.
OrdinaryDiffEqRosenbrock.Rodas5 — Type
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 notconcrete_jac: function of the formjac!(J, u, p, t)linsolve: custom solver for the inner linear systemsstep_limiter!: function of the formlimiter!(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.
OrdinaryDiffEqRosenbrock.Rodas5P — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRodas5P(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.step_limiter!: function of the formlimiter!(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
OrdinaryDiffEqRosenbrock.Rodas5Pe — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRodas5Pe(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.step_limiter!: function of the formlimiter!(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
OrdinaryDiffEqRosenbrock.Rodas5Pr — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRodas5Pr(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.step_limiter!: function of the formlimiter!(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
OrdinaryDiffEqRosenbrock.RosenbrockW6S4OS — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyRosenbrockW6S4OS(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
References
https://doi.org/10.1016/j.cam.2009.09.017
OrdinaryDiffEqRosenbrock.ROS2 — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.ROS2PR — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.ROS2S — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyROS2S(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
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
OrdinaryDiffEqRosenbrock.ROS3 — Type
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 notconcrete_jac: function of the formjac!(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)
OrdinaryDiffEqRosenbrock.ROS3PR — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.Scholz4_7 — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.ROS34PW1a — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyROS34PW1a(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
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.
OrdinaryDiffEqRosenbrock.ROS34PW1b — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyROS34PW1b(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
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.
OrdinaryDiffEqRosenbrock.ROS34PW2 — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyROS34PW2(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
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.
OrdinaryDiffEqRosenbrock.ROS34PW3 — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyROS34PW3(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
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.
OrdinaryDiffEqRosenbrock.ROS34PRw — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyROS34PRw(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
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
OrdinaryDiffEqRosenbrock.ROS3PRL — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.ROS3PRL2 — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.ROK4a — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyROK4a(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
References
- Tranquilli, Paul and Sandu, Adrian (2014): Rosenbrock–Krylov Methods for Large Systems of Differential Equations https://doi.org/10.1137/130923336
OrdinaryDiffEqRosenbrock.RosShamp4 — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.Veldd4 — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.Velds4 — Type
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 toAutoForwardDiff()for automatic differentiation, which by default useschunksize = 0, and thus uses the internal ForwardDiff.jl algorithm for the choice. To useFiniteDiff.jl, theAutoFiniteDiff()ADType can be used, which has a keyword argumentfdtypewith default valueVal{:forward}(), and alternativesVal{:central}()andVal{:complex}().concrete_jac: Specifies whether a Jacobian should be constructed. Defaults tonothing, which means it will be chosen true/false depending on circumstances of the solver, such as whether a Krylov subspace method is used forlinsolve.linsolve: Any LinearSolve.jl compatible linear solver. For example, to use KLU.jl, specifyVelds4(linsolve = KLUFactorization()). Whennothingis passed, usesDefaultLinearSolver.
References
- van Veldhuizen, D-stability and Kaps-Rentrop-methods, M. Computing (1984) 32: 229. doi:10.1007/BF02243574
OrdinaryDiffEqRosenbrock.GRK4T — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.GRK4A — Type
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 notconcrete_jac: function of the formjac!(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
OrdinaryDiffEqRosenbrock.Ros4LStab — Type
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 notconcrete_jac: function of the formjac!(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)