Implicit-Explicit (IMEX) Additively-Partitioned Runge-Kutta Timestepping

In this tutorial, we shall explore the use of IMplicit-EXplicit (IMEX) methods for the solution of nonautonomous (or non time-invariant) equations. For our model problem, we shall reuse the acoustic wave test in the GCM configuration. See its code for details on the model and parameters. For the purposes of this tutorial, we will only run the experiment for a total of 3600 simulation seconds. Details on this test case can be found in Sec. 4.3 of Francis X Giraldo , James F Kelly , Emil M Constantinescu (2013).

using ClimateMachine
const clima_dir = dirname(dirname(pathof(ClimateMachine)));
include(joinpath(
    clima_dir,
    "tutorials",
    "Numerics",
    "TimeStepping",
    "tutorial_acousticwave_config.jl",
));

The acoustic wave test case used in this tutorial represents a global-scale problem with inertia-gravity waves traveling around the entire planet. It has a hydrostatically balanced initial state that is given a pressure perturbation. This initial pressure perturbation causes an acoustic wave to travel to the antipode, coalesce, and return to the initial position. The exact solution of this test case is simple in that the (linear) acoustic theory allows one to verify the analytic speed of sound based on the thermodynamics variables. The initial condition is defined as a hydrostatically balanced atmosphere with background (reference) potential temperature.

To fully demonstrate the advantages of using an IMEX scheme over fully explicit schemes, we start here by going over a simple, fully explicit scheme. The reader can refer to the Single-rate Explicit Timestepping tutorial for detailes on such schemes. Here we use the the 14-stage LSRK method LSRK144NiegemannDiehlBusch, which contains the largest stability region of the low-storage methods available in ClimateMachine.jl.

FT = Float64
timeend = FT(100)

ode_solver = ClimateMachine.ExplicitSolverType(
    solver_method = LSRK144NiegemannDiehlBusch,
);

In the following example, the timestep calculation is based on the CFL condition for horizontally-propogating acoustic waves. We use a Courant number $C = 0.002$ (denoted by CFL in the code bellow) in the horizontal, which corresponds to a timestep size of approximately $1$ second.

CFL = FT(0.002)
cfl_direction = HorizontalDirection()
run_acousticwave(ode_solver, CFL, cfl_direction, timeend);
ClimateMachine.array_type() = Array
┌ Info: Model composition
│     physics = ClimateMachine.Atmos.AtmosPhysics{Float64,Main.##393.EarthParameterSet,ClimateMachine.Atmos.HydrostaticState{Thermodynamics.TemperatureProfiles.DecayingTemperatureProfile{Float64},Float64},ClimateMachine.Atmos.TotalEnergyModel,ClimateMachine.Atmos.DryModel,ClimateMachine.Atmos.Compressible,ClimateMachine.TurbulenceClosures.ConstantDynamicViscosity{Float64,ClimateMachine.TurbulenceClosures.WithoutDivergence},ClimateMachine.TurbulenceConvection.NoTurbConv,ClimateMachine.TurbulenceClosures.NoHyperDiffusion,ClimateMachine.TurbulenceClosures.NoViscousSponge,ClimateMachine.Atmos.NoPrecipitation,ClimateMachine.Atmos.NoRadiation,ClimateMachine.Atmos.NoTracers,ClimateMachine.Atmos.NoLSForcing}(Main.##393.EarthParameterSet(), ClimateMachine.Atmos.HydrostaticState{Thermodynamics.TemperatureProfiles.DecayingTemperatureProfile{Float64},Float64}(Thermodynamics.TemperatureProfiles.DecayingTemperatureProfile{Float64}(300.0, 300.0, 8776.832091407436), 0.0, true), ClimateMachine.Atmos.TotalEnergyModel(), ClimateMachine.Atmos.DryModel(), ClimateMachine.Atmos.Compressible(), ClimateMachine.TurbulenceClosures.ConstantDynamicViscosity{Float64,ClimateMachine.TurbulenceClosures.WithoutDivergence}(0.0, ClimateMachine.TurbulenceClosures.WithoutDivergence()), ClimateMachine.TurbulenceConvection.NoTurbConv(), ClimateMachine.TurbulenceClosures.NoHyperDiffusion(), ClimateMachine.TurbulenceClosures.NoViscousSponge(), ClimateMachine.Atmos.NoPrecipitation(), ClimateMachine.Atmos.NoRadiation(), ClimateMachine.Atmos.NoTracers(), ClimateMachine.Atmos.NoLSForcing())
│     problem = ClimateMachine.Atmos.AtmosProblem{Tuple{ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC},ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC}},Main.##393.AcousticWaveSetup{Float64},typeof(ClimateMachine.Atmos.atmos_problem_init_state_auxiliary)}((ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC}(ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip}(ClimateMachine.Atmos.FreeSlip()), ClimateMachine.Atmos.Insulating(), ClimateMachine.Atmos.Impermeable(), ClimateMachine.Atmos.OutflowPrecipitation(), ClimateMachine.Atmos.ImpermeableTracer(), ClimateMachine.TurbulenceConvection.NoTurbConvBC()), ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC}(ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip}(ClimateMachine.Atmos.FreeSlip()), ClimateMachine.Atmos.Insulating(), ClimateMachine.Atmos.Impermeable(), ClimateMachine.Atmos.OutflowPrecipitation(), ClimateMachine.Atmos.ImpermeableTracer(), ClimateMachine.TurbulenceConvection.NoTurbConvBC())), Main.##393.AcousticWaveSetup{Float64}(10000.0, 300.0, 3.0, 100.0, 1), ClimateMachine.Atmos.atmos_problem_init_state_auxiliary)
│     orientation = ClimateMachine.Orientations.SphericalOrientation()
│     source = DispatchedTuples.DispatchedTuple{Tuple{Tuple{ClimateMachine.Atmos.Momentum,ClimateMachine.Atmos.Gravity}},DispatchedTuples.NoDefaults} with 1 entries:
│   ClimateMachine.Atmos.Momentum() => ClimateMachine.Atmos.Gravity()
│   default => ()
│ 
│     data_config = nothing
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/driver_configs.jl:188

