API
Initial conditions
General
ClimaAtmos.InitialConditions.InitialCondition
— TypeInitialCondition
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
.
ClimaAtmos.InitialConditions.IsothermalProfile
— TypeIsothermalProfile(; temperature = 300)
An InitialCondition
with a uniform temperature profile.
ClimaAtmos.InitialConditions.DecayingProfile
— TypeDecayingProfile(; perturb = true)
An InitialCondition
with a decaying temperature profile, and with an optional perturbation to the temperature.
ClimaAtmos.InitialConditions.hydrostatic_pressure_profile
— Functionhydrostatic_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.
Plane / Box
ClimaAtmos.InitialConditions.ConstantBuoyancyFrequencyProfile
— TypeConstantBuoyancyFrequencyProfile()
An InitialCondition
with a constant Brunt-Vaisala frequency and constant wind velocity, where the pressure profile is hydrostatically balanced. This is currently the only InitialCondition
that supports the approximation of a steady-state solution.
ClimaAtmos.InitialConditions.DryDensityCurrentProfile
— TypeDryDensityCurrentProfile(; perturb = false)
An InitialCondition
with an isothermal background profile, with a negatively buoyant bubble, and with an optional perturbation to the temperature.
ClimaAtmos.InitialConditions.RisingThermalBubbleProfile
— TypeRisingThermalBubbleProfile(; perturb = false)
An InitialCondition
with an isothermal background profile, with a positively buoyant bubble, and with an optional perturbation to the temperature.
Sphere
ClimaAtmos.InitialConditions.DryBaroclinicWave
— TypeDryBaroclinicWave(; perturb = true, deep_atmosphere = false)
An InitialCondition
with a dry baroclinic wave, and with an optional perturbation to the horizontal velocity.
ClimaAtmos.InitialConditions.MoistBaroclinicWaveWithEDMF
— TypeMoistBaroclinicWaveWithEDMF(; perturb = true, deep_atmosphere = false)
The same InitialCondition
as MoistBaroclinicWave
, except with an initial TKE of 0 and an initial draft area fraction of 0.2.
ClimaAtmos.InitialConditions.MoistAdiabaticProfileEDMFX
— TypeMoistAdiabaticProfileEDMFX(; perturb = true)
An InitialCondition
with a moist adiabatic temperature profile, and with an optional perturbation to the temperature.
Cases from literature
ClimaAtmos.InitialConditions.Nieuwstadt
— TypeNieuwstadt
The InitialCondition
described in [7], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.GABLS
— TypeGABLS
The InitialCondition
described in [8], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.GATE_III
— TypeGATE_III
The InitialCondition
described in [9], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.ARM_SGP
— TypeARM_SGP
The InitialCondition
described in [10], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.DYCOMS_RF01
— TypeDYCOMS_RF01
The InitialCondition
described in [11], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.DYCOMS_RF02
— TypeDYCOMS_RF02
The InitialCondition
described in [12], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.Rico
— TypeRico
The InitialCondition
described in [13], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.TRMM_LBA
— TypeTRMM_LBA
The InitialCondition
described in [14], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.LifeCycleTan2018
— TypeLifeCycleTan2018
The InitialCondition
described in [15], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.Bomex
— TypeBomex
The InitialCondition
described in [16], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.Soares
— TypeSoares
The InitialCondition
described in [17], but with a hydrostatically balanced pressure profile.
Helper
ClimaAtmos.InitialConditions.ColumnInterpolatableField
— TypeColumnInterpolatableField(::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)
Jacobian
ClimaAtmos.Jacobian
— TypeJacobian(alg, Y, atmos; [verbose])
Wrapper for a JacobianAlgorithm
and its cache, which it uses to update and invert the Jacobian. The optional verbose
flag specifies whether debugging information should be printed during initialization.
ClimaAtmos.JacobianAlgorithm
— TypeJacobianAlgorithm
A description of how to compute the matrix $∂R/∂Y$, where $R(Y)$ denotes the residual of an implicit step with the state $Y$. Concrete implementations of this abstract type should define 3 methods:
jacobian_cache(alg::JacobianAlgorithm, Y, atmos; [verbose])
update_jacobian!(alg::JacobianAlgorithm, cache, Y, p, dtγ, t)
invert_jacobian!(alg::JacobianAlgorithm, cache, ΔY, R)
To facilitate debugging, concrete implementations should also define
first_column_block_arrays(alg::JacobianAlgorithm, Y, p, dtγ, t)
See Implicit Solver for additional background information.
ClimaAtmos.ManualSparseJacobian
— TypeManualSparseJacobian(
topography_flag,
diffusion_flag,
sgs_advection_flag,
sgs_entr_detr_flag,
sgs_mass_flux_flag,
sgs_nh_pressure_flag,
approximate_solve_iters,
)
A JacobianAlgorithm
that approximates the Jacobian using analytically derived tendency derivatives and inverts it using a specialized nested linear solver. Certain groups of derivatives can be toggled on or off by setting their DerivativeFlag
s to either UseDerivative
or IgnoreDerivative
.
Arguments
topography_flag::DerivativeFlag
: whether the derivative of vertical contravariant velocity with respect to horizontal covariant velocity should be computeddiffusion_flag::DerivativeFlag
: whether the derivatives of the grid-scale diffusion tendency should be computedsgs_advection_flag::DerivativeFlag
: whether the derivatives of the subgrid-scale advection tendency should be computedsgs_entr_detr_flag::DerivativeFlag
: whether the derivatives of the subgrid-scale entrainment and detrainment tendencies should be computedsgs_mass_flux_flag::DerivativeFlag
: whether the derivatives of the subgrid-scale mass flux tendency should be computedsgs_nh_pressure_flag::DerivativeFlag
: whether the derivatives of the subgrid-scale non-hydrostatic pressure drag tendency should be computedapproximate_solve_iters::Int
: number of iterations to take for the approximate linear solve required when thediffusion_flag
isUseDerivative
ClimaAtmos.AutoDenseJacobian
— TypeAutoDenseJacobian([max_simultaneous_derivatives])
A JacobianAlgorithm
that computes the Jacobian using forward-mode automatic differentiation, without making any assumptions about sparsity structure. After the dense matrix for each spatial column is updated, parallel_lu_factorize!
computes its LU factorization in parallel across all columns. The linear solver is also run in parallel with parallel_lu_solve!
.
To automatically compute the derivative of implicit_tendency!
with respect to Y
, we first create copies of Y
, p.precomputed
, and p.scratch
in which every floating-point number is replaced by a dual number from ForwardDiff.jl
. A dual number can be expressed as $Xᴰ = X + ε₁x₁ + ε₂x₂ + ... + εₙxₙ$, where $X$ and $xᵢ$ are floating-point numbers, and where $εᵢ$ is a hyperreal number that satisfies $εᵢεⱼ = 0$. If the $i$-th value in dual column state $Yᴰ$ is set to $Yᴰᵢ = Yᵢ + 1εᵢ$, where $Yᵢ$ is the $i$-th value in the column state $Y$, then evaluating the implicit tendency of the dual column state generates a dense representation of the Jacobian matrix $∂T/∂Y$. Specifically, the $i$-th value in the dual column tendency $Tᴰ = T(Yᴰ)$ is $Tᴰᵢ = Tᵢ + (∂Tᵢ/∂Y₁)ε₁ + ... + (∂Tᵢ/∂Yₙ)εₙ$, where $Tᵢ$ is the $i$-th value in the column tendency $T(Y)$, and where $n$ is the number of values in $Y$. In other words, the entry in the $i$-th row and $j$-th column of the matrix $∂T/∂Y$ is the coefficient of $εⱼ$ in $Tᴰᵢ$. The size of the dense matrix scales as $O(n^2)$, leading to very large memory requirements at higher vertical resolutions.
When the number of values in each column is very large, computing the entire dense matrix in a single evaluation of implicit_tendency!
can be too expensive to compile and run. So, the dual number components are split into partitions with a maximum size of max_simultaneous_derivatives
, and we call implicit_tendency!
once for each partition. That is, if the partition size is $s$, then the first partition evaluates the coefficients of $ε₁$ through $εₛ$, the second evaluates the coefficients of $εₛ₊₁$ through $ε₂ₛ$, and so on until $εₙ$. The default partition size is 32.
ClimaAtmos.AutoSparseJacobian
— TypeAutoSparseJacobian(sparse_jacobian_alg, [padding_bands_per_block])
A JacobianAlgorithm
that computes the Jacobian using forward-mode automatic differentiation, assuming that the Jacobian's sparsity structure is given by sparse_jacobian_alg
.
Only entries that are exptected to be nonzero according to the sparsity structure are updated, but any other entries that are nonzero can introduce errors to the updated entries. This issue can be avoided by adding padding bands to blocks that are likely to introduce errors. In cases where the default padding bands are insufficient, padding_bands_per_block
can be specified to add a fixed number of padding bands to every block.
For more information about this algorithm, see Implicit Solver.
Internals
ClimaAtmos.parallel_lu_factorize!
— Functionparallel_lu_factorize!(device, matrices, ::Val{N})
Runs a parallel LU factorization algorithm on the specified device
. If each slice matrices[1:N, 1:N, i]
represents a matrix $Mᵢ$, this function overwrites it with the lower triangular matrix $Lᵢ$ and the upper triangular matrix $Uᵢ$, where $Mᵢ = Lᵢ * Uᵢ$. The value of N
must be wrapped in a Val
to ensure that it is statically inferrable, which allows the LU factorization to avoid dynamic local memory allocations.
The runtime of this algorithm scales as $O(N^3)$.
ClimaAtmos.parallel_lu_solve!
— Functionparallel_lu_solve!(device, vectors, matrices, ::Val{N})
Runs a parallel LU solver algorithm on the specified device
. If each slice vectors[1:N, i]
represents a vector $vᵢ$, and if each slice matrices[1:N, 1:N, i]
represents a matrix $Lᵢ * Uᵢ$ that was factorized by parallel_lu_factorize!
, this function overwrites the slice vectors[1:N, i]
with $(Lᵢ * Uᵢ)⁻¹ * vᵢ$. The value of N
must be wrapped in a Val
to ensure that it is statically inferrable, which allows the LU solver to avoid dynamic local memory allocations.
The runtime of this algorithm scales as $O(N^2)$.