Driver

Solver types

ClimateMachine.HEVISplittingType
HEVISplitting

HEVI (horizontally explicit, vertically implicit) type method, where vertical acoustic waves are treated implicitly. All other dynamics are treated explicitly.

Note: Can potentially imagine several different types of HEVI splittings (for example, include vertical momentum and/or diffusion)

source
ClimateMachine.MISSolverTypeType

Description

MISSolverType(;
    splitting_type = SlowFastSplitting(),
    fast_model = AtmosAcousticGravityLinearModel,
    mis_method = MIS2,
    fast_method = LSRK54CarpenterKennedy,
    nsubsteps = 50,
)

This solver type constructs an ODE solver using a generalization of the split-explicit Runge-Kutta method. Known as the Multirate Infinitesimal Step (MIS) method, this solver solves ODEs with the partitioned form:

\[ \dot{Q} = f_{fast}(Q, t) + f_{slow}(Q, t)\]

where the right-hand-side functions f_fast and f_slow denote fast and slow dynamics respectively, depending on the state Q.

Arguments

  • splitting_type (DiscreteSplittingType): The type of discrete splitting to apply to the right-hand side. Default: SlowFastSplitting()
  • fast_model (Type): The model describing fast dynamics. Default: AtmosAcousticGravityLinearModel
  • mis_method (Function): Function defining the particular MIS method to be used. Default: MIS2
  • fast_method (Function): Function defining the fast solver. Default: LSRK54CarpenterKennedy
  • nsubsteps (Int): Integer denoting the total number of times to substep the fast process. Default: 50
  • discrete_splitting (Boolean): Boolean denoting whether a PDE level or discretized level splitting should be used. If true then the PDE is discretized in such a way that f_fast + f_slow is equivalent to discretizing the original PDE directly. Default: false

References

@article{KnothWensch2014,
    title={Generalized split-explicit Runge--Kutta methods for
        the compressible Euler equations},
    author={Knoth, Oswald and Wensch, Joerg},
    journal={Monthly Weather Review},
    volume={142},
    number={5},
    pages={2067--2081},
    year={2014}
}
source
ClimateMachine.MultirateSolverTypeType

Description

MultirateSolverType(;
    splitting_type = SlowFastSplitting(),
    fast_model = AtmosAcousticGravityLinearModel,
    implicit_solver = ManyColumnLU,
    implicit_solver_adjustable = false,
    slow_method = LSRK54CarpenterKennedy,
    fast_method = LSRK54CarpenterKennedy,
    timestep_ratio = 100,
)

This solver type constructs an ODE solver using a standard multirate Runge-Kutta implementation. This solver computes solutions to ODEs with the partitioned form:

\[ \dot{Q} = f_{fast}(Q, t) + f_{slow}(Q, t)\]

where the right-hand-side functions f_fast and f_slow denote fast and slow dynamics respectively, depending on the state Q.

Arguments

  • splitting_type (DiscreteSplittingType): The type of discrete splitting to apply to the right-hand side. Default: SlowFastSplitting()
  • fast_model (Type): The model describing fast dynamics. Default: AtmosAcousticGravityLinearModel
  • implicit_solver (Type): An implicit solver for inverting the implicit system of equations (if using HEVISplitting()). Default: ManyColumnLU
  • implicit_solver_adjustable (Bool): A flag identifying whether or not the implicit_solver can be updated as the time-step size changes. This is particularly important when using an implicit solver within a multirate scheme. Default: false
  • slow_method (Function): Function defining the particular explicit Runge-Kutta method to be used for the slow processes. Default: LSRK54CarpenterKennedy
  • fast_method (Function): Function defining the fast solver. Depending on the choice of splitting_type, this can be an explicit Runge Kutta method or a 1-D IMEX (additive Runge-Kutta) method. Default: LSRK54CarpenterKennedy
  • timestep_ratio (Int): Integer denoting the ratio between the slow and fast time-step sizes. Default: 100
  • discrete_splitting (Boolean): Boolean denoting whether a PDE level or discretized level splitting should be used. If true then the PDE is discretized in such a way that f_fast + f_slow is equivalent to discretizing the original PDE directly.

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
ClimateMachine.AbstractSolverTypeType
AbstractSolverType

