OrdinaryDiffEqBDF
Backward Differentiation Formula (BDF) methods are multistep implicit methods specifically designed for solving large stiff systems of differential equations. They are the preferred choice for very large systems (>1000 equations) where other implicit methods become computationally expensive.
Key Properties
BDF methods offer:
- Excellent efficiency for large systems (>1000 ODEs)
- L-stable behavior for orders 1 and 2 only
- Adaptive order and stepsize control for optimal performance
- Alpha-stability for higher orders (but less stable than L-stable methods for problems with large complex eigenvalues)
When to Use BDF Methods
BDF methods are recommended for:
- Large stiff systems with more than 1000 equations
- Very stiff problems where other implicit methods struggle
- Long-time integration of stiff systems
- Parabolic PDEs after spatial discretization
- Reaction-diffusion systems and chemical kinetics
- Circuit simulation and other engineering applications with large stiff systems
Solver Selection Guide
Recommended methods
- QNDF: Adaptive order quasi-constant timestep BDF, best general choice for large systems
- FBDF: Fixed-leading coefficient BDF, often more efficient than QNDF
Performance Characteristics
- Most efficient for systems with >1000 equations
- Outperform Runge-Kutta methods on very large stiff systems
- Memory efficient due to multistep structure
- Excel at very low tolerances (1e-9 and below)
- Particularly effective for problems arising from PDE discretizations
Comparison with Other Methods
Choose BDF methods over:
- Rosenbrock methods: When system size > 1000 equations
- SDIRK methods: For very large stiff systems where RK methods become expensive
- Explicit methods: For any stiff problem
Choose other methods over BDF when:
- System size < 100: Rosenbrock or SDIRK methods often more efficient
- Problems with large complex eigenvalues: Rosenbrock and L-stable SDIRK methods are more stable due to BDF methods only being alpha-stable
- Moderate stiffness: SDIRK methods may be more robust
- Non-stiff problems: Use explicit methods like Tsit5
Installation
To be able to access the solvers in OrdinaryDiffEqBDF, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqBDF")This will only install the solvers listed at the bottom of this page. If you want to explore other solvers for your problem, you will need to install some of the other libraries listed in the navigation bar on the left.
Example usage
using OrdinaryDiffEqBDF
function lorenz!(du, u, p, t)
    du[1] = 10.0 * (u[2] - u[1])
    du[2] = u[1] * (28.0 - u[3]) - u[2]
    du[3] = u[1] * u[2] - (8 / 3) * u[3]
end
u0 = [1.0; 0.0; 0.0]
tspan = (0.0, 100.0)
prob = ODEProblem(lorenz!, u0, tspan)
sol = solve(prob, QNDF())Full list of solvers
OrdinaryDiffEqBDF.ABDF2 — TypeABDF2(; chunk_size = Val{0}(),
        autodiff = AutoForwardDiff(),
        standardtag = Val{true}(),
        concrete_jac = nothing,
        linsolve = nothing,
        precs = DEFAULT_PRECS,
        κ = nothing,
        tol = nothing,
        nlsolve = NLNewton(),
        smooth_est = true,
        extrapolant = :linear,
        controller = :Standard,
        step_limiter! = trivial_limiter!)Multistep Method. An adaptive order 2 L-stable fixed leading coefficient multistep BDF 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- fdtypewith default value- Val{:forward}(), and alternatives- Val{:central}()and- Val{:complex}().
- standardtag: Specifies whether to use package-specific tags instead of the ForwardDiff default function-specific tags. For more information, see this blog post. Defaults to- Val{true}().
- 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- ABDF2(linsolve = KLUFactorization()). When- nothingis passed, uses- DefaultLinearSolver.
- precs: Any LinearSolve.jl-compatible preconditioner can be used as a left or right preconditioner. Preconditioners are specified by the- Pl,Pr = precs(W,du,u,p,t,newW,Plprev,Prprev,solverdata)function where the arguments are defined as:- W: the current Jacobian of the nonlinear system. Specified as either $I - \gamma J$ or $I/\gamma - J$ depending on the algorithm. This will commonly be a- WOperatortype defined by OrdinaryDiffEq.jl. It is a lazy representation of the operator. Users can construct the W-matrix on demand by calling- convert(AbstractMatrix,W)to receive an- AbstractMatrixmatching the- jac_prototype.
- du: the current ODE derivative
- u: the current ODE state
- p: the ODE parameters
- t: the current ODE time
- newW: a- Boolwhich specifies whether the- Wmatrix has been updated since the last call to- precs. It is recommended that this is checked to only update the preconditioner when- newW == true.
- Plprev: the previous- Pl.
- Prprev: the previous- Pr.
- solverdata: Optional extra data the solvers can give to the- precsfunction. Solver-dependent and subject to change.
 - (Pl,Pr)of the LinearSolve.jl-compatible preconditioners. To specify one-sided preconditioning, simply return- nothingfor the preconditioner which is not used. Additionally,- precsmust supply the dispatch:
 which is used in the solver setup phase to construct the integrator type with the preconditioners- Pl, Pr = precs(W, du, u, p, t, ::Nothing, ::Nothing, ::Nothing, solverdata)- (Pl,Pr). The default is- precs=DEFAULT_PRECSwhere the default preconditioner function is defined as:
 /n-- DEFAULT_PRECS(W, du, u, p, t, newW, Plprev, Prprev, solverdata) = nothing, nothing- κ: TBD
