Solver Functions
Prognostic Variables
ClimaLand.prognostic_vars
— Functionprognostic_vars(soil::RichardsModel)
A function which returns the names of the prognostic variables of RichardsModel
.
prognostic_vars(soil::EnergyHydrology)
A function which returns the names of the prognostic variables of EnergyHydrology
.
prognostic_vars(::SnowModel)
Returns the prognostic variable names of the snow model.
For this model, we track the snow water equivalent S in meters (liquid water volume per ground area) and the energy per unit ground area U [J/m^2] prognostically.
prognostic_vars(m::AbstractModel)
Returns the prognostic variable symbols for the model in the form of a tuple.
Note that this default suggests that a model has no prognostic variables, which is an invalid model setup. This function is meant to be extended for all models.
ClimaLand.prognostic_vars(::AbstractCanopyComponent)
Returns the prognostic vars of the canopy component passed in as an argument.
prognostic_vars(canopy::CanopyModel)
Returns the prognostic 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
.
prognostic_vars(model::PlantHydraulicsModel)
A function which returns the names of the prognostic variables of the PlantHydraulicsModel
.
ClimaLand.prognostic_types
— Functionprognostic_types(soil::EnergyHydrology{FT}) where {FT}
A function which returns the types of the prognostic variables of EnergyHydrology
.
prognostic_types(::SnowModel{FT})
Returns the prognostic variable types of the snow model; both snow water equivalent and energy per unit area are scalars.
prognostic_types(m::AbstractModel{FT}) where {FT}
Returns the prognostic 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.
Here, the coordinate points are those returned by coordinates(model).
Note that this default suggests that a model has no prognostic variables, which is an invalid model setup. This function is meant to be extended for all models.
ClimaLand.prognostic_types(::AbstractCanopyComponent)
Returns the prognostic types of the canopy component passed in as an argument.
prognostic_types(canopy::CanopyModel)
Returns the prognostic 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.prognostic_types(model::PlantHydraulicsModel{FT}) where {FT}
Defines the prognostic types for the PlantHydraulicsModel.
ClimaLand.prognostic_domain_names
— Functionprognostic_domain_names(::SnowModel)
Returns the prognostic variable domain names of the snow model; both snow water equivalent and energy per unit area are modeling only as a function of (x,y), and not as a function of depth. Therefore their domain name is ":surface".
prognostic_domain_names(m::AbstractModel)
Returns the domain names for the prognostic variables in the form of a tuple.
Examples: (:surface, :surface, :subsurface).
Note that this default suggests that a model has no prognostic variables, which is an invalid model setup. This function is meant to be extended for all models.
prognostic_domain_names(m::AbstractCanopyComponent)
Returns the domain names for the prognostic variables in the form of a tuple.
ClimaLand.initialize_prognostic
— Functioninitialize_prognostic(model::AbstractModel, state::NamedTuple)
Returns a FieldVector of prognostic variables for model
with the required structure, with values equal to similar(state)
. This assumes that all prognostic variables are defined over the entire domain, and that all prognostic variables have the same dimension and type.
If a model has no prognostic variables, the returned FieldVector 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 prognostic variables have different dimensions - require defining a new method.
initialize_prognostic(
component::AbstractCanopyComponent,
state,
)
Creates and returns a ClimaCore.Fields.FieldVector with the prognostic variables of the canopy component component
, stored using the name of the component.
The input state
is usually a ClimaCore Field object.
initialize_prognostic(
model::CanopyModel{FT},
coords,
) where {FT}
Creates the prognostic state vector of the CanopyModel
and returns it as a ClimaCore.Fields.FieldVector.
The input state
is usually a ClimaCore Field object.
This function loops over the components of the CanopyModel
and appends each component models prognostic state vector into a single state vector, structured by component name.
Sources
ClimaLand.AbstractSource
— TypeAbstractSource{FT <: AbstractFloat}
An abstract type for types of source terms.
ClimaLand.source!
— Functionsource!(dY::ClimaCore.Fields.FieldVector,
src::PhaseChange{FT},
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
model
)
Computes the source terms for phase change explicitly in time.
source!(dY::ClimaCore.Fields.FieldVector,
src::SoilSublimation{FT},
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
model
)
Updates dY.soil.θ_i in place with a term due to sublimation; this only affects the surface layer of soil.
source!(dY::ClimaCore.Fields.FieldVector,
src::AbstractSource,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple
)::ClimaCore.Fields.Field
A stub function, which is extended by ClimaLand.
source!(dY::ClimaCore.Fields.FieldVector,
src::SoilSublimationwithSnow{FT},
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
model
)
Updates dY.soil.θ_i
in place with a term due to sublimation; this only affects the surface layer of soil.
ClimaLand.source!(dY::ClimaCore.Fields.FieldVector,
src::RootExtraction,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple
model::EnergyHydrology)
An extension of the ClimaLand.source!
function, which computes source terms for the soil model; this method returns the water and energy loss/gain due to root extraction.
ClimaLand.source!(dY::ClimaCore.Fields.FieldVector,
src::MicrobeProduction,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
params)
A method which extends the ClimaLand source! function for the case of microbe production of CO2 in soil.
ClimaLand.source!(
dY::ClimaCore.Fields.FieldVector,
src::TOPMODELSubsurfaceRunoff,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
model::AbstractSoilModel{FT},
) where {FT}
Adjusts dY.soil.ϑ_l in place to account for the loss of water due to subsurface runoff.
The sink term is given by - Rss/h∇ H(twc - ν), where H is the Heaviside function, h∇ is the water table thickness (defined to be where twc>ν), where twc is the total water content, and Rss is the runoff as a flux(m/s).
Boundary Conditions
ClimaLand.AbstractBC
— TypeAbstractBC
An abstract type for types of boundary conditions, which will include prescribed functions of space and time as Dirichlet conditions or Neumann conditions, in addition to other convenient conditions.
ClimaLand.AbstractBoundary
— TypeAbstractBoundary
An abstract type to indicate which boundary we are doing calculations for. Currently, we support the top boundary (TopBoundary) and bottom boundary (BottomBoundary).
ClimaLand.TopBoundary
— TypeTopBoundary{} <: AbstractBoundary{}
A simple object which should be passed into a function to indicate that we are considering the top boundary.
ClimaLand.BottomBoundary
— TypeBottomBoundary{} <: AbstractBoundary{}
A simple object which should be passed into a function to indicate that we are considering the bottom boundary.
ClimaLand.boundary_vars
— Functionboundary_vars(::RichardsAtmosDrivenFluxBC{<:PrescribedPrecipitation,
<:Runoff.AbstractRunoffModel,
}, ::ClimaLand.TopBoundary)
An extension of the boundary_vars
method for RichardsAtmosDrivenFluxBC with runoff.
These variables are updated in place in boundary_flux!
.
boundary_vars(::AtmosDrivenFluxBC, ::ClimaLand.TopBoundary)
An extension of the boundary_vars
method for AtmosDrivenFluxBC. This adds the surface conditions (SHF, LHF, evaporation, and resistance) and the net radiation to the auxiliary variables.
These variables are updated in place in soil_boundary_fluxes!
.
boundary_vars(::MoistureStateBC, ::ClimaLand.TopBoundary)
An extension of the boundary_vars
method for MoistureStateBC at the top boundary.
These variables are updated in place in boundary_flux!
.
boundary_vars(::AbstractBC , ::ClimaLand.TopBoundary)
The list of symbols for additional variables to add to the model auxiliary state, for models solving PDEs, which defaults to adding storage for the top boundary flux fields, but which can be extended depending on the type of boundary condition used.
For the Soil and SoilCO2 models - which solve PDEs - the tendency functions and updateboundaryfluxes functions are coded to access the field :top_bc
to be present in the model cache, which is why this is the default. If this is not your (PDE) model's desired behavior, you can extend this function with a new method.
The field top_bc_wvec
is created to prevent allocations only; it is used in the tendency function only.
Use this function in the exact same way you would use auxiliary_vars
.
boundary_vars(::AbstractBC, ::ClimaLand.BottomBoundary)
The list of symbols for additional variables to add to the model auxiliary state, for models solving PDEs, which defaults to adding storage for the bottom boundary flux fields, but which can be extended depending on the type of boundary condition used.
For the Soil and SoilCO2 models - which solve PDEs - the tendency functions and updateboundaryfluxes functions are coded to access the field :bottom_bc
to be present in the model cache, which is why this is the default. If this is not your (PDE) model's desired behavior, you can extend this function with a new method.
The field bottom_bc_wvec
is created to prevent allocations only; it is used in the tendency function only.
Use this function in the exact same way you would use auxiliary_vars
.
ClimaLand.boundary_var_domain_names
— Functionboundary_var_domain_names(::RichardsAtmosDrivenFluxBC{<:PrescribedPrecipitation,
<:Runoff.AbstractRunoffModel,
},
::ClimaLand.TopBoundary)
An extension of the boundary_var_domain_names
method for RichardsAtmosDrivenFluxBC with runoff.
boundary_var_domain_names(::AtmosDrivenFluxBC,
::ClimaLand.TopBoundary)
An extension of the boundary_var_domain_names
method for AtmosDrivenFluxBC. This specifies the part of the domain on which the additional variables should be defined.
boundary_var_domain_names(::MoistureStateBC, ::ClimaLand.TopBoundary)
An extension of the boundary_var_domain_names
method for MoistureStateBC at the top boundary.
boundary_var_domain_names(::AbstractBC, ::ClimaLand.AbstractBoundary)
The list of domain names for additional variables to add to the model auxiliary state, for models solving PDEs, which defaults to adding storage on the surface domain for the top or bottom boundary flux fields, but which can be extended depending on the type of boundary condition used.
Use in conjunction with boundary_vars
, in the same way you would use auxiliary_var_domain_names
.
ClimaLand.boundary_var_types
— Functionboundary_var_types(::RichardsModel{FT},
::RichardsAtmosDrivenFluxBC{<:PrescribedPrecipitation,
<: Runoff.AbstractRunoffModel,
},
::ClimaLand.TopBoundary,
) where {FT}
An extension of the boundary_var_types
method for RichardsAtmosDrivenFluxBC with runoff.
boundary_var_types(::Soil.EnergyHydrology{FT}, ::AbstractEnergyHydrologyBC, ::ClimaLand.AbstractBoundary) where {FT}
The list of domain names for additional variables added to the EnergyHydrology model auxiliary state, which defaults to adding storage for the boundary flux field.
Because we supply boundary conditions for water and heat, we found it convenient to have these stored as a NamedTuple under the names top_bc
and bottom_bc
.
boundary_var_types(
::EnergyHydrology{FT},
::AtmosDrivenFluxBC,
::ClimaLand.TopBoundary,
) where {FT}
An extension of the boundary_var_types
method for AtmosDrivenFluxBC. This specifies the type of the additional variables.
boundary_var_types(
::EnergyHydrology{FT},
::AtmosDrivenFluxBC{<:CoupledAtmosphere, <:CoupledRadiativeFluxes},
::ClimaLand.TopBoundary,
) where {FT}
An extension of the boundary_var_types
method for AtmosDrivenFluxBC with coupled atmosphere and radiative fluxes. This specifies the type of the additional variables.
This method includes additional fluxes needed by the atmosphere: momentum fluxes (ρτxz
, ρτyz
) and the buoyancy flux (buoy_flux
). These are updated in place when the coupler computes turbulent fluxes, rather than in soil_boundary_fluxes!
.
Note that we currently store these in the land model because the coupler computes turbulent land/atmosphere fluxes using ClimaLand functions, and the land model needs to be able to store the fluxes as an intermediary. Once we compute fluxes entirely within the coupler, we can remove this.
boundary_var_types(::RichardsModel{FT},
::MoistureStateBC,
::ClimaLand.TopBoundary,
) where {FT}
An extension of the boundary_var_types
method for MoistureStateBC at the top boundary.
boundary_var_types(model::AbstractModel{FT}, ::AbstractBC, ::ClimaLand.AbstractBoundary) where {FT}
The list of types for additional variables to add to the model auxiliary state, for models solving PDEs, which defaults to adding a scalar variable on the surface domain for the top or bottom boundary flux fields, but which can be extended depending on the type of boundary condition used.
Use in conjunction with boundary_vars
, in the same way you would use auxiliary_var_types
. The use of a scalar is appropriate for models with a single PDE; models with multiple PDEs will need to supply multiple scalar fields.
Implicit Tendencies
ClimaLand.make_imp_tendency
— Functionmake_imp_tendency(model::AbstractImExModel)
Returns an imp_tendency
that updates auxiliary variables and updates the prognostic state of variables that are stepped implicitly.
compute_imp_tendency!
should be compatible with SciMLBase.jl solvers.
make_imp_tendency(model::AbstractModel)
Returns an imp_tendency
that does nothing. This model type is not stepped explicity.
ClimaLand.make_compute_imp_tendency
— Functionmake_compute_imp_tendency(model::RichardsModel)
An extension of the function make_compute_imp_tendency
, for the Richardson- Richards equation.
This function creates and returns a function which computes the entire right hand side of the PDE for ϑ_l
, and updates dY.soil.ϑ_l
in place with that value.
make_compute_imp_tendency(model::EnergyHydrology)
An extension of the function make_compute_imp_tendency
, for the integrated soil energy and heat equations, including phase change.
This version of this function computes the right hand side of the PDE for Y.soil.ϑ_l
, which is the only quantity we currently step implicitly.
This has been written so as to work with Differential Equations.jl.
make_compute_imp_tendency(model::AbstractModel)
Return a compute_imp_tendency!
function that updates state variables that we will be stepped implicitly. This fallback sets all tendencies of this model to zero, which is appropriate for models that do not have any implicit tendencies to update. Note that we cannot set dY .= 0
here because this would overwrite the tendencies of all models in the case of an integrated LSM.
compute_imp_tendency!
should be compatible with SciMLBase.jl solvers.
ClimaLand.make_compute_imp_tendency(component::AbstractCanopyComponent, canopy)
Creates the computeimptendency!(dY,Y,p,t) function for the canopy component
.
Since component models are not standalone models, other information may be needed and passed in (via the canopy
model itself). The right hand side for the entire canopy model can make use of these functions for the individual components.
make_compute_imp_tendency(canopy::CanopyModel)
Creates and returns the computeimptendency! for the CanopyModel
.
ClimaLand.make_compute_jacobian
— FunctionClimaLand.make_compute_jacobian(model::RichardsModel{FT}) where {FT}
Creates and returns the compute_jacobian! function for RichardsModel. This updates the contribution for the soil liquid water content.
Using this Jacobian with a backwards Euler timestepper is equivalent to using the modified Picard scheme of Celia et al. (1990).
ClimaLand.make_compute_jacobian(model::EnergyHydrology{FT}) where {FT}
Creates and returns the compute_jacobian! function for the EnergyHydrology model. This updates the contribution for the soil liquid water content only.
Using this Jacobian with a backwards Euler timestepper is equivalent to using the modified Picard scheme of Celia et al. (1990).
make_compute_jacobian(model::AbstractModel)
Creates and returns a function which computes the entries of the Jacobian matrix W
in place.
If the implicit tendency function is given by T!(dY, Y, p, t) = make_implicit_tendency(model)
, the Jacobian should be given by W_{i,j}! = ∂T!_i/∂Y_j
, where Y_j
is the j-th
state variable and T!_i
is the implicit tendency of the i-th
state variable.
The default is that no updates are required, but this function must be extended for models that use implicit timestepping.
ClimaLand.make_compute_jacobian(canopy::CanopyModel)
Creates and returns the compute_jacobian! for the CanopyModel
.
ClimaLand.make_jacobian
— Functionmake_jacobian(model::AbstractModel)
Creates and returns a function which updates the auxiliary variables p
in place and then updates the entries of the Jacobian matrix W
for the model
in place.
The default is that no updates are required, no implicit tendency is present, and hence the timestepping is entirely explicit.
Note that the returned function jacobian!
should be used as Wfact!
in ClimaTimeSteppers.jl
and SciMLBase.jl
.
ClimaLand.FieldMatrixWithSolver
— FunctionFieldMatrixWithSolver(Y::ClimaCore.Fields.FieldVector)
Outer constructor for the FieldMatrixWithSolver Jacobian matrix struct. This extends the constructor from ClimaCore.FieldMatrix, filling the object with ClimaLand-specific values.
For variables that will be stepped implicitly, the Jacobian matrix is a tridiagonal matrix. For variables that will be stepped explicitly, the Jacobian matrix is a negative identity matrix.
To run a model with one or more prognostic variables stepped implicitly, the Jacobian matrix must be constructed and passed to the solver. All implicitly-stepped variables of the model should be added to the implicit_names
tuple, and any explicitly-stepped variables should be added to the explicit_names
tuple.
Explicit Tendencies
ClimaLand.make_exp_tendency
— Functionmake_exp_tendency(model::AbstractModel)
Returns an exp_tendency
that updates auxiliary variables and updates the prognostic state of variables that are stepped explicitly.
compute_exp_tendency!
should be compatible with SciMLBase.jl solvers.
ClimaLand.make_compute_exp_tendency
— Functionmake_explicit_tendency(model::Soil.RichardsModel)
An extension of the function make_compute_imp_tendency
, for the Richardson- Richards equation.
Construct the tendency computation function for the explicit terms of the RHS, which are horizontal components and source/sink terms.
make_compute_exp_tendency(model::EnergyHydrology)
An extension of the function make_compute_exp_tendency
, for the integrated soil energy and heat equations, including phase change.
This function creates and returns a function which computes the entire right hand side of the PDE for Y.soil.ϑ_l, Y.soil.θ_i, Y.soil.ρe_int
, and updates dY.soil
in place with those values. All of these quantities will be stepped explicitly.
This has been written so as to work with Differential Equations.jl.
make_compute_exp_tendency(model::BucketModel{FT}) where {FT}
Creates the computeexptendency! function for the bucket model.
make_compute_exp_tendency(model::AbstractModel)
Return a compute_exp_tendency!
function that updates state variables that we will be stepped explicitly. This fallback sets all tendencies of this model to zero, which is appropriate for models that do not have any explicit tendencies to update. Note that we cannot set dY .= 0
here because this would overwrite the tendencies of all models in the case of an integrated LSM.
compute_exp_tendency!
should be compatible with SciMLBase.jl solvers.
ClimaLand.make_compute_exp_tendency(component::AbstractCanopyComponent, canopy)
Creates the computeexptendency!(dY,Y,p,t) function for the canopy component
.
Since component models are not standalone models, other information may be needed and passed in (via the canopy
model itself). The right hand side for the entire canopy model can make use of these functions for the individual components.
make_compute_exp_tendency(canopy::CanopyModel)
Creates and returns the computeexptendency! for the CanopyModel
.
make_compute_exp_tendency(model::SoilCO2Model)
An extension of the function make_compute_exp_tendency
, for the soilco2 equation. This function creates and returns a function which computes the entire right hand side of the PDE for C
, and updates dY.soil.C
in place with that value. These quantities will be stepped explicitly.
This has been written so as to work with Differential Equations.jl.
make_compute_exp_tendency(model::PlantHydraulicsModel, _)
A function which creates the computeexptendency! function for the PlantHydraulicsModel. The computeexptendency! function must comply with a rhs function of SciMLBase.jl.
Below, fa
denotes a flux multiplied by the relevant cross section (per unit area ground, or area index, AI). The tendency for the ith compartment can be written then as: ∂ϑ[i]/∂t = 1/(AI*dz)[fa[i]-fa[i+1]).
Note that if the area_index is zero because no plant is present, AIdz is zero, and the fluxes fa
appearing in the numerator are zero because they are scaled by AI.
To prevent dividing by zero, we change AI/(AI x dz)" to "AI/max(AI x dz, eps(FT))"