ClimaLand Auxiliary Variables Cache Constructors and Update Methods
ClimaLand.initialize_auxiliary
— Functioninitialize_auxiliary(model::AbstractModel, state::NamedTuple)
Returns a NamedTuple of auxiliary variables for model
with the required structure, with values equal to similar(state)
. This assumes that all auxiliary variables are defined over the entire domain, and that all auxiliary variables have the same dimension and type. The auxiliary variables NamedTuple can also hold preallocated objects which are not Fields.
If a model has no auxiliary variables, the returned NamedTuple contains only an empty array.
The input state
is an array-like object, usually a ClimaCore Field or a Vector{FT}.
Adjustments to this - for example because different auxiliary variables have different dimensions - require defining a new method.
initialize_auxiliary(
component::AbstractCanopyComponent,
state,
)
Creates and returns a ClimaCore.Fields.FieldVector with the auxiliary variables of the canopy component component
, stored using the name of the component.
The input state
is usually a ClimaCore Field object.
initialize_auxiliary(
model::CanopyModel{FT},
coords,
) where {FT}
Creates the auxiliary state vector of the CanopyModel
and returns it as a ClimaCore.Fields.FieldVector.
The input coords
is usually a ClimaCore Field object.
This function loops over the components of the CanopyModel
and appends each component models auxiliary state vector into a single state vector, structured by component name.
ClimaLand.make_set_initial_cache
— Functionmake_set_initial_cache(model::AbstractModel)
Returns the set_initial_cache!
function, which updates the auxiliary state p
in place with the initial values corresponding to Y(t=t0) = Y0.
In principle, this function is not needed, because in the very first evaluation of either explicit_tendency
or implicit_tendency
, at t=t0, the auxiliary state is updated using the initial conditions for Y=Y0. However, without setting the initial p
state prior to running the simulation, the value of p
in the saved output at t=t0 will be unset.
Furthermore, specific methods of this function may be useful for models which store time indepedent spatially varying parameter fields in the auxiliary state. In this case, update_aux!
does not need to do anything, but they do need to be set with the initial (constant) values before the simulation can be carried out.
make_set_initial_cache(model::Union{LandModel, SoilCanopyModel})
Creates the function with arguments (p,Y0,t0) that updates the cache p
with initial values corresponding to Y0 and t0.
We require a different method from the default for a model with a canopy, so this method is for any model with type ∈ Union{LandModel, SoilCanopyModel}. This is a close copy of the method for the CanopyModel, except unpacking model.canopy
rather than using model
directly.
ClimaLand.make_set_initial_cache(model::CanopyModel)
Set the initial cache p
for the canopy model. Note that if the photosynthesis model is the P-model, then set_initial_cache!
will also run set_historical_cache!
which sets the (t-1) values for Vcmax25opt, Jmax25opt, and ξ_opt.
ClimaLand.make_update_aux
— Functionmake_update_aux(model::RichardsModel)
An extension of the function make_update_aux
, for the Richardson- Richards equation.
This function creates and returns a function which updates the auxiliary variables p.soil.variable
in place.
This has been written so as to work with Differential Equations.jl.
make_update_aux(model::EnergyHydrology)
An extension of the function make_update_aux
, for the integrated soil hydrology and energy model.
This function creates and returns a function which updates the auxiliary variables p.soil.variable
in place.
This has been written so as to work with Differential Equations.jl.
make_update_aux(model::BucketModel{FT}) where {FT}
Creates the update_aux! function for the BucketModel.
make_update_aux(model::AbstractModel)
Return an update_aux!
function that updates auxiliary parameters p
.
ClimaLand.make_update_aux(canopy::CanopyModel)
Creates the update_aux!
function for the CanopyModel
Please note that the plant hydraulics model has auxiliary variables that are updated in its prognostic compute_exp_tendency!
function. While confusing, this is better for performance as it saves looping over the state vector multiple times.
The other sub-components rely heavily on each other, so the version of the CanopyModel
with these subcomponents has a single update_aux! function, given here.
make_update_aux(model::SoilCO2Model)
An extension of the function make_update_aux
, for the soilco2 equation. This function creates and returns a function which updates the auxiliary variables p.soil.variable
in place. This has been written so as to work with Differential Equations.jl.
ClimaLand.auxiliary_vars
— Functionauxiliary_vars(soil::RichardsModel)
A function which returns the names of the auxiliary variables of RichardsModel
.
auxiliary_vars(soil::EnergyHydrology)
A function which returns the names of the auxiliary variables of EnergyHydrology
.
auxiliary_vars(::SnowModel)
Returns the auxiliary variable names for the snow model. These include
- the specific humidity at the surface of the snow (
q_sfc
, unitless), - the mass fraction in liquid water (
q_l
, unitless), - the thermal conductivity (
κ
, W/m/K), - the bulk temperature (
T
, K), - the surface temperature (
T_sfc
, K), - the snow depth (
z_snow
, m), - the bulk snow density (
ρ_snow
, kg/m^3) - the SHF, LHF, and vapor flux (
turbulent_fluxes.shf
, etc), - the net radiation (
R_n, J/m^2/s)
, - the energy flux in liquid water runoff (
energy_runoff
, J/m^2/s), - the water volume in runoff (
water_runoff
, m/s),
and the total energy and water fluxes applied to the snowpack.
Since the snow can melt completely in one timestep, we clip the water and energy fluxes such that SWE cannot become negative and U cannot become unphysical. The clipped values are what are actually applied as boundary fluxes, and are stored in applied_
fluxes.
auxiliary_vars(m::AbstractModel)
Returns the auxiliary variable symbols for the model in the form of a tuple.
ClimaLand.auxiliary_vars(::AbstractCanopyComponent)
Returns the auxiliary types of the canopy component passed in as an argument.
ClimaLand.auxiliary_vars(model::PModel)
ClimaLand.auxiliary_types(model::PModel)
ClimaLand.auxiliary_domain_names(model::PModel)
Defines the auxiliary vars of the Pmode: canopy level net photosynthesis, canopy-level gross photosynthesis (GPP
), and dark respiration at the canopy level (Rd
), and
OptVars
: a NamedTuple with keys:ξ_opt
,:Vcmax25_opt
, and:Jmax25_opt
containing the acclimated optimal values of ξ, Vcmax25, and Jmax25, respectively. These are updated using an exponential moving average (EMA) at local noon.
auxiliary_vars(canopy::CanopyModel)
Returns the auxiliary variables for the canopy model by looping over each sub-component name in canopy_components
.
This relies on the propertynames of CanopyModel
being the same as those returned by canopy_components
.
auxiliary_vars(model::PlantHydraulicsModel)
A function which returns the names of the auxiliary variables of the PlantHydraulicsModel
, the water potential ψ
(m), the volume flux*cross section fa
(1/s), and the volume flux*root cross section in the roots fa_roots
(1/s), where the cross section can be represented by an area index.
ClimaLand.auxiliary_types
— Functionauxiliary_types(soil::RichardsModel)
A function which returns the names of the auxiliary types of RichardsModel
.
auxiliary_types(soil::EnergyHydrology{FT}) where {FT}
A function which returns the types of the auxiliary variables of EnergyHydrology
.
auxiliary_types(m::AbstractModel{FT}) where {FT}
Returns the auxiliary variable types for the model in the form of a tuple.
Types provided must have ClimaCore.RecursiveApply.rzero(T::DataType)
defined. Common examples include
Float64
,Float32
for scalar variables (a scalar value at each
coordinate point)
SVector{k,Float64}
for a mutable but statically sized array of
length k
at each coordinate point.
- Note that Arrays, MVectors are not isbits and cannot be used.
Here, the coordinate points are those returned by coordinates(model).
ClimaLand.auxiliary_types(::AbstractCanopyComponent)
Returns the auxiliary types of the canopy component passed in as an argument.
auxiliary_types(canopy::CanopyModel)
Returns the auxiliary types for the canopy model by looping over each sub-component name in canopy_components
.
This relies on the propertynames of CanopyModel
being the same as those returned by canopy_components
.
ClimaLand.auxiliary_types(model::PlantHydraulicsModel{FT}) where {FT}
Defines the auxiliary types for the PlantHydraulicsModel.
ClimaLand.auxiliary_domain_names
— Functionauxiliary_domain_names(soil::RichardsModel)
A function which returns the names of the auxiliary domain names of RichardsModel
.
auxiliary_domain_names(m::AbstractModel)
Returns the domain names for the auxiliary variables in the form of a tuple.
Examples: (:surface, :surface, :subsurface).
auxiliary_domain_names(m::AbstractCanopyComponent)
Returns the domain names for the auxiliary variables in the form of a tuple.
ClimaLand.make_update_cache
— Function make_update_cache(model::AbstractModel)
A helper function which updates all cache variables of a model; currently only used in set_initial_cache
since not all cache variables are updated at the same time.
ClimaLand.diffusive_flux
— Functiondiffusive_flux(K, x_2, x_1, Δz)
Calculates the diffusive flux of a quantity x (water content, temp, etc). Here, x2 = x(z + Δz) and x1 = x(z), so x_2 is at a larger z by convention.
ClimaLand.boundary_flux!
— Functionboundary_flux!(bc_field, bc::WaterFluxBC, _...)
A method of boundary fluxes which updates the desired flux.
boundary_flux!(bc_field, bc::RichardsAtmosDrivenFluxBC,
boundary::ClimaLand.AbstractBoundary,
model::RichardsModel{FT},
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
) where {FT}
A method of boundary fluxes which returns the desired water volume flux for the RichardsModel, at the top of the domain, in the case of a prescribed precipitation flux.
If model.runoff
is not of type NoRunoff
, surface runoff is accounted for when computing the infiltration.
boundary_flux!(bc_field, rre_bc::MoistureStateBC,
::ClimaLand.TopBoundary,
model::AbstractSoilModel,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of boundary fluxes which converts a state boundary condition on θ_l at the top of the domain into a flux of liquid water.
boundary_flux!(bc_field, rre_bc::MoistureStateBC,
::ClimaLand.BottomBoundary,
model::AbstractSoilModel,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of boundary fluxes which converts a state boundary condition on θ_l at the bottom of the domain into a flux of liquid water.
boundary_flux!(bc_field, bc::FreeDrainage,
boundary::ClimaLand.BottomBoundary,
model::AbstractSoilModel,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of boundary fluxes which enforces free drainage at the bottom of the domain.
boundary_flux!(bc_field, bc::HeatFluxBC, _...)
A method of boundary fluxes which updates the desired flux.
boundary_flux!(bc_field, heat_bc::TemperatureStateBC,
::ClimaLand.TopBoundary,
model::EnergyHydrology,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
):
A method of boundary fluxes which converts a state boundary condition on temperature at the top of the domain into a flux of energy.
boundary_flux!(bc_field, heat_bc::TemperatureStateBC,
::ClimaLand.BottomBoundary,
model::EnergyHydrology,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of boundary fluxes which converts a state boundary condition on temperature at the bottom of the domain into a flux of energy.
boundary_flux!(bc_field, bc::AbstractBC, bound_type::AbstractBoundary, Δz, _...)
A function which updates bc_field with the correct boundary flux given any boundary condition (BC).
function ClimaLand.boundary_flux!(bc_field,
bc::RunoffBC,
::TopBoundary,
model::Soil.RichardsModel,
Δz::FT,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
params,
)
Extension of the ClimaLand.boundary_flux
function, which returns the water volume boundary flux for the soil. At the top boundary, return the soil infiltration (computed each step and stored in p.soil_infiltration
).
ClimaLand.boundary_flux!(bc_field,
bc::SoilCO2FluxBC,
boundary::ClimaLand.AbstractBoundary,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of ClimaLand.boundary_flux which updates the soilco2 flux (kg CO2 /m^2/s) in the case of a prescribed flux BC at either the top or bottom of the domain.
ClimaLand.boundary_flux!(bc_field,
bc::SoilCO2StateBC,
boundary::ClimaLand.TopBoundary,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of ClimaLand.boundary_flux which returns the soilco2 flux in the case of a prescribed state BC at top of the domain.
ClimaLand.boundary_flux!(bc_field,
bc::SoilCO2StateBC,
boundary::ClimaLand.BottomBoundary,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of ClimaLand.boundary_flux which returns the soilco2 flux in the case of a prescribed state BC at bottom of the domain.
ClimaLand.boundary_flux!(bc_field,
bc::AtmosCO2StateBC,
boundary::ClimaLand.TopBoundary,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t,
)
A method of ClimaLand.boundary_flux which returns the soilco2 flux in the case when the atmospheric CO2 is ued at top of the domain.
ClimaLand.set_dfluxBCdY!
— FunctionClimaLand.set_dfluxBCdY!(
model::RichardsModel,
::MoistureStateBC,
boundary::ClimaLand.TopBoundary,
Δz,
Y,
p,
t,
)
Computes the derivative of the flux in the top layer (due to the boundary condition), with respect to the state variable in the top layer. This value is then updated in-place in the cache.
For Richards equation (a diffusion equation with a single state variable), this is given by ∂F_bc/∂Y_N= -K_N (∂ψ_bc/∂ϑ_N) / Δz
, where N
indicates the top layer cell index and ψ_bc
is the pressure head at the boundary condition.
set_dfluxBCdY!(::AbstractModel,
::AbstractBC,
::AbstractBoundary,
_...)::Union{ClimaCore.Fields.FieldVector, Nothing}
A function stub which returns the derivative of the implicit tendency term of the model
arising from the boundary condition, with respect to the state Y.
ClimaLand.make_update_boundary_fluxes
— Functionmake_update_boundary_fluxes(model::AbstractModel)
Return an update_boundary_fluxes!
function that updates the auxiliary parameters in p
corresponding to boundary fluxes or interactions between componets..
make_update_boundary_fluxes(
land::LandHydrology{FT, SM, SW},
) where {FT, SM <: Soil.RichardsModel{FT}, SW <: Pond.PondModel{FT}}
A method which makes a function; the returned function updates the auxiliary variable p.soil_infiltration
, which is needed for both the boundary condition for the soil model and the source term (runoff) for the surface water model.
This function is called each ode function evaluation.
make_update_boundary_fluxes(
land::SoilCanopyModel{FT, MM, SM, RM},
) where {
FT,
MM <: Soil.Biogeochemistry.SoilCO2Model{FT},
SM <: Soil.RichardsModel{FT},
RM <: Canopy.CanopyModel{FT}
}
A method which makes a function; the returned function updates the additional auxiliary variables for the integrated model, as well as updates the boundary auxiliary variables for all component models.
This function is called each ode function evaluation, prior to the tendency function evaluation.
make_update_boundary_fluxes(
land::SoilSnowModel{FT, SnM, SoM},
) where {
FT,
SnM <: Snow.SnowModel{FT},
SoM <: Soil.EnergyHydrology{FT},
}
A method which makes a function; the returned function updates the additional auxiliary variables for the integrated model, as well as updates the boundary auxiliary variables for all component models.
This function is called each ode function evaluation, prior to the tendency function evaluation.
In this method, we
- Compute the ground heat flux between soil and snow. This is required to update the snow and soil boundary fluxes
- Update the snow boundary fluxes, which also computes any excess flux of energy or water which occurs when the snow
completely melts in a step. In this case, that excess must go to the soil for conservation
- Update the soil boundary fluxes use precomputed ground heat flux and excess fluxes from snow.
- Compute the net flux for the atmosphere, which is useful for assessing conservation.
make_update_boundary_fluxes(
land::LandModel{FT, MM, SM, RM, SnM},
) where {
FT,
MM <: Soil.Biogeochemistry.SoilCO2Model{FT},
SM <: Soil.RichardsModel{FT},
RM <: Canopy.CanopyModel{FT}
SnM <: Snow.SnowModel{FT}
}
A method which makes a function; the returned function updates the additional auxiliary variables for the integrated model, as well as updates the boundary auxiliary variables for all component models.
This function is called each ode function evaluation, prior to the tendency function evaluation.