API

Simulation

ClimaAtmos.AtmosSimulationType
AtmosSimulation(config::AtmosConfig)

Construct a simulation from a YAML-based configuration. Construct an atmospheric simulation with the default floating-point type Float32. Equivalent to AtmosSimulation{Float32}(; kwargs...).

source
AtmosSimulation(; kwargs...)

Construct an atmospheric simulation with the default floating-point type Float32. Equivalent to AtmosSimulation{Float32}(; kwargs...).

source
AtmosSimulation{FT}(; kwargs...) where {FT}

Construct an atmospheric simulation with floating-point type FT (default: Float32).

Keyword Arguments

Model and domain

  • model::AtmosModel = AtmosModel(): Physics and parameterization configuration.
  • params::ClimaAtmosParameters = ClimaAtmosParameters(FT): Physical parameters.
  • grid::AbstractGrid = SphereGrid(FT; ...): Computational grid. Use ColumnGrid, BoxGrid, PlaneGrid, or SphereGrid.
  • setup = Setups.DecayingProfile(; perturb=true, params): Setup defining the initial state. See Setups for available options.

Time

  • dt = 600: Timestep in seconds.
  • t_start = 0: Start time in seconds.
  • t_end = 864000: End time in seconds (default: 10 days).
  • start_date = DateTime(2010, 1, 1): Calendar reference date.

Output

  • job_id::String = "atmos_sim": Run identifier, used in output directory naming.
  • output_dir = nothing: Output directory path. Auto-generated from job_id if nothing.
  • output_dir_style = "activelink": Output directory organization style.
  • checkpoint_frequency = Inf: How often to save restart checkpoints (seconds).
  • log_to_file::Bool = false: Write log output to a file in output_dir.

Diagnostics

  • diagnostics::DiagnosticsConfig = DiagnosticsConfig(): Specification of which diagnostics the simulation produces and how their NetCDF output is shaped. See DiagnosticsConfig.

Callbacks

  • default_callbacks::Bool = true: Enable common simulation callbacks.
  • callbacks = (): Additional user-provided callbacks.
  • callback_kwargs = (): Extra keyword arguments forwarded to default callbacks.

Restarts

  • restart_file = nothing: Path to a restart file to resume from.
  • detect_restart_file::Bool = false: Automatically detect the latest restart file in a structured output directory.

Numerics

  • ode_config: ODE solver algorithm. Default: IMEXAlgorithm(ARS343(), NewtonsMethod(...)).
  • jacobian::JacobianAlgorithm = ManualSparseJacobian(; approximate_solve_iters = 1): Jacobian algorithm for the implicit solve. Use ManualSparseJacobian, AutoSparseJacobian, or AutoDenseJacobian.
  • debug_jacobian::Bool = false: Enable Jacobian debugging output.
  • tracers = []: Additional tracer species.

Example

import ClimaAtmos as CA

# Minimal: 1-day global simulation with defaults
simulation = CA.AtmosSimulation{Float64}(; t_end = 86400)
CA.solve_atmos!(simulation)

# Single-column BOMEX case
simulation = CA.AtmosSimulation{Float64}(;
    grid = CA.ColumnGrid(Float64; z_elem = 60, z_max = 3000.0),
    setup = CA.Setups.Bomex(),
    dt = 5,
    t_end = 3600 * 6,
)
source

Presets

ClimaAtmos.Presets.baroclinic_waveFunction
baroclinic_wave([FT = Float32]; kwargs...)

Dry baroclinic-wave simulation preset: global SphereGrid, DryBaroclinicWave setup, and a dry model with disable_surface_flux_tendency = true. For the moist variant, pass setup = Setups.MoistBaroclinicWave() and model = Presets.equil_moist_0m(; disable_surface_flux_tendency = true). Keyword arguments are forwarded to AtmosSimulation.

source
ClimaAtmos.Presets.bomexFunction
bomex([FT = Float32]; kwargs...)