PDE: ∂_t Y_i + (∇•F_1(Y))_i + (∇•F_2(Y,G)))_i = (S(Y,G))_i
┌──────────┬────────────────────────────┬─────────────────────────────────┬───────────┐
│ Equation │           Flux{FirstOrder} │               Flux{SecondOrder} │    Source │
│    (Y_i) │                      (F_1) │                           (F_2) │       (S) │
├──────────┼────────────────────────────┼─────────────────────────────────┼───────────┤
│     Mass │                   (Advect) │                              () │        () │
│ Momentum │ (Advect, PressureGradient) │                 (ViscousStress) │ (Gravity) │
│   Energy │         (Advect, Pressure) │ (ViscousFlux, DiffEnthalpyFlux) │        () │
└──────────┴────────────────────────────┴─────────────────────────────────┴───────────┘

┌ Info: Establishing Atmos GCM configuration for GCM Driver: Acoustic wave test
│     precision               = Float64
│     horiz polynomial order  = 4
│     vert polynomial order   = 4
│     horiz cutoff order      = nothing
│     vert cutoff order       = nothing
│     # horiz elem            = 6
│     # vert elems            = 4
│     domain height           = 1.00e+04 m
│     MPI ranks               = 1
│     min(Δ_horz)             = 203711.48 m
│     min(Δ_vert)             = 431.68 m
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/driver_configs.jl:418
┌ Info: Initializing GCM Driver: Acoustic wave test
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/solver_configs.jl:185
┌ Info: Starting GCM Driver: Acoustic wave test
│     dt              = 1.16279e+00
│     timeend         =   100.00
│     number of steps = 86
│     norm(Q)         = 9.5073559986777266e+13
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/Driver.jl:802
┌ Info: Update
│     simtime =    15.12 /   100.00
│     wallclock = 00:01:00
│     efficiency (simtime / wallclock) =   0.2493
│     wallclock end (estimated) = 00:06:41
│     norm(Q) = 9.5073396077964141e+13
└ @ ClimateMachine.Callbacks /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/Callbacks/Callbacks.jl:75
┌ Info: Update
│     simtime =    80.23 /   100.00
│     wallclock = 00:02:00
│     efficiency (simtime / wallclock) =   0.6652
│     wallclock end (estimated) = 00:02:30
│     norm(Q) = 9.5073418293448797e+13
└ @ ClimateMachine.Callbacks /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/Callbacks/Callbacks.jl:75
┌ Info: Finished
│     norm(Q)            = 9.5073405874321281e+13
│     norm(Q) / norm(Q₀) = 9.9999837901877231e-01
│     norm(Q) - norm(Q₀) = -1.5411245598437500e+08
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/Driver.jl:853

