Balance Laws
The balance law
ClimateMachine.BalanceLaws.BalanceLaw — Typeabstract type BalanceLaw endAn abstract type representing a PDE balance law of the form:
where:
- $q$ is the prognostic state,
- $a$ is the auxiliary state,
- $g = G(q, a, t)$ is the gradient state (variables of which we compute the gradient),
- $h$ is the hyperdiffusive state.
Subtypes of BalanceLaw should define the following interfaces:
vars_stateto define the prognostic, auxiliary and intermediate variables.flux_first_order!to compute $F_1$flux_second_order!to compute $F_2$source!to compute $S$
If vars(bl, ::GradientFlux, FT) is non-empty, then the following should be defined:
compute_gradient_argument!to compute $G$compute_gradient_flux!is a linear transformation of $\nabla g$
If vars(bl, ::Hyperdiffusive, FT) is non-empty, then the following should be defined:
Additional functions:
wavespeedif using the Rusanov numerical flux.boundary_state!if using non-periodic boundary conditions.
State variable types
ClimateMachine.BalanceLaws.AbstractStateType — TypeAbstractStateTypeSubtypes of this describe the variables used by different parts of a BalanceLaw:
PrognosticAuxiliaryGradientGradientFluxGradientLaplacianHyperdiffusiveUpwardIntegralsDownwardIntegrals
See also vars_state.
ClimateMachine.BalanceLaws.Prognostic — TypePrognostic <: AbstractStateTypePrognostic variables in the PDE system, which are specified by the BalanceLaw, and solved for by the ODE solver.
ClimateMachine.BalanceLaws.Auxiliary — TypeAuxiliary <: AbstractStateTypeAuxiliary variables help serve several purposes:
- Pre-compute and store "expensive" variables, for example, quantities computed in vertical integrals.
- Diagnostic exports
ClimateMachine.BalanceLaws.Gradient — TypeGradient <: AbstractStateTypeVariables whose gradients must be computed.
ClimateMachine.BalanceLaws.GradientFlux — TypeGradientFlux <: AbstractStateTypeFlux variables, which are functions of gradients.
ClimateMachine.BalanceLaws.GradientLaplacian — TypeGradientLaplacian <: AbstractStateTypeGradient-Laplacian variables.
ClimateMachine.BalanceLaws.Hyperdiffusive — TypeHyperdiffusive <: AbstractStateTypeHyper-diffusive variables
ClimateMachine.BalanceLaws.UpwardIntegrals — TypeUpwardIntegrals <: AbstractStateTypeVariables computed in upward integrals
ClimateMachine.BalanceLaws.DownwardIntegrals — TypeDownwardIntegrals <: AbstractStateTypeVariables computed in downward integrals
Variable specification methods
ClimateMachine.BalanceLaws.vars_state — FunctionBalanceLaws.vars_state(::BL, ::AbstractStateType, FT)Defines the state variables of a BalanceLaw subtype BL with floating point type FT.
For each AbstractStateType, this should return a NamedTuple type, with element type either FT, an SArray with element type FT or another NamedTuple satisfying the same property.
For convenience, we recommend using the VariableTemplates.@vars macro.
Example
struct MyBalanceLaw <: BalanceLaw end
BalanceLaws.vars_state(::MyBalanceLaw, ::Prognostic, FT) =
@vars(x::FT, y::SVector{3, FT})
BalanceLaws.vars_state(::MyBalanceLaw, ::Auxiliary, FT) =
@vars(components::@vars(a::FT, b::FT))Initial condition methods
ClimateMachine.BalanceLaws.init_state_prognostic! — Functioninit_state_prognostic!(
::BL,
state_prognostic::Vars,
state_auxiliary::Vars,
coords,
args...,
)Sets the initial state of the prognostic variables state_prognostic at each node for a BalanceLaw subtype BL.
ClimateMachine.BalanceLaws.init_state_auxiliary! — Functioninit_state_auxiliary!(
bl::BalanceLaw,
f!,
statearray_auxiliary,
grid,
direction;
state_temporary = nothing
)Apply f!(bl, state_auxiliary, tmp, geom) at each node, storing the result in statearray_auxiliary, where tmp are the values at the corresponding node in state_temporary and geom contains the geometry information.
init_state_auxiliary!(
::BL,
statearray_auxiliary,
geom::LocalGeometry,
)Sets the initial state of the auxiliary variables state_auxiliary at each node for a BalanceLaw subtype BL. By default this calls nodal_init_state_auxiliary!.
init_state_auxiliary!(
m::AtmosModel,
aux::Vars,
grid,
direction
)Initialise auxiliary variables for each AtmosModel subcomponent. Store Cartesian coordinate information in aux.coord.
init_state_auxiliary!(::HBModel)sets the initial value for auxiliary variables (those that aren't related to vertical integrals) dispatches to oceaninitaux! which is defined in a problem file such as SimpleBoxProblem.jl
No need to init, initialize by full modelClimateMachine.BalanceLaws.nodal_init_state_auxiliary! — Functionnodal_init_state_auxiliary!(::BL, state_auxiliary, state_temporary, geom)Sets the initial state of the auxiliary variables state_auxiliary at each node for a BalanceLaw subtype BL.
See also init_state_auxiliary!.
Source term kernels
ClimateMachine.BalanceLaws.flux_first_order! — Functionflux_first_order!(
::BL,
flux::Grad,
state_prognostic::Vars,
state_auxiliary::Vars,
t::Real,
direction
)Sets the first-order (hyperbolic) flux terms for a BalanceLaw subtype BL.
ClimateMachine.BalanceLaws.flux_second_order! — Functionflux_second_order!(
::BL,
flux::Grad,
state_prognostic::Vars,
state_gradient_flux::Vars,
hyperdiffusive::Vars,
state_auxiliary::Vars,
t::Real
)Sets second-order (parabolic) flux terms for a BalanceLaw subtype BL.
ClimateMachine.BalanceLaws.source! — Functionsource!(
::BL,
source::Vars,
state_prognostic::Vars,
diffusive::Vars,
state_auxiliary::Vars,
t::Real
)Compute non-conservative source terms for a BalanceLaw subtype BL.
Integral kernels
ClimateMachine.BalanceLaws.indefinite_stack_integral! — Functionindefinite_stack_integral!Compute indefinite integral along stack.
ClimateMachine.BalanceLaws.reverse_indefinite_stack_integral! — Functionreverse_indefinite_stack_integral!Compute reverse indefinite integral along stack.
ClimateMachine.BalanceLaws.integral_load_auxiliary_state! — Functionintegral_load_auxiliary_state!(::BL, integrand, state_prognostic, state_aux)Specify variables integrand which will have their upward integrals computed.
See also UpwardIntegrals
ClimateMachine.BalanceLaws.integral_set_auxiliary_state! — Functionintegral_set_auxiliary_state!(::BL, state_aux, integral)Update auxiliary variables based on the upward integral integral defined in integral_load_auxiliary_state!.
ClimateMachine.BalanceLaws.reverse_integral_load_auxiliary_state! — Functionreverse_integral_load_auxiliary_state!(::BL, integrand, state_prognostic, state_aux)Specify variables integrand which will have their downward integrals computed.
See also DownwardIntegrals
ClimateMachine.BalanceLaws.reverse_integral_set_auxiliary_state! — Functionreverse_integral_set_auxiliary_state!(::BL, state_aux, integral)Update auxiliary variables based on the downward integral integral defined in reverse_integral_load_auxiliary_state!.
Gradient/Laplacian kernels
ClimateMachine.BalanceLaws.compute_gradient_flux! — Functioncompute_gradient_flux!(
::BL,
state_gradient_flux::Vars,
∇transformstate::Grad,
state_auxiliary::Vars,
t::Real
)Transformation of gradients to the diffusive variables for a BalanceLaw subtype BL. This should be a linear function of ∇transformstate
ClimateMachine.BalanceLaws.compute_gradient_argument! — Functioncompute_gradient_argument!(
::BL,
transformstate::Vars,
state_prognostic::Vars,
state_auxiliary::Vars,
t::Real
)Transformation of state variables state_prognostic to variables transformstate of which gradients are computed for a BalanceLaw subtype BL.
ClimateMachine.BalanceLaws.transform_post_gradient_laplacian! — Functiontransform_post_gradient_laplacian!(
::BL,
Qhypervisc_div::Vars,
∇Δtransformstate::Grad,
state_auxiliary::Vars,
t::Real
)Transformation of Laplacian gradients to the hyperdiffusive variables for a BalanceLaw subtype BL.
Auxiliary kernels
ClimateMachine.BalanceLaws.wavespeed — Functionwavespeed(
::BL,
n⁻,
state_prognostic::Vars,
state_auxiliary::Vars,
t::Real,
direction
)Wavespeed in the direction n⁻ for a BalanceLaw subtype BL. This is required to be defined if using a RusanovNumericalFlux numerical flux.
ClimateMachine.BalanceLaws.boundary_state! — Functionboundary_state!(
::NumericalFluxGradient,
::L,
state_prognostic⁺::Vars,
state_auxiliary⁺::Vars,
normal⁻,
state_prognostic⁻::Vars,
state_auxiliary⁻::Vars,
bctype,
t
)
boundary_state!(
::NumericalFluxFirstOrder,
::L,
state_prognostic⁺::Vars,
state_auxiliary⁺::Vars,
normal⁻,
state_prognostic⁻::Vars,
state_auxiliary⁻::Vars,
bctype,
t
)
boundary_state!(
::NumericalFluxSecondOrder,
::L,
state_prognostic⁺::Vars,
state_gradient_flux⁺::Vars,
state_auxiliary⁺:
Vars, normal⁻,
state_prognostic⁻::Vars,
state_gradient_flux⁻::Vars,
state_auxiliary⁻::Vars,
bctype,
t
)Apply boundary conditions for
NumericalFluxGradientnumerical flux (internal method)NumericalFluxFirstOrderfirst-order unknownsNumericalFluxSecondOrdersecond-order unknowns
ClimateMachine.BalanceLaws.update_auxiliary_state! — Functionupdate_auxiliary_state!(
dg::DGModel,
m::BalanceLaw,
statearray_aux,
t::Real,
elems::UnitRange,
[diffusive=false]
)Hook to update the auxiliary state variables before calling any other functions.
By default, this calls nodal_update_auxiliary_state! at each node.
If diffusive=true, then state_gradflux is also passed to nodal_update_auxiliary_state!.
update_auxiliary_state!(::HBModel)
applies the vertical filter to the zonal and meridional velocities to preserve numerical incompressibility
applies an exponential filter to θ to anti-alias the non-linear advective term
doesn't actually touch the aux variables any more, but we need a better filter interface than this anywaysClimateMachine.BalanceLaws.update_auxiliary_state_gradient! — Functionupdate_auxiliary_state_gradient!(
dg::DGModel,
m::BalanceLaw,
statearray_aux,
t::Real,
elems::UnitRange,
[diffusive=false]
)Hook to update the auxiliary state variables after the gradient computation.
By default, this calls nothing.
If diffusive=true, then state_gradflux is also passed to nodal_update_auxiliary_state!.
ClimateMachine.BalanceLaws.nodal_update_auxiliary_state! — Functionnodal_update_auxiliary_state!(::BL, state_prognostic, state_auxiliary, [state_gradflux,] t)Update the auxiliary state variables at each node for a BalanceLaw subtype BL. By default it does nothing.
Called by update_auxiliary_state!.