BOMEX shallow-cumulus single-column simulation preset: ColumnGrid (60 uniform levels, zmax = 3 km), Setups.Bomex setup, [`equilmoist0m](@ref) physics,dt = 10 s,tend = 6 h`.

No EDMF turbulence-convection scheme is enabled by default; pass model = Presets.diagnostic_edmf(FT) or Presets.prognostic_edmf(FT) to add one. Keyword arguments are forwarded to AtmosSimulation.

source
ClimaAtmos.Presets.dryFunction
dry(; kwargs...)

Dry atmosphere preset (microphysics_model = DryModel()). Keyword arguments are forwarded to AtmosModel.

source
ClimaAtmos.Presets.equil_moist_0mFunction
equil_moist_0m(; kwargs...)

Equilibrium-moisture preset with 0-moment microphysics, grid-scale cloud, prescribed zonally-symmetric SST, and idealized insolation. Keyword arguments are forwarded to AtmosModel.

source
ClimaAtmos.Presets.nonequil_moist_1mFunction
nonequil_moist_1m(; kwargs...)

Non-equilibrium-moisture preset with 1-moment microphysics, explicit microphysics tendency timestepping, grid-scale cloud, prescribed zonally-symmetric SST, and idealized insolation. Mirrors equil_moist_0m but with 1-moment non-equilibrium microphysics in place of 0-moment equilibrium. Keyword arguments are forwarded to AtmosModel.

source
ClimaAtmos.Presets.diagnostic_edmfFunction
diagnostic_edmf([FT = Float32]; area_fraction, n_updrafts, prognostic_tke, kwargs...)

Equilibrium-moist model with the DiagnosticEDMFX turbulence-convection scheme wired in with Generalized entrainment/detrainment, SGS mass & diffusive fluxes, and non-hydrostatic pressure drag (matches the canonical diagnostic_edmfx_* configs).

area_fraction defaults to 1e-5. Override microphysics_model to pair EDMF with a dry or non-equilibrium scheme. All remaining keyword arguments are forwarded to AtmosModel.

source
ClimaAtmos.Presets.prognostic_edmfFunction
prognostic_edmf([FT = Float32]; area_fraction, n_updrafts, prognostic_tke, kwargs...)

Equilibrium-moist model with the PrognosticEDMFX turbulence-convection scheme. In addition to the diagnostic_edmf EDMF settings, this also enables prognostic updraft vertical diffusion and the relaxation filter on negative updraft velocities (matches the canonical prognostic_edmfx_* configs).

All remaining keyword arguments are forwarded to AtmosModel.

source
ClimaAtmos.Presets.prognostic_edmf_1mFunction
prognostic_edmf_1m([FT = Float32]; kwargs...)

prognostic_edmf with 1-moment non-equilibrium microphysics and explicit microphysics tendency timestepping (matches the canonical prognostic_edmfx_* configs that use microphysics_model: "1M"). All keyword arguments are forwarded to prognostic_edmf and on to AtmosModel.

source

Grids

ClimaAtmos.ColumnGridFunction
ColumnGrid(::Type{FT}; kwargs...)

Create a ColumnGrid.

Arguments

  • FT: the floating-point type [Float32, Float64]

Keyword Arguments

  • context = ClimaComms.context(): the ClimaComms communications context
  • z_elem = 10: the number of z-points
  • z_max = 30000.0: the domain maximum along the z-direction
  • z_stretch = true: whether to use vertical stretching
  • dz_bottom = 500.0: bottom layer thickness for stretching
  • z_mesh: Optionally provide a custom z-mesh, instead of z_elem, z_max, z_stretch
source
ClimaAtmos.SphereGridFunction
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 context
  • z_elem = 10: the number of z-points
  • z_max = 30000.0: the domain maximum along the z-direction
  • z_stretch = true: whether to use vertical stretching
  • dz_bottom = 500.0: bottom layer thickness for stretching
  • z_mesh: Optionally provide a custom z-mesh, instead of z_elem, z_max, z_stretch
  • radius = 6.371229e6: the radius of the cubed sphere
  • h_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 then n_quad_points = nh_poly + 1
  • bubble = false: enables the "bubble correction" for more accurate element areas when computing the spectral element space
  • deep_atmosphere = true: use deep atmosphere equations and metric terms, otherwise assume columns are cylindrical (shallow atmosphere)
  • topography = NoTopography(): topography type
  • topography_damping_factor = 5.0: factor by which smallest resolved length-scale is to be damped
  • mesh_warp_type = SLEVEWarp{FT}(): mesh warping type (SLEVEWarp or LinearWarp)
  • topo_smoothing = false: apply topography smoothing
source
ClimaAtmos.PlaneGridFunction
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 context
  • x_elem = 6: the number of x-points
  • x_max = 300000.0: the domain maximum along the x-direction
  • z_elem = 10: the number of z-points
  • z_max = 30000.0: the domain maximum along the z-direction
  • z_mesh = DefaultZMesh(FT; z_min = 0, z_max, z_elem, z_stretch): the vertical mesh
  • nh_poly = 3: the polynomial order. Note: The number of quadrature points in 1D within each horizontal element is then n_quad_points = nh_poly + 1
  • z_stretch = true: whether to use vertical stretching
  • dz_bottom = 500.0: bottom layer thickness for stretching
  • periodic_x = true: use periodic domain along x-direction
  • topography = NoTopography(): topography type
  • topography_damping_factor = 5.0: factor by which smallest resolved length-scale is to be damped
  • mesh_warp_type = LinearWarp(): mesh warping type (SLEVEWarp or LinearWarp)
  • topo_smoothing = false: apply topography smoothing
source
ClimaAtmos.BoxGridFunction
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 context
  • x_elem = 6: the number of x-points
  • x_max = 300000.0: the domain maximum along the x-direction
  • y_elem = 6: the number of y-points
  • y_max = 300000.0: the domain maximum along the y-direction
  • z_elem = 10: the number of z-points
  • z_max = 30000.0: the domain maximum along the z-direction
  • nh_poly = 3: the polynomial order. Note: The number of quadrature points in 1D within each horizontal element is then n_quad_points = nh_poly + 1
  • z_stretch = true: whether to use vertical stretching
  • dz_bottom = 500.0: bottom layer thickness for vertical stretching
  • z_mesh: Optionally provide a custom z-mesh, instead of z_elem, z_max, z_stretch
  • bubble = false: enables the "bubble correction" for more accurate element areas when computing the spectral element space.
  • periodic_x = true: use periodic domain along x-direction
  • periodic_y = true: use periodic domain along y-direction
  • topography = NoTopography(): topography type
  • topography_damping_factor = 5.0: factor by which smallest resolved length-scale is to be damped
  • mesh_warp_type = LinearWarp(): mesh warping type (SLEVEWarp or LinearWarp)
  • topo_smoothing = false: apply topography smoothing
source

Jacobian

ClimaAtmos.JacobianType
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.

source
ClimaAtmos.JacobianAlgorithmType
JacobianAlgorithm

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.

source
ClimaAtmos.ManualSparseJacobianType
ManualSparseJacobian(; approximate_solve_iters = 1)

A JacobianAlgorithm that approximates the Jacobian using analytically derived tendency derivatives and inverts it using a specialized nested linear solver.

Which derivative blocks are computed is determined automatically from the AtmosModel (topography, diffusion mode, EDMF modes) when the cache is built — users do not configure them directly.

Arguments

  • approximate_solve_iters::Int = 1: number of iterations to take for the approximate linear solve required when grid-scale diffusion is treated implicitly.
source
ClimaAtmos.AutoDenseJacobianType
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.

source
ClimaAtmos.AutoSparseJacobianType
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.

source

Diagnostics

ClimaAtmos.Diagnostics.DiagnosticsConfigType
DiagnosticsConfig(; default = true, additional = (), interpolation_num_points = nothing, output_at_levels = true)

Specification of which diagnostics a simulation produces and how their NetCDF output is shaped. A single DiagnosticsConfig value is passed to AtmosSimulation via the diagnostics keyword argument.

Fields

  • default::Bool = true: include the built-in ClimaAtmos diagnostic set for the chosen AtmosModel.

  • additional = (): extra user-supplied diagnostics. Each entry can be:

    • a ClimaDiagnostics.ScheduledDiagnostic (full control);
    • a Pair like "ua" => (; period = "30mins", reduction = "average") (short_name => options as a NamedTuple);
    • a NamedTuple with at least short_name and period, e.g. (; short_name = "ts", period = "1hours");
    • a YAML-style Dict{String,Any} (the same shape produced by the diagnostics: YAML key).

    Mixed lists are allowed.

  • interpolation_num_points = nothing: override the NetCDF remap grid (e.g. (180, 90, 10)). When nothing, falls back to the default chosen from the underlying space.

  • output_at_levels::Bool = true: write at model levels (no vertical interpolation). Set false to interpolate to pressure levels.

A simulation produces no diagnostics when default = false and additional is empty.

source

Topography

ClimaAtmos.CosineTopographyType
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)
source
ClimaAtmos.AgnesiTopographyType
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)
source
ClimaAtmos.ScharTopographyType
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)
source
ClimaAtmos.SLEVEWarpType
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.7
  • s: 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.