However, as it is imaginable, for real-world climate processes a time step of 1 second would lead to extemely long time-to-solution simulations. How can we do better? To be able to take larger time step, we can treat the most restrictive wave speeds (vertical acoustic) implicitly rather than explicitly. This motivates the use of an IMplicit-EXplicit (IMEX) methods.

In general, a single step of an $s$-stage, $N$-part additive RK method (ARK_N) is defined by its generalized Butcher tableau:

\[\begin{align} \begin{array}{c|c|c|c} \boldsymbol{c} &\boldsymbol{A}_{1} & \cdots & \boldsymbol{A}_{N}\\ \hline & \boldsymbol{b}_1^T & \cdots & \boldsymbol{b}_N^T\\ \hline & \widehat{\boldsymbol{b}}_1^T & \cdots & \widehat{\boldsymbol{b}}_N^T \end{array} = \begin{array}{c|c c c | c | c c c } c_1 & a^{[ 1 ]}_{1,1} & \cdots & a^{[ 1 ]}_{1,s} & \cdots & a^{[ \nu ]}_{1,1} & \cdots & a^{[ \nu ]}_{1,s}\\ \vdots & \vdots & \ddots & \vdots & \cdots & \vdots & \ddots & \vdots \\ c_s & a^{[ 1 ]}_{s,1} & \cdots & a^{[ 1 ]}_{s,s} & \cdots & a^{[ \nu ]}_{s,1} & \cdots & a^{[ \nu ]}_{s,s}\\ \hline & b^{[ 1 ]}_1 & \cdots & b^{[ 1 ]}_s & \cdots & b^{[ \nu ]}_1 & \cdots & b^{[ \nu ]}_s\\ \hline & \widehat{b}^{[ 1 ]}_1 & \cdots & \widehat{b}^{[ 1 ]}_s & & \widehat{b}^{[ \nu ]}_1 & \cdots & \widehat{b}^{[ \nu ]}_s \end{array} \end{align}\]

and is given by

$\boldsymbol{q}^{n+1} = \boldsymbol{q}^n + \Delta t \left( \underbrace{\sum_{i=1}^{s}}_{\textrm{Stages}} \underbrace{\sum_{\nu=1}^{N}}_{\textrm{Components}} b_i^{[ \nu ]} {\mathcal{T}}^{[ \nu ]}(\boldsymbol{Q}^i)) \right)$

where $s$ denotes the stages and $N$ the components, and where the stage values are given by:

$\boldsymbol{Q}^i = \boldsymbol{q}^n + \Delta t \sum_{j=1}^{s} \sum_{\nu = 1}^{N} a_{i,j}^{[ \nu ]} {\mathcal{T}}^{[ \nu]}(\boldsymbol{Q}^j).$

Similar to standard RK methods, the stage vectors are approximations to the state at each stage of the ARK method. Moreover, the temporal coefficients $c_i$ satisfy a similar row-sum condition, holding for all $\nu = 1, \cdots, N$:

$c_i = \sum_{j=1}^{s} a_{i, j}^{[ \nu ]}, \quad \forall \nu = 1, \cdots, N.$

The Butcher coefficients $\boldsymbol{c}$, $\boldsymbol{b}_{\nu}$, $\boldsymbol{A}_{\nu}$, and $\widehat{\boldsymbol{b}}_{\nu}$ are constrained by certain accuracy and stability requirements, which are summarized in Christopher Alan Kennedy (2001).

A common setting is the case $N = 2$. This gives the typical context for Implicit-Explicit (IMEX) splitting methods, where the tendency ${\mathcal{T}}$ is assumed to have the decomposition:

$\dot{\boldsymbol{q}} = \mathcal{T}(\boldsymbol{q}) \equiv {\mathcal{T}}_{s}(\boldsymbol{q}) + {\mathcal{T}}_{ns}(\boldsymbol{q}),$ where the right-hand side has been split into a "stiff" component ${\mathcal{T}}_{s}$, to be treated implicitly, and a non-stiff part ${\mathcal{T}}_{ns}$ to be treated explicitly.

