OrdinaryDiffEqPDIRK
PDIRK methods are parallel DIRK methods. SDIRK methods, or singly-diagonally implicit methods, have to build and solve a factorize a Jacobian of the form W = I-gammaJ
where gamma
is dependent on the chosen method. PDIRK methods use multiple different choices of gamma
, i.e. W_i = I-gamma_iJ
, which are all used in the update process. There are some advantages to this, as no SDIRK method can be a higher order than 5, while DIRK methods generally can have arbitrarily high order and lower error coefficients, leading to lower errors at larger dt sizes. With the right construction of the tableau, these matrices can be factorized and the underlying steps can be computed in parallel, which is why these are the parallel DIRK methods.
OrdinaryDiffEqPDIRK
is experimental, as there are no parallel DIRK tableaus that achieve good performance in the literature.
Installation
To be able to access the solvers in OrdinaryDiffEqPDIRK
, you must first install them use the Julia package manager:
using Pkg
Pkg.add("OrdinaryDiffEqPDIRK")
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 OrdinaryDiffEqPDIRK
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, PDIRK44())
Full list of solvers
OrdinaryDiffEqPDIRK.PDIRK44
— TypePDIRK44(; chunk_size = Val{0}(),
autodiff = AutoForwardDiff(),
standardtag = Val{true}(),
concrete_jac = nothing,
diff_type = Val{:forward},
linsolve = nothing,
precs = DEFAULT_PRECS,
nlsolve = NLNewton(),
extrapolant = :constant,
thread = OrdinaryDiffEq.True())
Parallel Diagonally Implicit Runge-Kutta Method. A 2 processor 4th order diagonally non-adaptive implicit 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 argumentfdtype
with default valueVal{:forward}()
, and alternativesVal{:central}()
andVal{: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 toVal{true}()
.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, specifyPDIRK44(linsolve = KLUFactorization()
). Whennothing
is passed, usesDefaultLinearSolver
.precs
: Any LinearSolve.jl-compatible preconditioner can be used as a left or right preconditioner. Preconditioners are specified by thePl,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 aWOperator
type defined by OrdinaryDiffEq.jl. It is a lazy representation of the operator. Users can construct the W-matrix on demand by callingconvert(AbstractMatrix,W)
to receive anAbstractMatrix
matching thejac_prototype
.du
: the current ODE derivativeu
: the current ODE statep
: the ODE parameterst
: the current ODE timenewW
: aBool
which specifies whether theW
matrix has been updated since the last call toprecs
. It is recommended that this is checked to only update the preconditioner whennewW == true
.Plprev
: the previousPl
.Prprev
: the previousPr
.solverdata
: Optional extra data the solvers can give to theprecs
function. Solver-dependent and subject to change.
(Pl,Pr)
of the LinearSolve.jl-compatible preconditioners. To specify one-sided preconditioning, simply returnnothing
for the preconditioner which is not used. Additionally,precs
must supply the dispatch:
which is used in the solver setup phase to construct the integrator type with the preconditionersPl, Pr = precs(W, du, u, p, t, ::Nothing, ::Nothing, ::Nothing, solverdata)
(Pl,Pr)
. The default isprecs=DEFAULT_PRECS
where the default preconditioner function is defined as:DEFAULT_PRECS(W, du, u, p, t, newW, Plprev, Prprev, solverdata) = nothing, nothing
nlsolve
: TBD,extrapolant
: TBD,thread
: determines whether internal broadcasting on appropriate CPU arrays should be serial (thread = OrdinaryDiffEq.False()
) or use multiple threads (thread = OrdinaryDiffEq.True()
) when Julia is started with multiple threads.
References
"@article{iserles1990theory, title={On the theory of parallel Runge—Kutta methods}, author={Iserles, Arieh and Norrsett, SP}, journal={IMA Journal of numerical Analysis}, volume={10}, number={4}, pages={463–488}, year={1990}, publisher={Oxford University Press}}