source
ClimaAtmos.LinearWarpType
LinearWarp()

Linear mesh warping that uniformly distributes vertical levels between the surface and top of the domain.

source

Surface

See the Surface Conditions page for a guide to these types and how to choose among them.

ClimaAtmos.AtmosSurfaceType
AtmosSurface

Groups surface-related models and types.

source

Flux schemes

ClimaAtmos.SurfaceConditions.MoninObukhovType
MoninObukhov(; z0, z0m, z0b, fluxes, shf, lhf, θ_flux, q_flux, ustar)

Container for storing values used to calculate surface conditions using Monin-Obukhov Similarity Theory. See the SurfaceFluxes.jl MOST documentation for more information.

Roughness (required)

  • z0: Roughness (sets both z0m and z0b)
  • z0m, z0b: Roughness for momentum and scalars (specify both, or use z0)

Prescribed fluxes (optional) — specify via one of:

  • fluxes: A HeatFluxes/θAndQFluxes struct, or a callable (t, FT) -> HeatFluxes/θAndQFluxes for time-varying fluxes (resolved once per surface update by resolve_flux_scheme, before the per-cell broadcast)
  • shf, lhf: Sensible/latent heat fluxes (W/m²) — constructs HeatFluxes
  • θ_flux, q_flux: θ and q fluxes (K·m/s, kg/kg·m/s) — constructs θAndQFluxes