Referencing the canonical form introduced in Time integration we have that in this particular forumlation $\mathcal{T}_{ns}(t, \boldsymbol{q}) \equiv \mathcal{G}(t, \boldsymbol{q})$ and $\mathcal{T}_{s}(t, \boldsymbol{q}) \equiv \mathcal{F}(t, \boldsymbol{q})$.

Two different RK methods are applied to ${\mathcal{T}}_{s}$ and ${\mathcal{T}}_{ns}$ separately, which have been specifically designed and coupled. Examples can be found in Francis X Giraldo , James F Kelly , Emil M Constantinescu (2013). The Butcher Tableau for an ARK_2 method will have the form

\[\begin{align} \begin{array}{c|c|c} \boldsymbol{c} &\boldsymbol{A}_E &\boldsymbol{A}_I\\ \hline & \boldsymbol{b}_E^T & \boldsymbol{b}_I^T \\ \hline & \widehat{\boldsymbol{b}}_E^T & \widehat{\boldsymbol{b}}_I^T \end{array}, \end{align}\]

with

$\boldsymbol{A}_O = \left\lbrace a_{i, j}^O \right\rbrace, \quad \boldsymbol{b}_O = \left\lbrace b_{i}^O \right\rbrace, \quad \widehat{\boldsymbol{b}}_O = \left\lbrace \widehat{b}_{i}^O \right\rbrace,$

where $O$ denotes the label (either $E$ for explicit or $I$ for implicit).

For the acoustic wave example used here, we use 4th order polynomials in our discontinuous Galerkin approximation, with 6 elements in each horizontal direction and 4 elements in the vertical direction, on the cubed-sphere. This gives an effective minimal node-distance (distance between LGL nodes) of roughly 203000 m. As in the previous tutorial, we can determine our $\Delta t$ by specifying our desired horizontal Courant number $C$ (the timestep calculation is based on the CFL condition for horizontally-propogating acoustic waves). In this very simple test case, we can use a value of 0.5, which corresponds to a time-step size of around 257 seconds. But for this particular example, even higher values might work.

timeend = FT(3600)
ode_solver = ClimateMachine.IMEXSolverType(
    solver_method = ARK2GiraldoKellyConstantinescu,
)
CFL = FT(0.5)
cfl_direction = HorizontalDirection()
run_acousticwave(ode_solver, CFL, cfl_direction, timeend);
ClimateMachine.array_type() = Array
┌ Info: Model composition
│     physics = ClimateMachine.Atmos.AtmosPhysics{Float64,Main.##393.EarthParameterSet,ClimateMachine.Atmos.HydrostaticState{Thermodynamics.TemperatureProfiles.DecayingTemperatureProfile{Float64},Float64},ClimateMachine.Atmos.TotalEnergyModel,ClimateMachine.Atmos.DryModel,ClimateMachine.Atmos.Compressible,ClimateMachine.TurbulenceClosures.ConstantDynamicViscosity{Float64,ClimateMachine.TurbulenceClosures.WithoutDivergence},ClimateMachine.TurbulenceConvection.NoTurbConv,ClimateMachine.TurbulenceClosures.NoHyperDiffusion,ClimateMachine.TurbulenceClosures.NoViscousSponge,ClimateMachine.Atmos.NoPrecipitation,ClimateMachine.Atmos.NoRadiation,ClimateMachine.Atmos.NoTracers,ClimateMachine.Atmos.NoLSForcing}(Main.##393.EarthParameterSet(), ClimateMachine.Atmos.HydrostaticState{Thermodynamics.TemperatureProfiles.DecayingTemperatureProfile{Float64},Float64}(Thermodynamics.TemperatureProfiles.DecayingTemperatureProfile{Float64}(300.0, 300.0, 8776.832091407436), 0.0, true), ClimateMachine.Atmos.TotalEnergyModel(), ClimateMachine.Atmos.DryModel(), ClimateMachine.Atmos.Compressible(), ClimateMachine.TurbulenceClosures.ConstantDynamicViscosity{Float64,ClimateMachine.TurbulenceClosures.WithoutDivergence}(0.0, ClimateMachine.TurbulenceClosures.WithoutDivergence()), ClimateMachine.TurbulenceConvection.NoTurbConv(), ClimateMachine.TurbulenceClosures.NoHyperDiffusion(), ClimateMachine.TurbulenceClosures.NoViscousSponge(), ClimateMachine.Atmos.NoPrecipitation(), ClimateMachine.Atmos.NoRadiation(), ClimateMachine.Atmos.NoTracers(), ClimateMachine.Atmos.NoLSForcing())
│     problem = ClimateMachine.Atmos.AtmosProblem{Tuple{ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC},ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC}},Main.##393.AcousticWaveSetup{Float64},typeof(ClimateMachine.Atmos.atmos_problem_init_state_auxiliary)}((ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC}(ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip}(ClimateMachine.Atmos.FreeSlip()), ClimateMachine.Atmos.Insulating(), ClimateMachine.Atmos.Impermeable(), ClimateMachine.Atmos.OutflowPrecipitation(), ClimateMachine.Atmos.ImpermeableTracer(), ClimateMachine.TurbulenceConvection.NoTurbConvBC()), ClimateMachine.Atmos.AtmosBC{ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip},ClimateMachine.Atmos.Insulating,ClimateMachine.Atmos.Impermeable,ClimateMachine.Atmos.OutflowPrecipitation,ClimateMachine.Atmos.ImpermeableTracer,ClimateMachine.TurbulenceConvection.NoTurbConvBC}(ClimateMachine.Atmos.Impenetrable{ClimateMachine.Atmos.FreeSlip}(ClimateMachine.Atmos.FreeSlip()), ClimateMachine.Atmos.Insulating(), ClimateMachine.Atmos.Impermeable(), ClimateMachine.Atmos.OutflowPrecipitation(), ClimateMachine.Atmos.ImpermeableTracer(), ClimateMachine.TurbulenceConvection.NoTurbConvBC())), Main.##393.AcousticWaveSetup{Float64}(10000.0, 300.0, 3.0, 100.0, 1), ClimateMachine.Atmos.atmos_problem_init_state_auxiliary)
│     orientation = ClimateMachine.Orientations.SphericalOrientation()
│     source = DispatchedTuples.DispatchedTuple{Tuple{Tuple{ClimateMachine.Atmos.Momentum,ClimateMachine.Atmos.Gravity}},DispatchedTuples.NoDefaults} with 1 entries:
│   ClimateMachine.Atmos.Momentum() => ClimateMachine.Atmos.Gravity()
│   default => ()
│ 
│     data_config = nothing
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/driver_configs.jl:188

