ODESolvers

Low Storage Runge Kutta methods

ClimateMachine.ODESolvers.LowStorageRungeKutta2NType
LowStorageRungeKutta2N(f, RKA, RKB, RKC, Q; dt, t0 = 0)

This is a time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

The constructor builds a low-storage Runge-Kutta scheme using 2N storage based on the provided RKA, RKB and RKC coefficient arrays.

The available concrete implementations are:

source
ClimateMachine.ODESolvers.LSRK54CarpenterKennedyFunction
LSRK54CarpenterKennedy(f, Q; dt, t0 = 0)

This function returns a LowStorageRungeKutta2N time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This uses the fourth-order, low-storage, Runge–Kutta scheme of Carpenter and Kennedy (1994) (in their notation (5,4) 2N-Storage RK scheme).

References

@TECHREPORT{CarpenterKennedy1994,
  author = {M.~H. Carpenter and C.~A. Kennedy},
  title = {Fourth-order {2N-storage} {Runge-Kutta} schemes},
  institution = {National Aeronautics and Space Administration},
  year = {1994},
  number = {NASA TM-109112},
  address = {Langley Research Center, Hampton, VA},
}
source
ClimateMachine.ODESolvers.LSRK144NiegemannDiehlBuschFunction
LSRK144NiegemannDiehlBusch((f, Q; dt, t0 = 0)

This function returns a LowStorageRungeKutta2N time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This uses the fourth-order, 14-stage, low-storage, Runge–Kutta scheme of Niegemann, Diehl, and Busch (2012) with optimized stability region

References

@article{niegemann2012efficient,
  title={Efficient low-storage Runge--Kutta schemes with optimized stability regions},
  author={Niegemann, Jens and Diehl, Richard and Busch, Kurt},
  journal={Journal of Computational Physics},
  volume={231},
  number={2},
  pages={364--372},
  year={2012},
  publisher={Elsevier}
}
source

Low Storage (3N) Runge Kutta methods

ClimateMachine.ODESolvers.LowStorageRungeKutta3NType
LowStorageRungeKutta3N(f, RKA, RKB, RKC, RKW, Q; dt, t0 = 0)

This is a time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

The constructor builds a low-storage Runge–Kutta scheme using 3N storage based on the provided RKA, RKB and RKC coefficient arrays. RKC (vector of length the number of stages ns) set nodal points position; RKA and RKB (size: ns x 2) set weight for tendency and stage-state; RKW (unused) provides RK weight (last row in Butcher's tableau).

The 3-N storage formulation from Fyfe (1966) is applicable to any 4-stage, fourth-order RK scheme. It is implemented here as:

\[\hspace{-20mm} for ~~ j ~~ in ~ [1:ns]: \hspace{10mm} t_j = t^n + \Delta t ~ rkC_j\]
\[ dQ_j = dQ^*_j + f(Q_j,t_j)\]
\[ Q_{j+1} = Q_{j} + \Delta t \{ rkB_{j,1} ~ dQ_j + rkB_{j,2} ~ dR_j \}\]
\[ dR_{j+1} = dR_j + rkA_{j+1,2} ~ dQ_j\]
\[ dQ^*_{j+1} = rkA_{j+1,1} ~ dQ_j\]

The available concrete implementations are:

References

@article{Fyfe1966,
   title = {Economical Evaluation of Runge-Kutta Formulae},
   author = {Fyfe, David J.},
   journal = {Mathematics of Computation},
   volume = {20},
   pages = {392--398},
   year = {1966}
}
source
ClimateMachine.ODESolvers.LS3NRK44ClassicFunction
LS3NRK44Classic(f, Q; dt, t0 = 0)

This function returns a LowStorageRungeKutta3N time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This uses the classic 4-stage, fourth-order Runge–Kutta scheme in the low-storage implementation of Blum (1962)

References

@article {Blum1962,
   title = {A Modification of the Runge-Kutta Fourth-Order Method}
   author = {Blum, E. K.},
   journal = {Mathematics of Computation},
   volume = {16},
   pages = {176-187},
   year = {1962}
}
source
ClimateMachine.ODESolvers.LS3NRK33HeunsFunction
LS3NRK33Heuns(f, Q; dt, t0 = 0)

This function returns a LowStorageRungeKutta3N time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This method uses the 3-stage, third-order Heun's Runge–Kutta scheme.

References

@article {Heun1900,
   title = {Neue Methoden zur approximativen Integration der
   Differentialgleichungen einer unabh"{a}ngigen Ver"{a}nderlichen}
   author = {Heun, Karl},
   journal = {Z. Math. Phys},
   volume = {45},
   pages = {23--38},
   year = {1900}
}
source

Strong Stability Preserving RungeKutta methods

ClimateMachine.ODESolvers.StrongStabilityPreservingRungeKuttaType
StrongStabilityPreservingRungeKutta(f, RKA, RKB, RKC, Q; dt, t0 = 0)

This is a time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

The constructor builds a strong-stability-preserving Runge–Kutta scheme based on the provided RKA, RKB and RKC coefficient arrays.

The available concrete implementations are:

source
ClimateMachine.ODESolvers.SSPRK33ShuOsherFunction
SSPRK33ShuOsher(f, Q; dt, t0 = 0)

This function returns a StrongStabilityPreservingRungeKutta time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This uses the third-order, 3-stage, strong-stability-preserving, Runge–Kutta scheme of Shu and Osher (1988)

References

@article{shu1988efficient,
  title={Efficient implementation of essentially non-oscillatory shock-capturing schemes},
  author={Shu, Chi-Wang and Osher, Stanley},
  journal={Journal of computational physics},
  volume={77},
  number={2},
  pages={439--471},
  year={1988},
  publisher={Elsevier}
}
source
ClimateMachine.ODESolvers.SSPRK34SpiteriRuuthFunction
SSPRK34SpiteriRuuth(f, Q; dt, t0 = 0)

This function returns a StrongStabilityPreservingRungeKutta time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This uses the third-order, 4-stage, strong-stability-preserving, Runge–Kutta scheme of Spiteri and Ruuth (1988)

References

@article{spiteri2002new,
  title={A new class of optimal high-order strong-stability-preserving time discretization methods},
  author={Spiteri, Raymond J and Ruuth, Steven J},
  journal={SIAM Journal on Numerical Analysis},
  volume={40},
  number={2},
  pages={469--491},
  year={2002},
  publisher={SIAM}
}
source

Additive Runge Kutta methods

ClimateMachine.ODESolvers.AdditiveRungeKuttaType
AdditiveRungeKutta(f, l, backward_euler_solver, RKAe, RKAi, RKB, RKC, Q;
                   split_explicit_implicit, variant, dt, t0 = 0)

This is a time stepping object for implicit-explicit time stepping of a decomposed differential equation. When split_explicit_implicit == false the equation is assumed to be decomposed as

\[ \dot{Q} = [l(Q, t)] + [f(Q, t) - l(Q, t)]\]

where Q is the state, f is the full tendency and l is the chosen implicit operator. When split_explicit_implicit == true the assumed decomposition is

\[ \dot{Q} = [l(Q, t)] + [f(Q, t)]\]

where f is now only the nonlinear tendency. For both decompositions the implicit operator l is integrated implicitly whereas the remaining part is integrated explicitly. Other arguments are the required time step size dt and the optional initial time t0. The resulting backward Euler type systems are solved using the provided backward_euler_solver. This time stepping object is intended to be passed to the solve! command.

The constructor builds an additive Runge–Kutta scheme based on the provided RKAe, RKAi, RKB and RKC coefficient arrays. Additionally variant specifies which of the analytically equivalent but numerically different formulations of the scheme is used.

The available concrete implementations are:

source
ClimateMachine.ODESolvers.ARK1ForwardBackwardEulerFunction
ARK1ForwardBackwardEuler(f, l, backward_euler_solver, Q; dt, t0,
                         split_explicit_implicit, variant)

This function returns an AdditiveRungeKutta time stepping object, see the documentation of AdditiveRungeKutta for arguments definitions. This time stepping object is intended to be passed to the solve! command.

This uses a first-order-accurate two-stage additive Runge–Kutta scheme by combining a forward Euler explicit step with a backward Euler implicit correction.

References

@article{Ascher1997,
  title = {Implicit-explicit Runge-Kutta methods for time-dependent
           partial differential equations},
  author = {Uri M. Ascher and Steven J. Ruuth and Raymond J. Spiteri},
  volume = {25},
  number = {2-3},
  pages = {151--167},
  year = {1997},
  journal = {Applied Numerical Mathematics},
  publisher = {Elsevier {BV}}
}
source
ClimateMachine.ODESolvers.ARK2ImplicitExplicitMidpointFunction
ARK2ImplicitExplicitMidpoint(f, l, backward_euler_solver, Q; dt, t0,
                             split_explicit_implicit, variant)

This function returns an AdditiveRungeKutta time stepping object, see the documentation of AdditiveRungeKutta for arguments definitions. This time stepping object is intended to be passed to the solve! command.

This uses a second-order-accurate two-stage additive Runge–Kutta scheme by combining the implicit and explicit midpoint methods.

References

@article{Ascher1997,
  title = {Implicit-explicit Runge-Kutta methods for time-dependent
           partial differential equations},
  author = {Uri M. Ascher and Steven J. Ruuth and Raymond J. Spiteri},
  volume = {25},
  number = {2-3},
  pages = {151--167},
  year = {1997},
  journal = {Applied Numerical Mathematics},
  publisher = {Elsevier {BV}}
}
source
ClimateMachine.ODESolvers.ARK2GiraldoKellyConstantinescuFunction
ARK2GiraldoKellyConstantinescu(f, l, backward_euler_solver, Q; dt, t0,
                               split_explicit_implicit, variant, paperversion)

This function returns an AdditiveRungeKutta time stepping object, see the documentation of AdditiveRungeKutta for arguments definitions. This time stepping object is intended to be passed to the solve! command.

paperversion=true uses the coefficients from the paper, paperversion=false uses coefficients that make the scheme (much) more stable but less accurate

This uses the second-order-accurate 3-stage additive Runge–Kutta scheme of Giraldo, Kelly and Constantinescu (2013).

References

@article{giraldo2013implicit,
  title={Implicit-explicit formulations of a three-dimensional
         nonhydrostatic unified model of the atmosphere ({NUMA})},
  author={Giraldo, Francis X and Kelly, James F and Constantinescu, Emil M},
  journal={SIAM Journal on Scientific Computing},
  volume={35},
  number={5},
  pages={B1162--B1194},
  year={2013},
  publisher={SIAM}
}
source
ClimateMachine.ODESolvers.ARK548L2SA2KennedyCarpenterFunction
ARK548L2SA2KennedyCarpenter(f, l, backward_euler_solver, Q; dt, t0,
                            split_explicit_implicit, variant)

This function returns an AdditiveRungeKutta time stepping object, see the documentation of AdditiveRungeKutta for arguments definitions. This time stepping object is intended to be passed to the solve! command.

This uses the fifth-order-accurate 8-stage additive Runge–Kutta scheme of Kennedy and Carpenter (2013).

References

@article{kennedy2019higher,
  title={Higher-order additive Runge--Kutta schemes for ordinary
         differential equations},
  author={Kennedy, Christopher A and Carpenter, Mark H},
  journal={Applied Numerical Mathematics},
  volume={136},
  pages={183--205},
  year={2019},
  publisher={Elsevier}
}
source
ClimateMachine.ODESolvers.ARK437L2SA1KennedyCarpenterFunction
ARK437L2SA1KennedyCarpenter(f, l, backward_euler_solver, Q; dt, t0,
                            split_explicit_implicit, variant)

This function returns an AdditiveRungeKutta time stepping object, see the documentation of AdditiveRungeKutta for arguments definitions. This time stepping object is intended to be passed to the solve! command.

This uses the fourth-order-accurate 7-stage additive Runge–Kutta scheme of Kennedy and Carpenter (2013).

References

@article{kennedy2019higher,
  title={Higher-order additive Runge--Kutta schemes for ordinary
         differential equations},
  author={Kennedy, Christopher A and Carpenter, Mark H},
  journal={Applied Numerical Mathematics},
  volume={136},
  pages={183--205},
  year={2019},
  publisher={Elsevier}
}
source
ClimateMachine.ODESolvers.SSPRK22RalstonsFunction
SSPRK22Ralstons(f, Q; dt, t0 = 0)

This function returns a StrongStabilityPreservingRungeKutta time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This uses the second-order, 2-stage, strong-stability-preserving, Runge–Kutta scheme of Shu and Osher (1988) (also known as Ralstons's method.) Exact choice of coefficients from wikipedia page for Heun's method :)

References

@article{shu1988efficient,
  title={Efficient implementation of essentially non-oscillatory shock-capturing schemes},
  author={Shu, Chi-Wang and Osher, Stanley},
  journal={Journal of computational physics},
  volume={77},
  number={2},
  pages={439--471},
  year={1988},
  publisher={Elsevier}
}
@article{ralston1962runge,
  title={Runge-Kutta methods with minimum error bounds},
  author={Ralston, Anthony},
  journal={Mathematics of computation},
  volume={16},
  number={80},
  pages={431--437},
  year={1962},
  doi={10.1090/S0025-5718-1962-0150954-0}
}
source
ClimateMachine.ODESolvers.SSPRK22HeunsFunction
SSPRK22Heuns(f, Q; dt, t0 = 0)

This function returns a StrongStabilityPreservingRungeKutta time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This uses the second-order, 2-stage, strong-stability-preserving, Runge–Kutta scheme of Shu and Osher (1988) (also known as Heun's method.) Exact choice of coefficients from wikipedia page for Heun's method :)

References

@article{shu1988efficient,
  title={Efficient implementation of essentially non-oscillatory shock-capturing schemes},
  author={Shu, Chi-Wang and Osher, Stanley},
  journal={Journal of computational physics},
  volume={77},
  number={2},
  pages={439--471},
  year={1988},
  publisher={Elsevier}
}
@article {Heun1900,
   title = {Neue Methoden zur approximativen Integration der
   Differentialgleichungen einer unabh"{a}ngigen Ver"{a}nderlichen}
   author = {Heun, Karl},
   journal = {Z. Math. Phys},
   volume = {45},
   pages = {23--38},
   year = {1900}
}
source
ClimateMachine.ODESolvers.LSRKEulerMethodFunction
LSRKEulerMethod(f, Q; dt, t0 = 0)

This function returns a LowStorageRungeKutta2N time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This method uses the LSRK2N framework to implement a simple Eulerian forward time stepping scheme for the use of debugging.

References

source

Multi-rate Runge Kutta Methods

ClimateMachine.ODESolvers.MultirateRungeKuttaType
MultirateRungeKutta(slow_solver, fast_solver; dt, t0 = 0)

This is a time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f_fast(Q, t) + f_slow(Q, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

The constructor builds a multirate Runge-Kutta scheme using two different RK solvers. This is based on

Currently only the low storage RK methods can be used as slow solvers

References

@article{SchlegelKnothArnoldWolke2012,
  title={Implementation of multirate time integration methods for air
         pollution modelling},
  author={Schlegel, M and Knoth, O and Arnold, M and Wolke, R},
  journal={Geoscientific Model Development},
  volume={5},
  number={6},
  pages={1395--1405},
  year={2012},
  publisher={Copernicus GmbH}
}
source

Split-explicit methods

ClimateMachine.ODESolvers.SplitExplicitSolverType
SplitExplicitSolver(slow_solver, fast_solver; dt, t0 = 0, coupled = true)

This is a time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q_fast} = f_fast(Q_fast, Q_slow, t) \dot{Q_slow} = f_slow(Q_slow, Q_fast, t)\]

with the required time step size dt and optional initial time t0. This time stepping object is intended to be passed to the solve! command.

This method performs an operator splitting to timestep the vertical average of the model at a faster rate than the full model. This results in a first- order time stepper.

source

GARK methods

ClimateMachine.ODESolvers.MRIGARKESDIRK24LSAFunction
MRIGARKESDIRK24LSA(f!,
                   fastsolver,
                   Q;
                   dt,
                   t0 = 0,
                   γ = 0.2,
                   c3 = (2γ + 1) / 2,
                   a32 = 0.2,
                   α = -0.1,
                   β1 = c3 / 10,
                   β2 = c3 / 10,
                   )

A 2nd order, 4 stage decoupled implicit scheme. It is based on an L-Stable, stiffly-accurate ESDIRK.

source
ClimateMachine.ODESolvers.MRIGARKExplicitType
MRIGARKExplicit(f!, fastsolver, Γs, γ̂s, Q, Δt, t0)

Construct an explicit MultiRate Infinitesimal General-structure Additive Runge–Kutta (MRI-GARK) scheme to solve

\[ \dot{y} = f(y, t) + g(y, t)\]

where f is the slow tendency function and g is the fast tendency function; see Sandu (2019).

The fast tendency is integrated using the fastsolver and the slow tendency using the MRI-GARK scheme. Namely, at each stage the scheme solves

\[ v(T_i) &= Y_i \\ \dot{v} &= f(v, t) + \sum_{j=1}^{i} \bar{γ}_{ij}(t) R_j \\ \bar{γ}_{ijk}(t) &= \sum_{k=0}^{NΓ-1} γ_{ijk} τ(t)^k / Δc_s \\ τ(t) &= (t - t_s) / Δt \\ Y_{i+1} &= v(T_i + c_s * Δt)\]

where $Y_1 = y_n$ and $y_{n+1} = Y_{Nstages+1}$.

Here $R_j = g(Y_j, t_0 + c_j * Δt)$ is the tendency for stage $j$, $γ_{ijk}$ are the GARK coupling coefficients, $NΓ$ is the number of sets of GARK coupling coefficients there are $Δc_s = \sum_{j=1}^{Nstages} γ_{sj1} = c_{s+1} - c_s$ is the scaling increment between stage times. The ODE for $v(t)$ is solved using the fastsolver. Note that this form of the scheme is based on Definition 2.2 of Sandu (2019), but ODE for $v(t)$ is written to go from $t_s$ to $T_i + c_s * Δt$ as opposed to $0$ to $1$.

Currently only LowStorageRungeKutta2N schemes are supported for fastsolver

The coefficients defined by γ̂s can be used for an embedded scheme (only the last stage is different).

The available concrete implementations are:

References

@article{Sandu2019,
    title={A class of multirate infinitesimal gark methods},
    author={Sandu, Adrian},
    journal={SIAM Journal on Numerical Analysis},
    volume={57},
    number={5},
    pages={2300--2327},
    year={2019},
    publisher={SIAM},
    doi={10.1137/18M1205492}
}
source
ClimateMachine.ODESolvers.MRIGARKESDIRK23LSAFunction
MRIGARKESDIRK23LSA(f!, fastsolver, Q; dt, t0 = 0, δ = 0

A 2nd order, 3 stage decoupled implicit scheme. It is based on L-Stable, stiffly-accurate ESDIRK scheme of Bank et al (1985); see also Kennedy and Carpenter (2016).

The free parameter δ can take any values for accuracy.

References

@article{Bank1985,
    title={Transient simulation of silicon devices and circuits},
    author={R. E. Bank and W. M. Coughran and W. Fichtner and
            E. H. Grosse and D. J. Rose and R. K. Smith},
    journal={IEEE Transactions on Computer-Aided Design of Integrated
             Circuits and Systems},
    volume={4},
    number={4},
    pages={436-451},
    year={1985},
    publisher={IEEE},
    doi={10.1109/TCAD.1985.1270142}
}

@techreport{KennedyCarpenter2016,
    title = {Diagonally implicit Runge-Kutta methods for ordinary
             differential equations. A review},
             author = {C. A. Kennedy and M. H. Carpenter},
    institution = {National Aeronautics and Space Administration},
    year = {2016},
    number = {NASA/TM–2016–219173},
    address = {Langley Research Center, Hampton, VA}
}
source
ClimateMachine.ODESolvers.MRIGARKERK33aSanduFunction
MRIGARKERK33aSandu(f!, fastsolver, Q; dt, t0 = 0, δ = -1 // 2)

The 3rd order, 3 stage scheme from Sandu (2019). The parameter δ defaults to the value suggested by Sandu, but can be varied.

source
ClimateMachine.ODESolvers.MRIGARKDecoupledImplicitType
MRIGARKDecoupledImplicit(f!, backward_euler_solver, fastsolver, Γs, γ̂s, Q,
                         Δt, t0)

Construct a decoupled implicit MultiRate Infinitesimal General-structure Additive Runge–Kutta (MRI-GARK) scheme to solve

\[ \dot{y} = f(y, t) + g(y, t)\]

where f is the slow tendency function and g is the fast tendency function; see Sandu (2019).

The fast tendency is integrated using the fastsolver and the slow tendency using the MRI-GARK scheme. Since this is a decoupled, implicit MRI-GARK there is no implicit coupling between the fast and slow tendencies.

The backward_euler_solver should be of type AbstractBackwardEulerSolver or LinearBackwardEulerSolver, and is used to perform the backward Euler solves for y given the slow tendency function, namely

\[ y = z + α f(y, t; p)\]

Currently only LowStorageRungeKutta2N schemes are supported for fastsolver

The coefficients defined by γ̂s can be used for an embedded scheme (only the last stage is different).

The available concrete implementations are:

References

@article{Sandu2019,
    title={A class of multirate infinitesimal gark methods},
    author={Sandu, Adrian},
    journal={SIAM Journal on Numerical Analysis},
    volume={57},
    number={5},
    pages={2300--2327},
    year={2019},
    publisher={SIAM},
    doi={10.1137/18M1205492}
}
source

Euler methods

ClimateMachine.ODESolvers.LinearBackwardEulerSolverType
LinearBackwardEulerSolver(::AbstractSystemSolver; isadjustable = false)

Helper type for specifying building a backward Euler solver with a linear solver. If isadjustable == true then the solver can be updated with a new time step size.

source
ClimateMachine.ODESolvers.NonLinearBackwardEulerSolverType
struct NonLinearBackwardEulerSolver{NLS}
    nlsolver::NLS
    isadjustable::Bool
    preconditioner_update_freq::Int64
end

Helper type for specifying building a nonlinear backward Euler solver with a nonlinear solver.

Arguments

  • nlsolver: iterative nonlinear solver, i.e., JacobianFreeNewtonKrylovSolver
  • isadjustable: TODO not used, might use for updating preconditioner
  • preconditioner_update_freq: relavent to Jacobian free -1: no preconditioner; positive number, update every freq times
source

Differential Equations

ClimateMachine.ODESolvers.DiffEqJLIMEXSolverType
DiffEqJLSolver(f, RKA, RKB, RKC, Q; dt, t0 = 0)

This is a time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f_I(Q, t) + f_E(Q, t)\]

via a DifferentialEquations.jl DEAlgorithm, which includes support for OrdinaryDiffEq.jl, Sundials.jl, and more.

source
ClimateMachine.ODESolvers.DiffEqJLSolverType
DiffEqJLSolver(f, RKA, RKB, RKC, Q; dt, t0 = 0)

This is a time stepping object for explicitly time stepping the differential equation given by the right-hand-side function f with the state Q, i.e.,

\[ \dot{Q} = f(Q, t)\]

via a DifferentialEquations.jl DEAlgorithm, which includes support for OrdinaryDiffEq.jl, Sundials.jl, and more.

source

ODE Solvers

ClimateMachine.ODESolvers.solve!Function
solve!(Q, solver::AbstractODESolver; timeend,
       stopaftertimeend=true, numberofsteps, callbacks)

Solves an ODE using the solver starting from a state Q. The state Q is updated inplace. The final time timeend or numberofsteps must be specified.

A series of optional callback functions can be specified using the tuple callbacks; see the GenericCallbacks module.

source

Generic Callbacks

ClimateMachine.GenericCallbacksModule
GenericCallbacks

This module defines interfaces and wrappers for callbacks to be used with an AbstractODESolver.

A callback cb defines three methods:

  • GenericCallbacks.init!(cb, solver, Q, param, t), to be called at solver initialization.

  • GenericCallbacks.call!(cb, solver, Q, param, t), to be called after each time step: the return value dictates what action should be taken:

    • 0 or nothing: continue time stepping as usual
    • 1: stop time stepping after all callbacks have been executed
    • 2: stop time stepping immediately
  • GenericCallbacks.fini!(cb, solver, Q, param, t), to be called at solver finish.

Additionally, wrapper callbacks can be used to execute the callbacks under certain conditions:

For convenience, the following objects can also be used as callbacks:

  • A Function object f, init! and fini! are no-ops, and call! will call f(), and ignore the return value.
  • A Tuple object will call init!, call! and fini! on each element of the tuple.
source