API
Initial conditions
General
ClimaAtmos.InitialConditions.InitialCondition — Type
InitialConditionA 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 — Type
IsothermalProfile(; temperature = 300)An InitialCondition with a uniform temperature profile.
ClimaAtmos.InitialConditions.DecayingProfile — Type
DecayingProfile(; perturb = true)An InitialCondition with a decaying temperature profile, and with an optional perturbation to the temperature.
ClimaAtmos.InitialConditions.hydrostatic_pressure_profile — Function
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.
Plane / Box
ClimaAtmos.InitialConditions.ConstantBuoyancyFrequencyProfile — Type
ConstantBuoyancyFrequencyProfile()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 — Type
DryDensityCurrentProfile(; 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 — Type
RisingThermalBubbleProfile(; 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 — Type
DryBaroclinicWave(; perturb = true, deep_atmosphere = false)An InitialCondition with a dry baroclinic wave, and with an optional perturbation to the horizontal velocity.
ClimaAtmos.InitialConditions.MoistBaroclinicWaveWithEDMF — Type
MoistBaroclinicWaveWithEDMF(; 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 — Type
MoistAdiabaticProfileEDMFX(; perturb = true)An InitialCondition with a moist adiabatic temperature profile, and with an optional perturbation to the temperature.
Cases from literature
ClimaAtmos.InitialConditions.Nieuwstadt — Type
NieuwstadtThe InitialCondition described in [9], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.GABLS — Type
GABLSThe InitialCondition described in [10], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.GATE_III — Type
GATE_IIIThe InitialCondition described in [11], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.ARM_SGP — Type
ARM_SGPThe InitialCondition described in [12], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.DYCOMS_RF01 — Type
DYCOMS_RF01The InitialCondition described in [13], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.DYCOMS_RF02 — Type
DYCOMS_RF02The InitialCondition described in [14], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.Rico — Type
RicoThe InitialCondition described in [15], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.TRMM_LBA — Type
TRMM_LBAThe InitialCondition described in [16], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.LifeCycleTan2018 — Type
LifeCycleTan2018The InitialCondition described in [17], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.Bomex — Type
BomexThe InitialCondition described in [18], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.Soares — Type
SoaresThe InitialCondition described in [19], but with a hydrostatically balanced pressure profile.
ClimaAtmos.InitialConditions.RCEMIPIIProfile — Type
RCEMIPIIProfile(temperature, humidity)An InitialCondition following the sounding to initialize simulations for RCEMIPII as described by Wing et. al. (2018) (https://doi.org/10.5194/gmd-11-793-2018). There are three input profiles: RCEMIPIIProfile295, RCEMIPIIProfile300, and RCEMIPIIProfile305, that specify three different SST temperatures and different initial specific humidity profiles. Note: this should be used for RCEsmall and NOT RCElarge - RCElarge must be initialized with the final state of RCE_small.
Helper
ClimaAtmos.InitialConditions.ColumnInterpolatableField — Type
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)Grids
ClimaAtmos.ColumnGrid — Function
ColumnGrid(::Type{FT}; kwargs...)Create a ColumnGrid.
Arguments
FT: the floating-point type [Float32,Float64]
Keyword Arguments
context = ClimaComms.context(): the ClimaComms communications contextz_elem = 10: the number of z-pointsz_max = 30000.0: the domain maximum along the z-directionz_stretch = true: whether to use vertical stretchingdz_bottom = 500.0: bottom layer thickness for stretching
ClimaAtmos.SphereGrid — Function
SphereGrid(::Type{FT}; kwargs...)Create an ExtrudedCubedSphereGrid with topography support.
Arguments
FT: the floating-point type [Float32,Float64]
Keyword Arguments
context = ClimaComms.context(): the ClimaComms communications contextz_elem = 10: the number of z-pointsz_max = 30000.0: the domain maximum along the z-directionz_stretch = true: whether to use vertical stretchingdz_bottom = 500.0: bottom layer thickness for stretchingradius = 6.371229e6: the radius of the cubed sphereh_elem = 6: the number of horizontal elements per side of every panel (6 panels in total)nh_poly = 3: the polynomial order. Note: The number of quadrature points in 1D within each horizontal element is thenn_quad_points = nh_poly + 1bubble = false: enables the "bubble correction" for more accurate element areas when computing the spectral element spacedeep_atmosphere = true: use deep atmosphere equations and metric terms, otherwise assume columns are cylindrical (shallow atmosphere)topography = NoTopography(): topography typetopography_damping_factor = 5.0: factor by which smallest resolved length-scale is to be dampedmesh_warp_type = SLEVEWarp{FT}(): mesh warping type (SLEVEWarporLinearWarp)topo_smoothing = false: apply topography smoothing
ClimaAtmos.PlaneGrid — Function
PlaneGrid(::Type{FT}; kwargs...)Create a SliceXZGrid with topography support.
Arguments
FT: the floating-point type [Float32,Float64]
Keyword Arguments
context = ClimaComms.context(): the ClimaComms communications contextx_elem = 6: the number of x-pointsx_max = 300000.0: the domain maximum along the x-directionz_elem = 10: the number of z-pointsz_max = 30000.0: the domain maximum along the z-directionnh_poly = 3: the polynomial order. Note: The number of quadrature points in 1D within each horizontal element is thenn_quad_points = nh_poly + 1z_stretch = true: whether to use vertical stretchingdz_bottom = 500.0: bottom layer thickness for stretchingbubble = false: enables the "bubble correction" for more accurate element areas when computing the spectral element space. Note: Currently not supported by SliceXZGrid in ClimaCore.periodic_x = true: use periodic domain along x-directiontopography = NoTopography(): topography typetopography_damping_factor = 5.0: factor by which smallest resolved length-scale is to be dampedmesh_warp_type = LinearWarp(): mesh warping type (SLEVEWarporLinearWarp)topo_smoothing = false: apply topography smoothing
ClimaAtmos.BoxGrid — Function
BoxGrid(::Type{FT}; kwargs...)Create a Box3DGrid with topography support.
Arguments
FT: the floating-point type [Float32,Float64]
Keyword Arguments
context = ClimaComms.context(): the ClimaComms communications contextx_elem = 6: the number of x-pointsx_max = 300000.0: the domain maximum along the x-directiony_elem = 6: the number of y-pointsy_max = 300000.0: the domain maximum along the y-directionz_elem = 10: the number of z-pointsz_max = 30000.0: the domain maximum along the z-directionnh_poly = 3: the polynomial order. Note: The number of quadrature points in 1D within each horizontal element is thenn_quad_points = nh_poly + 1z_stretch = true: whether to use vertical stretchingdz_bottom = 500.0: bottom layer thickness for stretchingbubble = false: enables the "bubble correction" for more accurate element areas when computing the spectral element space.periodic_x = true: use periodic domain along x-directionperiodic_y = true: use periodic domain along y-directiontopography = NoTopography(): topography typetopography_damping_factor = 5.0: factor by which smallest resolved length-scale is to be dampedmesh_warp_type = LinearWarp(): mesh warping type (SLEVEWarporLinearWarp)topo_smoothing = false: apply topography smoothing
Jacobian
ClimaAtmos.Jacobian — Type
Jacobian(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 — Type
JacobianAlgorithmA 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 — Type
ManualSparseJacobian(
topography_flag,
diffusion_flag,
sgs_advection_flag,
sgs_entr_detr_flag,
sgs_mass_flux_flag,
sgs_nh_pressure_flag,
sgs_vertdiff_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 DerivativeFlags 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 computedsgs_vertdiff_flag::DerivativeFlag: whether the derivatives of the subgrid-scale vertical diffusion tendency should be computedapproximate_solve_iters::Int: number of iterations to take for the approximate linear solve required when thediffusion_flagisUseDerivative
ClimaAtmos.AutoDenseJacobian — Type
AutoDenseJacobian([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 — Type
AutoSparseJacobian(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.
Topography
ClimaAtmos.CosineTopography — Type
CosineTopography{D, FT}(; h_max = 25, λ = 25e3)Cosine hill topography in 2D or 3D.
Arguments
h_max::FT: Maximum elevation (m)λ::FT: Wavelength of the cosine hills (m)
ClimaAtmos.AgnesiTopography — Type
AgnesiTopography{FT}(; h_max = 25, x_center = 50e3, a = 5e3)Witch of Agnesi mountain topography for 2D simulations.
Arguments
h_max: Maximum elevation (m)x_center: Center position (m)a: Mountain width parameter (m)
ClimaAtmos.ScharTopography — Type
ScharTopography{FT}(; h_max = 25, x_center = 50e3, λ = 4e3, a = 5e3)Schar mountain topography for 2D simulations.
Arguments
h_max: Maximum elevation (m)x_center: Center position (m)λ: Wavelength parameter (m)a: Mountain width parameter (m)
ClimaAtmos.EarthTopography — Type
EarthTopography()Earth topography from ETOPO2022 data files.
ClimaAtmos.DCMIP200Topography — Type
DCMIP200Topography()Surface elevation for the DCMIP-2-0-0 test problem.
ClimaAtmos.Hughes2023Topography — Type
Hughes2023Topography()Surface elevation for baroclinic wave test from Hughes and Jablonowski (2023).
ClimaAtmos.SLEVEWarp — Type
SLEVEWarp(; eta = 0.7, s = 10.0)Smooth Level Vertical (SLEVE) coordinate warping for terrain-following meshes.
Arguments
eta: Threshold parameter (if z/z_top > eta, no warping is applied). Default: 0.7s: Decay scale parameter controlling how quickly the warping decays with height. Default: 10.0
References
Schär et al. (2002), "A new terrain-following vertical coordinate formulation for atmospheric prediction models", Mon. Wea. Rev.
ClimaAtmos.LinearWarp — Type
LinearWarp()Linear mesh warping that uniformly distributes vertical levels between the surface and top of the domain.
Internals
ClimaAtmos.parallel_lu_factorize! — Function
parallel_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! — Function
parallel_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)$.