AtmosModel
AtmosProblem
ClimateMachine.Atmos.AtmosProblem
— TypeAtmosProblem
The default problem definition (initial and boundary conditions) for AtmosModel
.
AtmosModel balance law
ClimateMachine.Atmos.AtmosPhysics
— TypeAtmosPhysics
An AtmosPhysics
for atmospheric physics
Usage
AtmosPhysics(
param_set,
ref_state,
energy,
moisture,
compressibility,
turbulence,
turbconv,
hyperdiffusion,
precipitation,
radiation,
tracers,
lsforcing,
)
Fields
param_set
Parameter Set (type to dispatch on, e.g., planet parameters. See CLIMAParameters.jl package)
ref_state
Reference State (For initial conditions, or for linearisation when using implicit solvers)
energy
Energy sub-model, can be energy-based or θliqice-based
moisture
Moisture Model (Equations for dynamics of moist variables)
compressibility
Compressibility switch
turbulence
Turbulence Closure (Equations for dynamics of under-resolved turbulent flows)
turbconv
Turbulence Convection Closure (e.g., EDMF)
hyperdiffusion
Hyperdiffusion Model (Equations for dynamics of high-order spatial wave attenuation)
viscoussponge
Viscous sponge layers
precipitation
Precipitation Model (Equations for dynamics of precipitating species)
radiation
Radiation Model (Equations for radiative fluxes)
tracers
Tracer Terms (Equations for dynamics of active and passive tracers)
lsforcing
Large-scale forcing (Forcing information from GCMs, reanalyses, or observations)
ClimateMachine.Atmos.AtmosModel
— TypeAtmosModel <: BalanceLaw
A BalanceLaw
for atmosphere modeling. Users may over-ride prescribed default values for each field.
Usage
AtmosModel(
physics,
problem,
orientation,
source,
data_config,
)
Fields
physics
Atmospheric physics
problem
Problem (initial and boundary conditions)
orientation
An orientation model
source
Source Terms (Problem specific source terms)
data_config
Data Configuration (Helper field for experiment configuration)
AtmosModel methods
ClimateMachine.BalanceLaws.flux_first_order!
— Methodflux_first_order!(
bl::BalanceLaw,
flux::Grad,
state::Vars,
aux::Vars,
t::Real
)
Computes (and assembles) flux terms F¹(Y)
in:
∂Y
-- + ∇ • F¹(Y) + ∇ • F²(Y,G) = S(Y, G), G = ∇Y
∂t
Computes and assembles non-diffusive fluxes in the model equations.
For this fallback to work, several methods must be defined:
optionally,
and individual flux
kernels that are defined for each type that eq_tends
returns.
ClimateMachine.BalanceLaws.flux_second_order!
— Methodflux_second_order!(
bl::BalanceLaw,
flux::Grad,
state::Vars,
diffusive::Vars,
hyperdiffusive::Vars,
aux::Vars,
t::Real
)
Computes (and assembles) flux terms F²(Y, G)
in:
∂Y
-- + ∇ • F¹(Y) + ∇ • F²(Y,G) = S(Y, G), G = ∇Y
∂t
Diffusive fluxes in BalanceLaw. Viscosity, diffusivity are calculated in the turbulence subcomponent and accessed within the diffusive flux function. Contributions from subcomponents are then assembled (pointwise).
For this fallback to work, several methods must be defined:
optionally,
and individual flux
kernels that are defined for each type that eq_tends
returns.
ClimateMachine.BalanceLaws.init_state_auxiliary!
— Methodinit_state_auxiliary!(
m::AtmosModel,
aux::Vars,
grid,
direction
)
Initialise auxiliary variables for each AtmosModel subcomponent. Store Cartesian coordinate information in aux.coord
.
ClimateMachine.BalanceLaws.source!
— Methodsource!(
bl::BalanceLaw,
source::Vars,
state::Vars,
diffusive::Vars,
aux::Vars,
t::Real,
direction::Direction,
)
Computes (and assembles) source terms S(Y)
in:
∂Y
-- + ∇ • F¹(Y) + ∇ • F²(Y,G) = S(Y, G), G = ∇Y
∂t
For this fallback to work, several methods must be defined:
optionally,
and individual source
kernels that are defined for each type that eq_tends
returns.
ClimateMachine.BalanceLaws.init_state_prognostic!
— Methodinit_state_prognostic!(
m::AtmosModel,
state::Vars,
aux::Vars,
localgeo,
t,
args...,
)
Initialise state variables. args...
provides an option to include configuration data (current use cases include problem constants, spline-interpolants).
Compressibility
ClimateMachine.Atmos.Compressible
— TypeCompressible <: Compressibilty
Dispatch on compressible model (default)
- Density is prognostic
ClimateMachine.Atmos.Anelastic1D
— TypeAnelastic1D <: Compressibilty
Dispatch on Anelastic1D model
- The state density is taken constant in time and equal to the reference density. This constant density profile is used in all equations and conversions from conservative to specific variables per unit mass. The density can be accessed using the dispatch function
density(atmos, state, aux)
. - The thermodynamic state is constructed from the reference pressure (constant in time), and the internal energy (which evolves in time).
- The state density is not consistent with the thermodynamic state, since we neglect buoyancy perturbations on all equations except in the vertical buoyancy flux.
- The density obtained from the thermodynamic state,
air_density(ts)
, recovers the full density, which should only be used to compute buoyancy and buoyancy fluxes, and in the FV reconstruction. - Removes momentum z-component tendencies, assuming balance between the pressure gradient and buoyancy forces.
Reference states
ClimateMachine.Atmos.HydrostaticState
— TypeHydrostaticState{P,T} <: ReferenceState
A hydrostatic state specified by a virtual temperature profile and relative humidity.
By default, this is a dry hydrostatic reference state.
ClimateMachine.Atmos.InitStateBC
— TypeInitStateBC
Set the value at the boundary to match the init_state_prognostic!
function. This is mainly useful for cases where the problem has an explicit solution.
TODO: This should be fixed later once BCs are figured out (likely want
different things here?)
ClimateMachine.Atmos.ReferenceState
— TypeReferenceState
Hydrostatic reference state, for example, used as initial condition or for linearization.
ClimateMachine.Atmos.NoReferenceState
— TypeNoReferenceState <: ReferenceState
No reference state used
Thermodynamics
ClimateMachine.Atmos.recover_thermo_state
— Functionrecover_thermo_state(atmos::AtmosModel, state::Vars, aux::Vars)
An atmospheric thermodynamic state.
For now, we are directly calling newthermostate to avoid inconsistent aux states in kernels where the aux states are out of sync with the boundary state.
TODO: Define/call recover_thermo_state
when it's safely implemented
(see https://github.com/CliMA/ClimateMachine.jl/issues/1648)
ClimateMachine.Atmos.new_thermo_state
— Functionnew_thermo_state(atmos::AtmosModel, state::Vars, aux::Vars)
Create a new thermodynamic state, based on the state
, and not the aux
state.
This method calls the iterative saturation adjustment procedure for EquilMoist models.
ClimateMachine.Atmos.recover_thermo_state_anelastic
— Functionrecover_thermo_state_anelastic(atmos::AtmosModel, state::Vars, aux::Vars)
An atmospheric thermodynamic state.
For now, we are directly calling newthermostate_anelastic to avoid inconsistent aux states in kernels where the aux states are out of sync with the boundary state.
TODO: Define/call recover_thermo_state_anelastic
when it's safely implemented
(see https://github.com/CliMA/ClimateMachine.jl/issues/1648)
ClimateMachine.Atmos.new_thermo_state_anelastic
— Functionnew_thermo_state_anelastic(atmos::AtmosModel, state::Vars, aux::Vars)
Create a new thermodynamic state, based on the state
, and not the aux
state.
This method calls the iterative saturation adjustment procedure for EquilMoist models.
Moisture and Precipitation
ClimateMachine.Atmos.DryModel
— TypeDryModel
Assumes the moisture components is in the dry limit.
ClimateMachine.Atmos.EquilMoist
— TypeEquilMoist
Assumes the moisture components are computed via thermodynamic equilibrium.
ClimateMachine.Atmos.NonEquilMoist
— TypeNonEquilMoist
Does not assume that the moisture components are in equilibrium.
ClimateMachine.Atmos.NoPrecipitation
— TypeNoPrecipitation <: PrecipitationModel
No precipitation.
ClimateMachine.Atmos.RainModel
— TypeRainModel <: PrecipitationModel
Precipitation model with rain.
ClimateMachine.Atmos.RainSnowModel
— TypeRainSnowModel <: PrecipitationModel
Precipitation model with rain and snow.
Stabilization
ClimateMachine.Atmos.RayleighSponge
— TypeRayleighSponge{FT} <: TendencyDef{Source}
Rayleigh Damping (Linear Relaxation) for top wall momentum components Assumes laterally periodic boundary conditions for LES flows. Momentum components are relaxed to reference values (zero velocities) at the top boundary.
BCs
ClimateMachine.Atmos.AtmosBC
— TypeAtmosBC(momentum = Impenetrable(FreeSlip())
energy = Insulating()
moisture = Impermeable()
precipitation = OutflowPrecipitation()
tracer = ImpermeableTracer())
The standard boundary condition for AtmosModel
. The default options imply a "no flux" boundary condition.
ClimateMachine.Atmos.DragLaw
— TypeDragLaw(fn) :: MomentumDragBC
Drag law for momentum parallel to the boundary. The drag coefficient is C = fn(state, aux, t, normu_int_tan)
, where normu_int_tan
is the internal speed parallel to the boundary. _int
refers to the first interior node.
ClimateMachine.Atmos.Impermeable
— TypeImpermeable() :: MoistureBC
No moisture flux.
ClimateMachine.Atmos.PrescribedMoistureFlux
— TypePrescribedMoistureFlux(fn) :: MoistureBC
Prescribe the net inward moisture flux across the boundary by fn
, a function with signature fn(state, aux, t)
, returning the flux (in kg/m^2).
ClimateMachine.Atmos.BulkFormulaMoisture
— TypeBulkFormulaMoisture(fn) :: MoistureBC
Calculate the net inward moisture flux across the boundary using the bulk formula. The drag coefficient is C_q = fn_C_q(state, aux, t, normu_int_tan)
. The surface qtot at the boundary is `qtot = fnqtot(state, aux, t)`.
Return the flux (in kg m^-2 s^-1).
ClimateMachine.Atmos.FreeSlip
— TypeFreeSlip() :: MomentumDragBC
No surface drag on momentum parallel to the boundary.
ClimateMachine.Atmos.PrescribedTemperature
— TypePrescribedTemperature(fn) :: EnergyBC
Prescribe the temperature at the boundary by fn
, a function with signature fn(state, aux, t)
returning the temperature (in K).
ClimateMachine.Atmos.PrescribedEnergyFlux
— TypePrescribedEnergyFlux(fn) :: EnergyBC
Prescribe the net inward energy flux across the boundary by fn
, a function with signature fn(state, aux, t)
, returning the flux (in W/m^2).
ClimateMachine.Atmos.NishizawaEnergyFlux
— TypeNishizawaEnergyFlux(fn) :: EnergyBC
Calculate the net inward energy flux across the boundary following Nishizawa and Kitamura (2018). Return the flux (in W m^-2).
ClimateMachine.Atmos.Adiabaticθ
— TypeAdiabaticθ(fn) :: EnergyBC
Prescribe the net inward potential temperature flux across the boundary by fn
, a function with signature fn(state, aux, t)
, returning the flux (in kgK/m^2).
ClimateMachine.Atmos.BulkFormulaEnergy
— TypeBulkFormulaEnergy(fn) :: EnergyBC
Calculate the net inward energy flux across the boundary. The drag coefficient is C_h = fn_C_h(atmos, state, aux, t, normu_int_tan)
. The surface temperature and qtot are `T, qtot = fnTandqtot(atmos, state, aux, t)`. Return the flux (in W m^-2).
ClimateMachine.Atmos.OutflowPrecipitation
— TypeOutflowPrecipitation() :: PrecipitationBC
Free flux out of the domain.
ClimateMachine.Atmos.ImpermeableTracer
— TypeImpermeableTracer() :: TracerBC
No tracer diffusion across boundary
ClimateMachine.Atmos.Impenetrable
— TypeImpenetrable(drag::MomentumDragBC) :: MomentumBC
Defines an impenetrable wall model for momentum. This implies:
- no flow in the direction normal to the boundary, and
- flow parallel to the boundary is subject to the
drag
condition.
ClimateMachine.Atmos.Insulating
— TypeInsulating() :: EnergyBC
No energy flux across the boundary.
ClimateMachine.Atmos.NoSlip
— TypeNoSlip() :: MomentumDragBC
Zero momentum at the boundary.
ClimateMachine.Atmos.average_density
— Functionaverage_density(ρ_sfc, ρ_int)
Average density between the surface and the interior point, given
ρ_sfc
density at the surfaceρ_int
density at the interior point
Sources
ClimateMachine.Atmos.RemovePrecipitation
— TypeRemovePrecipitation <: TendencyDef{Source}
A sink to q_tot
when cloud condensate is exceeding a threshold. The threshold is defined either in terms of condensate or supersaturation. The removal rate is implemented as a relaxation term in the CloudMicrophysics.jl Microphysics_0M module. The default thresholds and timescale are defined in CLIMAParameters.jl.
ClimateMachine.Atmos.CreateClouds
— TypeCreateClouds <: TendencyDef{Source}
A source/sink to q_liq
and q_ice
implemented as a relaxation towards equilibrium in the Microphysics module. The default relaxation timescales are defined in CLIMAParameters.jl.
ClimateMachine.Atmos.WarmRain_1M
— TypeWarmRain_1M <: TendencyDef{Source}
A collection of source/sink terms related to 1-moment warm rain microphysics. The microphysics process rates are implemented in the CloudMicrophysics.jl Microphysics_1M module.
ClimateMachine.Atmos.RainSnow_1M
— TypeRainSnow_1M <: TendencyDef{Source}
A collection of source/sink terms related to 1-moment rain and snow microphysics The microphysics process rates are implemented in the CloudMicrophysics.jl Microphysics_1M module.
Large-scale forcing
ClimateMachine.Atmos.NoLSForcing
— TypeNoLSForcing <: LSForcingModel
No large-scale forcing
ClimateMachine.Atmos.HadGEMVertical
— TypeContainer for GCM variables from HadGEM2-A forcing,
used in the AMIP experiments