This is an abstract type representing a generic solver. By a "solver," we mean an ODE solver together with any potential implicit solver (linear solvers).

source
ClimateMachine.DiscreteSplittingTypeType
DiscreteSplittingType

This is an abstract type representing a temporal splitting in the discrete equations. For example, HEVI (horizontally explicit, vertically implicit) type methods.

source
ClimateMachine.ExplicitSolverTypeType

Description

ExplicitSolverType(;
    solver_method = LSRK54CarpenterKennedy,
)

This solver type constructs an ODE solver using an explicit Runge-Kutta method.

Arguments

  • solver_method (Function): Function defining the explicit Runge-Kutta solver. Default: LSRK54CarpenterKennedy
source
ClimateMachine.ImplicitSolverTypeType

Description

ImplicitSolverType(;
    solver_method = KenCarp4,
)

This solver type constructs an ODE solver using a fully implicit method.

Arguments

  • solver_method (Function): Function defining the implicit solver. Default: KenCarp4
source
ClimateMachine.SplitExplicitSolverTypeType

Description

SplitExplicitSolverType

This solver type constructs an ODE solver using the SplitExplicitLSRK2nSolver.

Arguments

  • dt_slow (AbstractFloat): Time step for the slow solver
  • dt_fast (AbstractFloat): Time step for the fast solver
  • slow_method (Function): Function defining the explicit Runge-Kutta solver for the slow model. Default: LSRK54CarpenterKennedy
  • fast_method (Function): Function defining the explicit Runge-Kutta solver for the fast model. Default: LSRK54CarpenterKennedy
source
ClimateMachine.IMEXSolverTypeType

Description

IMEXSolverType(;
    splitting_type = HEVISplitting(),
    implicit_model = AtmosAcousticGravityLinearModel,
    implicit_solver = ManyColumnLU,
    implicit_solver_adjustable = false,
    solver_method = ARK2GiraldoKellyConstantinescu,
    solver_storage_variant = LowStorageVariant(),
    split_explicit_implicit = false,
    discrete_splitting = true,
)

This solver type constructs a solver for ODEs with the additively-partitioned form. 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)] + [n(Q, t)]\]

where n is now only the nonlinear tendency.

Arguments

  • splitting_type (DiscreteSplittingType): The type of discrete splitting to apply to the right-hand side. Default: HEVISplitting()
  • implicit_model (Type): The model describing dynamics to be treated implicitly. Default: AtmosAcousticGravityLinearModel
  • implicit_solver (Type): A solver for inverting the implicit system of equations. Default: ManyColumnLU
  • implicit_solver_adjustable (Bool): A flag identifying whether or not the implicit_solver can be updated as the time-step size changes. Default: false
  • solver_method (Function): Function defining the particular additive Runge-Kutta method to be used for the IMEX method. Default: ARK2GiraldoKellyConstantinescu
  • solver_storage_variant (Type): Storage type for the additive Runge-Kutta method. Default: LowStorageVariant()
  • split_explicit_implicit (Boolean): Whether the tendency is split in explicit and implicit parts or not.
  • discrete_splitting (Boolean): Boolean denoting whether a PDE level or discretized level splitting should be used. If true then the PDE is discretized in such a way that f_fast + f_slow is equivalent to discretizing the original PDE directly.

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

Configurations

ClimateMachine.InterpolationConfigurationFunction
InterpolationConfiguration(
    driver_config::DriverConfiguration,
    boundaries::Array,
    resolution::Tuple,
)

Creates an InterpolationTopology (either an InterpolationBrick or an InterpolationCubedSphere) to be used with a DiagnosticsGroup.

source
ClimateMachine.ConservationCheckType
ClimateMachine.ConservationCheck

