Balance Laws

The balance law

ClimateMachine.BalanceLaws.BalanceLawType
abstract type BalanceLaw end

An abstract type representing a PDE balance law of the form:

\[\frac{dq}{dt} = \nabla \cdot F_1(q, a, t) + \nabla \cdot F_2(q, \nabla g, h, a, t) + S(q, \nabla g, a, t)\]

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:

If vars(bl, ::GradientFlux, FT) is non-empty, then the following should be defined:

If vars(bl, ::Hyperdiffusive, FT) is non-empty, then the following should be defined:

Additional functions:

source

State variable types

ClimateMachine.BalanceLaws.AuxiliaryType
Auxiliary <: AbstractStateType

Auxiliary variables help serve several purposes:

  • Pre-compute and store "expensive" variables, for example, quantities computed in vertical integrals.
  • Diagnostic exports
source

Variable specification methods

ClimateMachine.BalanceLaws.vars_stateFunction
BalanceLaws.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))
source

Initial condition methods

ClimateMachine.BalanceLaws.init_state_auxiliary!Function
init_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.

source
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!.

source
init_state_auxiliary!(
    m::AtmosModel,
    aux::Vars,
    grid,
    direction
)

Initialise auxiliary variables for each AtmosModel subcomponent. Store Cartesian coordinate information in aux.coord.

source
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

source
No need to init, initialize by full model
source

Source term kernels

Integral kernels

Gradient/Laplacian kernels

Auxiliary kernels

ClimateMachine.BalanceLaws.wavespeedFunction
wavespeed(
    ::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.

source
ClimateMachine.BalanceLaws.boundary_state!Function
boundary_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

  • NumericalFluxGradient numerical flux (internal method)
  • NumericalFluxFirstOrder first-order unknowns
  • NumericalFluxSecondOrder second-order unknowns
source
ClimateMachine.BalanceLaws.update_auxiliary_state!Function
update_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!.

source
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 anyways
source
ClimateMachine.BalanceLaws.update_auxiliary_state_gradient!Function
update_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!.

source