PDE: ∂_t Y_i + (∇•F_1(Y))_i + (∇•F_2(Y,G)))_i = (S(Y,G))_i
┌──────────┬────────────────────────────┬─────────────────────────────────┬───────────┐
│ Equation │           Flux{FirstOrder} │               Flux{SecondOrder} │    Source │
│    (Y_i) │                      (F_1) │                           (F_2) │       (S) │
├──────────┼────────────────────────────┼─────────────────────────────────┼───────────┤
│     Mass │                   (Advect) │                              () │        () │
│ Momentum │ (Advect, PressureGradient) │                 (ViscousStress) │ (Gravity) │
│   Energy │         (Advect, Pressure) │ (ViscousFlux, DiffEnthalpyFlux) │        () │
└──────────┴────────────────────────────┴─────────────────────────────────┴───────────┘

┌ Info: Establishing Atmos GCM configuration for GCM Driver: Acoustic wave test
│     precision               = Float64
│     horiz polynomial order  = 4
│     vert polynomial order   = 4
│     horiz cutoff order      = nothing
│     vert cutoff order       = nothing
│     # horiz elem            = 6
│     # vert elems            = 4
│     domain height           = 1.00e+04 m
│     MPI ranks               = 1
│     min(Δ_horz)             = 203711.48 m
│     min(Δ_vert)             = 431.68 m
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/driver_configs.jl:418
┌ Info: Initializing GCM Driver: Acoustic wave test
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/solver_configs.jl:185
┌ Info: Starting GCM Driver: Acoustic wave test
│     dt              = 2.76923e+02
│     timeend         =  3600.00
│     number of steps = 13
│     norm(Q)         = 9.5073559986777266e+13
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/Driver.jl:802
┌ Info: Finished
│     norm(Q)            = 9.5073459060339328e+13
│     norm(Q) / norm(Q₀) = 9.9999893843842647e-01
│     norm(Q) - norm(Q₀) = -1.0092643793750000e+08
└ @ ClimateMachine /central/scratch/climaci/climatemachine-docs/1430/climatemachine-docs/src/Driver/Driver.jl:853

References


This page was generated using Literate.jl.