Other (optional)

  • ustar: Friction velocity (m/s)

Valid combinations:

  • roughness
  • roughness and fluxes
  • roughness and ustar
  • roughness and fluxes and ustar
source
ClimaAtmos.SurfaceConditions.ExchangeCoefficientsType
ExchangeCoefficients(; Cd, Ch)
ExchangeCoefficients(C)

Bulk-aerodynamic surface flux closure with fixed, dimensionless exchange coefficients — a SurfaceParameterization alternative to MoninObukhov in which the turbulent fluxes scale linearly with the near-surface wind speed and the air–surface differences (rather than being derived from Monin–Obukhov stability).

  • Cd: momentum (drag) exchange coefficient.
  • Ch: thermal/scalar (heat and moisture) exchange coefficient.

The single-argument form ExchangeCoefficients(C) sets Cd = Ch = C.

source
ClimaAtmos.SurfaceConditions.HeatFluxesType
HeatFluxes(; shf, lhf = nothing)

Prescribed surface turbulent energy fluxes, used as the fluxes field of a MoninObukhov closure. Both use the sign convention that positive is upward (directed from the surface into the atmosphere).

  • shf: sensible heat flux (W/m²).
  • lhf: latent heat flux (W/m²). Optional — nothing is treated as zero, and lhf must be left unset for a DryModel (specifying it with a dry model is an error).
source
ClimaAtmos.SurfaceConditions.θAndQFluxesType
θAndQFluxes(; θ_flux, q_flux = nothing)

