API

Initial conditions

General

ClimaAtmos.InitialConditions.InitialConditionType
InitialCondition

A mechanism for specifying the LocalState of an AtmosModel at every point in the domain. Given some initial_condition, calling initial_condition(params) returns a function of the form local_state(local_geometry)::LocalState.

source
ClimaAtmos.InitialConditions.hydrostatic_pressure_profileFunction
hydrostatic_pressure_profile(; thermo_params, p_0, [T, θ, q_tot, z_max])

Solves the initial value problem p'(z) = -g * ρ(z) for all z ∈ [0, z_max], given p(0), either T(z) or θ(z), and optionally also q_tot(z). If q_tot(z) is not given, it is assumed to be 0. If z_max is not given, it is assumed to be 30 km. Note that z_max should be the maximum elevation to which the specified profiles T(z), θ(z), and/or q_tot(z) are valid.

source

Plane / Box

Sphere

Cases from literature

Implicit Solver

ClimaAtmos.ImplicitEquationJacobianType
ImplicitEquationJacobian(
    Y, atmos;
    approximate_solve_iters, diffusion_flag, topography_flag, sgs_advection_flag, transform_flag
)

A wrapper for the matrix $∂E/∂Y$, where $E(Y)$ is the "error" of the implicit step with the state $Y$.

Background

When we use an implicit or split implicit-explicit (IMEX) timestepping scheme, we end up with a nonlinear equation of the form $E(Y) = 0$, where

\[ E(Y) = Y_{imp}(Y) - Y = \hat{Y} + Δt * T_{imp}(Y) - Y.\]

In this expression, $Y_{imp}(Y)$ denotes the state at some time $t + Δt$. This can be expressed as the sum of $\hat{Y}$, the contribution from the state at time $t$ (and possibly also at earlier times, depending on the order of the timestepping scheme), and $Δt * T_{imp}(Y)$, the contribution from the implicit tendency $T_{imp}$ between times $t$ and $t + Δt$. The new state at the end of each implicit step in the timestepping scheme is the value of $Y$ that solves this equation, i.e., the value of $Y$ that is consistent with the state $Y_{imp}(Y)$ predicted by the implicit step.

Note: When we use a higher-order timestepping scheme, the full step $Δt$ is divided into several sub-steps or "stages", where the duration of stage $i$ is $Δt * γ_i$ for some constant $γ_i$ between 0 and 1.

In order to solve this equation using Newton's method, we must specify the derivative $∂E/∂Y$. Since $\hat{Y}$ does not depend on $Y$ (it is only a function of the state at or before time $t$), this derivative is

\[ E'(Y) = Δt * T_{imp}'(Y) - I.\]

In addition, we must specify how to divide $E(Y)$ by this derivative, i.e., how to solve the linear equation

\[ E'(Y) * ΔY = E(Y).\]

Note: This equation comes from assuming that there is some $ΔY$ such that $E(Y - ΔY) = 0$ and making the first-order approximation

\[ E(Y - ΔY) \approx E(Y) - E'(Y) * ΔY.\]

After initializing $Y$ to $Y[0] = \hat{Y}$, Newton's method executes the following steps:

  • Compute the derivative $E'(Y[0])$.
  • Compute the implicit tendency $T_{imp}(Y[0])$ and use it to get $E(Y[0])$.
  • Solve the linear equation $E'(Y[0]) * ΔY[0] = E(Y[0])$ for $ΔY[0]$.
  • Update $Y$ to $Y[1] = Y[0] - ΔY[0]$.

If the number of Newton iterations is limited to 1, this new value of $Y$ is taken to be the solution of the implicit equation. Otherwise, this sequence of steps is repeated, i.e., $ΔY[1]$ is computed and used to update $Y$ to $Y[2] = Y[1] - ΔY[1]$, then $ΔY[2]$ is computed and used to update $Y$ to $Y[3] = Y[2] - ΔY[2]$, and so on. The iterative process is terminated either when the error $E(Y)$ is sufficiently close to 0 (according to the convergence condition passed to Newton's method), or when the maximum number of iterations is reached.

Arguments

  • Y::FieldVector: the state of the simulation
  • atmos::AtmosModel: the model configuration
  • approximate_solve_iters::Int: number of iterations to take for the approximate linear solve required when diffusion_flag = UseDerivative()
  • diffusion_flag::DerivativeFlag: whether the derivative of the diffusion tendency with respect to the quantities being diffused should be computed or approximated as 0; must be either UseDerivative() or IgnoreDerivative() instead of a Bool to ensure type-stability
  • topography_flag::DerivativeFlag: whether the derivative of vertical contravariant velocity with respect to horizontal covariant velocity should be computed or approximated as 0; must be either UseDerivative() or IgnoreDerivative() instead of a Bool to ensure type-stability
  • sgs_advection_flag::DerivativeFlag: whether the derivative of the subgrid-scale advection tendency with respect to the updraft quantities should be computed or approximated as 0; must be either UseDerivative() or IgnoreDerivative() instead of a Bool to ensure type-stability
  • transform_flag::Bool: whether the error should be transformed from $E(Y)$ to $E(Y)/Δt$, which is required for non-Rosenbrock timestepping schemes from OrdinaryDiffEq.jl
source

Helper

ClimaAtmos.InitialConditions.ColumnInterpolatableFieldType
ColumnInterpolatableField(::Fields.ColumnField)

A column field object that can be interpolated in the z-coordinate. For example:

cif = ColumnInterpolatableField(column_field)
z = 1.0
column_field_at_z = cif(z)
Warn

This function allocates and is not GPU-compatible so please avoid using this inside step! only use this for initialization.

source