Pass a tuple of these to ClimateMachine.invoke! to perform a conservation check of each varname at the specified interval. This computes Σv = weightedsum(Q.varname) and δv = (Σv - Σv₀) / Σv. invoke! throws an error if abs(δv) exceeds error_threshold. Ifshow,δv` is displayed.

source

Initialize / solve

ClimateMachine.array_typeFunction
ClimateMachine.array_type()

Return the array type used by ClimateMachine. This defaults to (CPU-based) Array and is only correctly set (based on choice from the command line, from an environment variable, or from experiment code) after ClimateMachine.init() is called.

source
ClimateMachine.initFunction
ClimateMachine.init(;
    parse_clargs::Bool = false,
    custom_clargs::Union{Nothing, ArgParseSettings} = nothing,
    init_driver::Bool = true,
    keyword_args...,
)

Initialize the ClimateMachine. If parse_clargs is set, parse command line arguments (additional driver-specific arguments can be added by specifying custom_clargs).

Setting init_driver = false will set up the ClimateMachine.Settings singleton values without initializing the ClimateMachine runtime. Otherwise, the runtime will be initialized (see init_runtime()).

Finally, key-value pairs can be supplied to ClimateMachine.init() to set system default settings – the final settings are decided as follows (in order of precedence):

  1. Command line arguments (if parse_clargs = true).
  2. Environment variables.
  3. Keyword arguments to init().
  4. Defaults (in ClimateMachine_Settings).

Recognized keyword arguments are:

  • disable_gpu::Bool = false: do not use the GPU
  • show_updates::String = "60secs": interval at which to show simulation updates
  • diagnostics::String = "never": interval at which to collect diagnostics"
  • vtk::String = "never": inteverval at which to write simulation vtk output
  • vtk-number-sample-points::Int = 0: the number of sampling points in each element for VTK output
  • monitor_timestep_duration::String = "never": interval in time-steps at which to output wall-clock time per time-step
  • monitor_courant_numbers::String = "never": interval at which to output acoustic, advective, and diffusive Courant numbers"
  • checkpoint::String = "never": interval at which to output a checkpoint
  • checkpoint_keep_one::Bool = true: (interval) keep all checkpoints (instead of just the most recent)"
  • checkpoint_at_end::Bool = false: create a checkpoint at the end of the simulation"
  • checkpoint_dir::String = "checkpoint": absolute or relative path to checkpoint directory
  • restart_from_num::Int = -1: checkpoint number from which to restart (in checkpoint_dir)
  • fix_rng_seed::Bool = false: set RNG seed to a fixed value for reproducibility
  • log_level::String = "INFO": log level for ClimateMachine global default runtime logger
  • disable_custom_logger::String = false: disable using a global custom logger for ClimateMachine
  • output_dir::String = "output": (path) absolute or relative path to output data directory
  • debug_init::Bool = false: fill state arrays with NaNs and dump them post-initialization
  • integration_testing::Bool = false: enable integration_testing

Returns nothing, or if parse_clargs = true, returns parsed command line arguments.

source
ClimateMachine.invoke!Function
ClimateMachine.invoke!(
    solver_config::SolverConfiguration;
    adjustfinalstep = false,
    diagnostics_config = nothing,
    user_callbacks = (),
    user_info_callback = () -> nothing,
    check_cons = (),
    check_euclidean_distance = false,
)

Run the simulation defined by solver_config.

Keyword Arguments:

The value of 'adjustfinalstepis passed to the ODE solver; see [solve!`](@ref ODESolvers.solve!).

The user_callbacks are passed to the ODE solver as callback functions; see solve!.

The function user_info_callback is called after the default info callback (which is called every Settings.show_updates interval). The single input argument init is true when the callback is called for initialization (before time stepping begins) and false when called during the actual ODE solve; see GenericCallbacks and solve!.

If conservation checks are to be performed, check_cons must be a tuple of ConservationCheck.

If check_euclidean_distance is true, then the Euclidean distance between the final solution and initial condition function evaluated withsolver_config.timeend` is reported.

source