Prescribed surface kinematic fluxes of potential temperature and total specific humidity, used as the fluxes field of a MoninObukhov closure. They are converted per surface point into the sensible/latent heat fluxes actually applied, via shf = θ_flux * ρ_sfc * cp_m and lhf = q_flux * ρ_sfc * Lᵥ. Positive is upward (surface into atmosphere).

  • θ_flux: potential-temperature flux (K·m/s).
  • q_flux: total-specific-humidity flux (kg/kg·m/s). Optional — nothing is treated as zero, and q_flux must be left unset for a DryModel.
source

Surface temperature

ClimaAtmos.SurfaceConditions.SurfaceTemperatureType
SurfaceTemperature

Abstract supertype for ways of obtaining the surface temperature T_sfc used when computing surface conditions. Concrete subtypes:

surface_temperature(t, Y, p, t_time) produces the value(s) that update_surface_conditions! will broadcast across the surface.

source
ClimaAtmos.SurfaceConditions.AnalyticTemperatureType
AnalyticTemperature(f)

A surface temperature given by f(coordinates, params, t). Used for the analytic SST formulas (zonally-symmetric, RCEMIPII), time-varying setups (e.g. GABLS), and spatially-uniform constants (AnalyticTemperature(Returns(T))). f is broadcast per-coordinate inside the surface update; if the formula does not depend on time, simply ignore the t argument.

source
ClimaAtmos.SurfaceConditions.ExternalTemperatureType
ExternalTemperature()

A surface temperature read from a time-varying external input. surface_temperature(::ExternalTemperature, Y, p, t) evaluates the field from p.external_forcing.surface_inputs.ts (driven by surface_timevaryinginputs.ts), so this temperature is only valid when the setup populates external_forcing.surface_inputs (e.g. InterpolatedColumnProfile).

source

Boundary overrides

ClimaAtmos.SurfaceConditions.SurfaceBoundaryOverridesType
SurfaceBoundaryOverrides(; p, q_vap, u, v, gustiness, beta)

Per-point overrides for surface boundary values used by SurfaceFluxes. Fields default to nothing, in which case sensible defaults are used:

  • p: surface pressure (default: hydrostatic from interior)
  • q_vap: surface specific humidity (default: q_vap_sat at T_sfc)
  • u, v: surface horizontal winds (default: 0)
  • gustiness: turbulent gustiness (default: 1)
  • beta: moisture availability (default: 1)

For the coupler use case, a Fields.Field{<:SurfaceBoundaryOverrides} may be stored on the cache so that the coupler can override per-cell values.

source

Surface albedo

ClimaAtmos.ConstantAlbedoType
struct ConstantAlbedo{FT} <: SurfaceAlbedoModel

A constant surface albedo model. The default value is α = 0.38. It is used purely for idealized experiments.

source
ClimaAtmos.CouplerAlbedoType
CouplerAlbedo()

Surface albedo supplied by an external driver (the coupler), which writes the direct/diffuse shortwave albedos into the radiation cache. ClimaAtmos performs no albedo computation of its own in this mode.

source

Core functions

ClimaAtmos.SurfaceConditions.update_surface_conditions!Function
update_surface_conditions!(Y, p, t)

Updates p.precomputed.sfc_conditions based on the current state Y and time t. Skips work if the surface model has no flux parameterization (isnothing(atmos.surface.flux_scheme)), which is the coupler-handoff case.

source
ClimaAtmos.SurfaceConditions.surface_state_to_conditionsFunction
surface_state_to_conditions(
    overrides, flux_scheme, T_sfc_in,
    surface_local_geometry,
    T_int, ρ_int, q_tot_int, q_liq_int, q_ice_int, u_int, v_int, z_int,
    thermo_params, surface_fluxes_params, surface_temp_params,
    atmos,
)

Compute the surface conditions at one point. T_sfc_in is either a scalar, the resolved temperature field value, or an AnalyticTemperature to evaluate against the local coordinates.

source
ClimaAtmos.SurfaceConditions.atmos_surface_conditionsFunction
atmos_surface_conditions(surface_conditions, ρ_sfc, surface_local_geometry)

Adds local geometry information to the SurfaceFluxes.SurfaceFluxConditions struct. The resulting values are the ones actually used by ClimaAtmos operator boundary conditions.

source

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)$.

source
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)$.

source