CFTimeSchemes
Documentation for CFTimeSchemes.
CFTimeSchemes.ARK_TRBDF2CFTimeSchemes.BackwardEulerCFTimeSchemes.ForwardBackwardEulerCFTimeSchemes.IVPSolverCFTimeSchemes.KinnmarkGrayCFTimeSchemes.MidpointCFTimeSchemes.RungeKutta4CFTimeSchemes.TRBDF2CFTimeSchemes.Update.manageCFTimeSchemes.advance!CFTimeSchemes.advance!CFTimeSchemes.max_time_stepCFTimeSchemes.model_dstateCFTimeSchemes.scratch_spaceCFTimeSchemes.tendencies!CFTimeSchemes.update!
CFTimeSchemes.ARK_TRBDF2 — Type
scheme = ARK_TRBDF2(model)Second-order IMEX scheme from Giraldo et al. 2013. Implicit terms of model are integrated with a TRBDF2 scheme and explicit terms with a 3-stage second-order RK scheme. Pass scheme to IVPSolver.
CFTimeSchemes.BackwardEuler — Type
scheme = BackwardEuler(model)First-order backward Euler scheme for model. Pass scheme to IVPSolver.
CFTimeSchemes.ForwardBackwardEuler — Type
scheme = ForwardBackwardEuler(model)First-order IMEX scheme for model: backward Euler for implicit terms followed by forward Euler for explicit terms. Pass scheme to IVPSolver.
CFTimeSchemes.IVPSolver — Type
solver = IVPSolver(scheme, dt) # non-mutating
solver = IVPSolver(scheme, dt, state0, time0) # mutatingReturn initial-value problem solver solver, to be used with advance!.
The first syntax returns a non-mutating solver, known to work with Zygote but which allocates. When state0 and time0 are provided, scratch spaces are allocated and a mutating solver is returned. advance should then not allocate and have better performance due to memory reuse.
state0 and time0 are used only to allocate scratch spaces. Provided zero(time0) is implemented., a special value of time0 (with user-defined type) may be passed to avoid the needless computation of tendencies, see scratch_space and model_dstate.
CFTimeSchemes.KinnmarkGray — Type
scheme = KinnmarkGray{2, 5}(model)second-order, 5-stage Runge Kutta scheme for model. Pass scheme to IVPSolver. This scheme has a maximum Courant number of 4 for imaginary eigenvalues, which is the best that can be achieved with a 5-stage RK scheme (Kinnmark & Gray, 1984a, https://doi.org/10.1016/0378-4754(84)90039-9).
scheme = KinnmarkGray{3, 5}(model) third-order, 5-stage Runge Kutta scheme for model. Pass scheme to IVPSolver. This scheme has a maximum Courant number of √15 for imaginary eigenvalues. (Kinnmark & Gray, 1984b, https://doi.org/10.1016/0378-4754(84)90056-9).
See also Guba et al. 2020 https://doi.org/10.5194/gmd-13-6467-2020 .
CFTimeSchemes.Midpoint — Type
scheme = Midpoint(model)Midpoint rule, second-order: forward Euler scheme for 1/2 time step followed by backward Euler scheme for 1/2 time step. Pass scheme to IVPSolver.
CFTimeSchemes.RungeKutta4 — Type
scheme = RungeKutta4(model)4th-order Runge Kutta scheme for model. Pass scheme to IVPSolver.
CFTimeSchemes.TRBDF2 — Type
scheme = TRBDF2(model)Three-stage, second-order, L-stable scheme with two implicit stages. Pass scheme to IVPSolver.
CFTimeSchemes.advance! — Function
future, t = advance!(future, scratch, scheme, present, t, dt) # mutating
future, t = advance!(void, void, scheme, present, t, dt) # non-mutatingIntegrate in time by one time step, respectively mutating (non-allocating) and non-mutating
CFTimeSchemes.advance! — Method
future, t = advance!(future, solver, present, t, N)
future, t = advance!(void, solver, present, t, N)CFTimeSchemes.max_time_step — Function
dt = max_time_step(scheme, dt_lim)Return the max time step dt for time-stepping scheme, assuming that stability is limited by imaginary eigenvalues. dt_lim is the smallest time scale (inverse of largest pulsation for a linear system) in the model.
CFTimeSchemes.model_dstate — Method
k = model_dstate(model, state, time)Returns an object k in which we can store tendencies. Argument state is the model state. Its type may contain information needed to allocate k, e.g. arrays may be of eltype ForwardDiff.Dual.
The fallback implementation looks like: k, _ = tendencies(void, void, model, state, time) which incurs the needless computation of tendencies. If this behavior is undesirable, one may pass a special value for time such as nothing and implement a specialized version of tendencies! which only allocates and skips computations.
CFTimeSchemes.scratch_space — Method
scratch = scratch_space(model, state, time)
scratch = scratch_space(scheme, state, time)Return scratch space scratch, to be used later to compute tendencies for model or to hold sub-stages of Runge-Kutta scheme scheme.
Returns an object k in which we can store tendencies. Argument state is the model state. Its type may contain information needed to allocate k, e.g. arrays may be of eltype ForwardDiff.Dual.
For a model, the fallback implementation looks like: _, scratch = tendencies(void, void, model, state, time) which incurs the needless computation of tendencies. If this behavior is undesirable, one may pass a special value for time such as nothing and implement a specialized version of tendencies! which only allocates and skips computations.
CFTimeSchemes.tendencies! — Function
# This variant is called by explicit time schemes.
dstate, scratch = tendencies!(dstate, scratch, model, state, time) # Mutating
dstate, scratch = tendencies!(void, void, model, state, time) # Non-mutatingReturn tendencies dstate and scratch space scratch for a certain model, state and time. Pass void as output arguments for non-mutating variant.
# This variant is called by diagonally-implicit time schemes.
dstate, scratch = tendencies!(dstate, scratch, model, state, time, tau) # Mutating
dstate, scratch = tendencies!(void, void, model, state, time, tau) # Non-mutatingPerform a backward Euler time step from time time to time+tau then return tendencies dstate evaluated at time+tau and scratch space scratch for a certain model and state. tau must be non-negative and may be zero. Pass void as output arguments for non-mutating variant.
# This variant is called by IMEX time schemes.
dstate_exp, dstate_imp, scratch = tendencies!(dstate_exp, dstate_imp, scratch, model, state, time, tau) # Mutating
dstate_exp, dstate_imp, scratch = tendencies!(void, void, void, model, state, time, tau) # Non-mutatingPerform a backward Euler time step of length tau for implicit tendencies and return explicit tendencies dstate_exp and implicit tendencies dtstate_imp, all evaluated at time time+tau. Also return scratch space scratch. Pass void as output arguments for non-mutating variant.
This function is not implemented. It is meant to be implemented by the user for user-defined type Model.
CFTimeSchemes.update! — Method
new = update!(new, model, new, factor1, increment1, [factor2, increment2, ...])
new = update!(new, model, old, factor1, increment1, [factor2, increment2, ...])
new = update!(void, model, old, factor1, increment1, [factor2, increment2, ...])Return newmodel state, obtained by updating theoldstate byfactor1*increment1. Extra argumentsfactor2, increment2, ...can be appended. Allocatesnewifnew::Void. Acts in-place ifnew==old`.
The provided implementation ignores argument model and calls Update.update!(new, nothing, old, factor, increment) which operates recursively on nested tuples / named tuples and looks like new = @. new = old + factor*increment for arrays.
If a different behavior is desired, one may specialize update! for a specific type of model!.
See also Update.manage
CFTimeSchemes.Update.manage — Method
managed_x = manage(x, mgr)Return a decorated object managed_x to be used as the l.h.s of a broadcasted assignment in place of x. This function is called by update!(mgr, ...) and implemented for mgr==nothing as: managed(x, ::Nothing) = x
More interesting behavior can be obtained by (i) passing a user-defined mgr to Update.update!, (ii) implementing a specialized method for Update.manage and (iii) implement special broadcasting behavior for the object managed_x.
For instance: CFTimesSchemes.Update.manage(x, mgr::LoopManager) = mgr[x] where LoopManager is a type provided by package LoopManagers, which also implements the broadcasting machinery.