- tol: TBD
- nlsolve: TBD
- smooth_est: TBD
- extrapolant: TBD
- controller: TBD
- step_limiter!: function of the form- limiter!(u, integrator, p, t)
References
E. Alberdi Celayaa, J. J. Anza Aguirrezabalab, P. Chatzipantelidisc. Implementation of an Adaptive BDF2 Formula and Comparison with The MATLAB Ode15s. Procedia Computer Science, 29, pp 1014-1026, 2014. doi: https://doi.org/10.1016/j.procs.2014.05.091
OrdinaryDiffEqBDF.QNDF — TypeQNDF(; chunk_size = Val{0}(),
       autodiff = AutoForwardDiff(),
       standardtag = Val{true}(),
       concrete_jac = nothing,
       linsolve = nothing,
       precs = DEFAULT_PRECS,
       κ = nothing,
       tol = nothing,
       nlsolve = NLNewton(),
       extrapolant = :linear,
       kappa =  promote(-0.1850, -1 // 9, -0.0823, -0.0415, 0),
       controller = :Standard,
       step_limiter! = trivial_limiter!)Multistep Method. An adaptive order quasi-constant timestep NDF method. Similar to MATLAB's ode15s. Uses Shampine's accuracy-optimal coefficients. Performance improves with larger, more complex ODEs. Good for medium to highly stiff problems. Recommended for large systems (>1000 ODEs).
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- fdtypewith default value- Val{:forward}(), and alternatives- Val{:central}()and- Val{:complex}().
- standardtag: Specifies whether to use package-specific tags instead of the ForwardDiff default function-specific tags. For more information, see this blog post. Defaults to- Val{true}().
- 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- QNDF(linsolve = KLUFactorization()). When- nothingis passed, uses- DefaultLinearSolver.
- precs: Any LinearSolve.jl-compatible preconditioner can be used as a left or right preconditioner. Preconditioners are specified by the- Pl,Pr = precs(W,du,u,p,t,newW,Plprev,Prprev,solverdata)function where the arguments are defined as:- W: the current Jacobian of the nonlinear system. Specified as either $I - \gamma J$ or $I/\gamma - J$ depending on the algorithm. This will commonly be a- WOperatortype defined by OrdinaryDiffEq.jl. It is a lazy representation of the operator. Users can construct the W-matrix on demand by calling- convert(AbstractMatrix,W)to receive an- AbstractMatrixmatching the- jac_prototype.
- du: the current ODE derivative
- u: the current ODE state
- p: the ODE parameters
- t: the current ODE time
- newW: a- Boolwhich specifies whether the- Wmatrix has been updated since the last call to- precs. It is recommended that this is checked to only update the preconditioner when- newW == true.
- Plprev: the previous- Pl.
- Prprev: the previous- Pr.
- solverdata: Optional extra data the solvers can give to the- precsfunction. Solver-dependent and subject to change.
 - (Pl,Pr)of the LinearSolve.jl-compatible preconditioners. To specify one-sided preconditioning, simply return- nothingfor the preconditioner which is not used. Additionally,- precsmust supply the dispatch:
 which is used in the solver setup phase to construct the integrator type with the preconditioners- Pl, Pr = precs(W, du, u, p, t, ::Nothing, ::Nothing, ::Nothing, solverdata)- (Pl,Pr). The default is- precs=DEFAULT_PRECSwhere the default preconditioner function is defined as:
 /n-- DEFAULT_PRECS(W, du, u, p, t, newW, Plprev, Prprev, solverdata) = nothing, nothing- κ: TBD
- tol: TBD
- nlsolve: TBD
- extrapolant: TBD
- kappa: TBD
- controller: TBD
- step_limiter!: function of the form- limiter!(u, integrator, p, t)
References
@article{shampine1997matlab, title={The matlab ode suite}, author={Shampine, Lawrence F and Reichelt, Mark W}, journal={SIAM journal on scientific computing}, volume={18}, number={1}, pages={1–22}, year={1997}, publisher={SIAM} }
OrdinaryDiffEqBDF.QNDF1 — TypeQNDF1(; chunk_size = Val{0}(),
        autodiff = AutoForwardDiff(),
        standardtag = Val{true}(),
        concrete_jac = nothing,
        linsolve = nothing,
        precs = DEFAULT_PRECS,
        nlsolve = NLNewton(),
        extrapolant = :linear,
        kappa = -0.1850,
        controller = :Standard,
        step_limiter! = trivial_limiter!)Multistep Method. An adaptive order 1 quasi-constant timestep L-stable numerical differentiation function 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- fdtypewith default value- Val{:forward}(), and alternatives- Val{:central}()and- Val{:complex}().
- standardtag: Specifies whether to use package-specific tags instead of the ForwardDiff default function-specific tags. For more information, see this blog post. Defaults to- Val{true}().
- 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- QNDF1(linsolve = KLUFactorization()). When- nothingis passed, uses- DefaultLinearSolver.
- precs: Any LinearSolve.jl-compatible preconditioner can be used as a left or right preconditioner. Preconditioners are specified by the- Pl,Pr = precs(W,du,u,p,t,newW,Plprev,Prprev,solverdata)function where the arguments are defined as:- W: the current Jacobian of the nonlinear system. Specified as either $I - \gamma J$ or $I/\gamma - J$ depending on the algorithm. This will commonly be a- WOperatortype defined by OrdinaryDiffEq.jl. It is a lazy representation of the operator. Users can construct the W-matrix on demand by calling- convert(AbstractMatrix,W)to receive an- AbstractMatrixmatching the- jac_prototype.
- du: the current ODE derivative
- u: the current ODE state
- p: the ODE parameters
- t: the current ODE time
- newW: a- Boolwhich specifies whether the- Wmatrix has been updated since the last call to- precs. It is recommended that this is checked to only update the preconditioner when- newW == true.
- Plprev: the previous- Pl.
- Prprev: the previous- Pr.
- solverdata: Optional extra data the solvers can give to the- precsfunction. Solver-dependent and subject to change.
 - (Pl,Pr)of the LinearSolve.jl-compatible preconditioners. To specify one-sided preconditioning, simply return- nothingfor the preconditioner which is not used. Additionally,- precsmust supply the dispatch:
 which is used in the solver setup phase to construct the integrator type with the preconditioners- Pl, Pr = precs(W, du, u, p, t, ::Nothing, ::Nothing, ::Nothing, solverdata)- (Pl,Pr). The default is- precs=DEFAULT_PRECSwhere the default preconditioner function is defined as:
 /n-- DEFAULT_PRECS(W, du, u, p, t, newW, Plprev, Prprev, solverdata) = nothing, nothing- nlsolve: TBD
- extrapolant: TBD
- kappa: TBD
- controller: TBD
- step_limiter!: function of the form- limiter!(u, integrator, p, t)
References
@article{shampine1997matlab, title={The matlab ode suite}, author={Shampine, Lawrence F and Reichelt, Mark W}, journal={SIAM journal on scientific computing}, volume={18}, number={1}, pages={1–22}, year={1997}, publisher={SIAM} }
OrdinaryDiffEqBDF.QNDF2 — TypeQNDF2(; chunk_size = Val{0}(),
        autodiff = AutoForwardDiff(),
        standardtag = Val{true}(),
        concrete_jac = nothing,
        linsolve = nothing,
        precs = DEFAULT_PRECS,
        nlsolve = NLNewton(),
        extrapolant = :linear,
        kappa =  -1 // 9,
        controller = :Standard,
        step_limiter! = trivial_limiter!)Multistep Method. An adaptive order 2 quasi-constant timestep L-stable numerical differentiation function (NDF) 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- fdtypewith default value- Val{:forward}(), and alternatives- Val{:central}()and- Val{:complex}().
- standardtag: Specifies whether to use package-specific tags instead of the ForwardDiff default function-specific tags. For more information, see this blog post. Defaults to- Val{true}().
- 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- QNDF2(linsolve = KLUFactorization()). When- nothingis passed, uses- DefaultLinearSolver.
- precs: Any LinearSolve.jl-compatible preconditioner can be used as a left or right preconditioner. Preconditioners are specified by the- Pl,Pr = precs(W,du,u,p,t,newW,Plprev,Prprev,solverdata)function where the arguments are defined as:- W: the current Jacobian of the nonlinear system. Specified as either $I - \gamma J$ or $I/\gamma - J$ depending on the algorithm. This will commonly be a- WOperatortype defined by OrdinaryDiffEq.jl. It is a lazy representation of the operator. Users can construct the W-matrix on demand by calling- convert(AbstractMatrix,W)to receive an- AbstractMatrixmatching the- jac_prototype.
- du: the current ODE derivative
- u: the current ODE state
- p: the ODE parameters
- t: the current ODE time
- newW: a- Boolwhich specifies whether the- Wmatrix has been updated since the last call to- precs. It is recommended that this is checked to only update the preconditioner when- newW == true.
- Plprev: the previous- Pl.
- Prprev: the previous- Pr.
- solverdata: Optional extra data the solvers can give to the- precsfunction. Solver-dependent and subject to change.
 - (Pl,Pr)of the LinearSolve.jl-compatible preconditioners. To specify one-sided preconditioning, simply return- nothingfor the preconditioner which is not used. Additionally,- precsmust supply the dispatch:
 which is used in the solver setup phase to construct the integrator type with the preconditioners- Pl, Pr = precs(W, du, u, p, t, ::Nothing, ::Nothing, ::Nothing, solverdata)- (Pl,Pr). The default is- precs=DEFAULT_PRECSwhere the default preconditioner function is defined as:
 /n-- DEFAULT_PRECS(W, du, u, p, t, newW, Plprev, Prprev, solverdata) = nothing, nothing- nlsolve: TBD
- extrapolant: TBD
- kappa: TBD
- controller: TBD
- step_limiter!: function of the form- limiter!(u, integrator, p, t)
References
@article{shampine1997matlab, title={The matlab ode suite}, author={Shampine, Lawrence F and Reichelt, Mark W}, journal={SIAM journal on scientific computing}, volume={18}, number={1}, pages={1–22}, year={1997}, publisher={SIAM} }
OrdinaryDiffEqBDF.QBDF — FunctionQBDF: Multistep Method
An alias of QNDF with κ=0.
OrdinaryDiffEqBDF.QBDF1 — FunctionQBDF1: Multistep Method
An alias of QNDF1 with κ=0.
OrdinaryDiffEqBDF.QBDF2 — FunctionQBDF2: Multistep Method
An alias of QNDF2 with κ=0.
OrdinaryDiffEqBDF.MEBDF2 — TypeMEBDF2(; chunk_size = Val{0}(),
         autodiff = AutoForwardDiff(),
         standardtag = Val{true}(),
         concrete_jac = nothing,
         linsolve = nothing,
         precs = DEFAULT_PRECS,
         nlsolve = NLNewton(),
         extrapolant = :constant)Multistep Method. The second order Modified Extended BDF method, which has improved stability properties over the standard BDF. Fixed timestep 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- fdtypewith default value- Val{:forward}(), and alternatives- Val{:central}()and- Val{:complex}().
- standardtag: Specifies whether to use package-specific tags instead of the ForwardDiff default function-specific tags. For more information, see this blog post. Defaults to- Val{true}().
- 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- MEBDF2(linsolve = KLUFactorization()). When- nothingis passed, uses- DefaultLinearSolver.
- precs: Any LinearSolve.jl-compatible preconditioner can be used as a left or right preconditioner. Preconditioners are specified by the- Pl,Pr = precs(W,du,u,p,t,newW,Plprev,Prprev,solverdata)function where the arguments are defined as:- W: the current Jacobian of the nonlinear system. Specified as either $I - \gamma J$ or $I/\gamma - J$ depending on the algorithm. This will commonly be a- WOperatortype defined by OrdinaryDiffEq.jl. It is a lazy representation of the operator. Users can construct the W-matrix on demand by calling- convert(AbstractMatrix,W)to receive an- AbstractMatrixmatching the- jac_prototype.
- du: the current ODE derivative
- u: the current ODE state
- p: the ODE parameters
- t: the current ODE time
- newW: a- Boolwhich specifies whether the- Wmatrix has been updated since the last call to- precs. It is recommended that this is checked to only update the preconditioner when- newW == true.
- Plprev: the previous- Pl.
- Prprev: the previous- Pr.
- solverdata: Optional extra data the solvers can give to the- precsfunction. Solver-dependent and subject to change.
 - (Pl,Pr)of the LinearSolve.jl-compatible preconditioners. To specify one-sided preconditioning, simply return- nothingfor the preconditioner which is not used. Additionally,- precsmust supply the dispatch:
 which is used in the solver setup phase to construct the integrator type with the preconditioners- Pl, Pr = precs(W, du, u, p, t, ::Nothing, ::Nothing, ::Nothing, solverdata)- (Pl,Pr). The default is- precs=DEFAULT_PRECSwhere the default preconditioner function is defined as:
 /n-- DEFAULT_PRECS(W, du, u, p, t, newW, Plprev, Prprev, solverdata) = nothing, nothing- nlsolve: TBD
- extrapolant: TBD
References
@article{cash2000modified, title={Modified extended backward differentiation formulae for the numerical solution of stiff initial value problems in ODEs and DAEs}, author={Cash, JR}, journal={Journal of Computational and Applied Mathematics}, volume={125}, number={1-2}, pages={117–130}, year={2000}, publisher={Elsevier}}
OrdinaryDiffEqBDF.FBDF — TypeFBDF(; chunk_size = Val{0}(),
       autodiff = AutoForwardDiff(),
       standardtag = Val{true}(),
       concrete_jac = nothing,
       linsolve = nothing,
       precs = DEFAULT_PRECS,
       κ = nothing,
       tol = nothing,
       nlsolve = NLNewton(),
       extrapolant = :linear,
       controller = :Standard,
       step_limiter! = trivial_limiter!,
       max_order::Val{MO} = Val{5}())Multistep Method. An adaptive order quasi-constant timestep NDF method. Fixed leading coefficient BDF. Utilizes Shampine's accuracy-optimal kappa values as defaults (has a keyword argument for a tuple of kappa coefficients).
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- fdtypewith default value- Val{:forward}(), and alternatives- Val{:central}()and- Val{:complex}().
- standardtag: Specifies whether to use package-specific tags instead of the ForwardDiff default function-specific tags. For more information, see this blog post. Defaults to- Val{true}().
- 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- FBDF(linsolve = KLUFactorization()). When- nothingis passed, uses- DefaultLinearSolver.
- precs: Any LinearSolve.jl-compatible preconditioner can be used as a left or right preconditioner. Preconditioners are specified by the- Pl,Pr = precs(W,du,u,p,t,newW,Plprev,Prprev,solverdata)function where the arguments are defined as:- W: the current Jacobian of the nonlinear system. Specified as either $I - \gamma J$ or $I/\gamma - J$ depending on the algorithm. This will commonly be a- WOperatortype defined by OrdinaryDiffEq.jl. It is a lazy representation of the operator. Users can construct the W-matrix on demand by calling- convert(AbstractMatrix,W)to receive an- AbstractMatrixmatching the- jac_prototype.
- du: the current ODE derivative
- u: the current ODE state
- p: the ODE parameters
- t: the current ODE time
- newW: a- Boolwhich specifies whether the- Wmatrix has been updated since the last call to- precs. It is recommended that this is checked to only update the preconditioner when- newW == true.
- Plprev: the previous- Pl.
- Prprev: the previous- Pr.
- solverdata: Optional extra data the solvers can give to the- precsfunction. Solver-dependent and subject to change.
 - (Pl,Pr)of the LinearSolve.jl-compatible preconditioners. To specify one-sided preconditioning, simply return- nothingfor the preconditioner which is not used. Additionally,- precsmust supply the dispatch:
 which is used in the solver setup phase to construct the integrator type with the preconditioners- Pl, Pr = precs(W, du, u, p, t, ::Nothing, ::Nothing, ::Nothing, solverdata)- (Pl,Pr). The default is- precs=DEFAULT_PRECSwhere the default preconditioner function is defined as:
 /n-- DEFAULT_PRECS(W, du, u, p, t, newW, Plprev, Prprev, solverdata) = nothing, nothing- κ: TBD
- tol: TBD
- nlsolve: TBD
- extrapolant: TBD
- controller: TBD
- step_limiter!: function of the form- limiter!(u, integrator, p, t)
- max_order: TBD
References
@article{shampine2002solving, title={Solving 0= F (t, y (t), y′(t)) in Matlab}, author={Shampine, Lawrence F}, year={2002}, publisher={Walter de Gruyter GmbH \& Co. KG}}