Shared Utilities

Domains

ClimaLSM.Domains.SphericalShellType
struct SphericalShell{FT} <: AbstractDomain{FT}
    radius::FT
    depth::FT
    dz_tuple::Union{Tuple{FT, FT}, Nothing}
    nelements::Tuple{Int, Int}
    npolynomial::Int
end

A struct holding the necessary information to construct a domain, a mesh, a 2d spectral element space (non-radial directions) x a 1d finite difference space (radial direction), and the resulting coordinate field.

Fields

  • radius: The radius of the shell

  • depth: The radial extent of the shell

  • dz_tuple: Tuple for mesh stretching specifying target (dzbottom, dztop) (m). If nothing, no stretching is applied.

  • nelements: The number of elements to be used in the non-radial and radial directions

  • npolynomial: The polynomial order to be used in the non-radial directions

  • space: The associated ClimaCore Space

source
ClimaLSM.Domains.SphericalSurfaceType
struct SphericalSurface{FT} <: AbstractDomain{FT}
    radius::FT
    nelements::Tuple{Int, Int}
    npolynomial::Int
end

A struct holding the necessary information to construct a domain, a mesh, a 2d spectral element space (non-radial directions) and the resulting coordinate field.

Fields

  • radius: The radius of the surface

  • nelements: The number of elements to be used in the non-radial directions

  • npolynomial: The polynomial order to be used in the non-radial directions

  • space: The associated ClimaCore Space

source
ClimaLSM.Domains.HybridBoxType
struct HybridBox{FT} <: AbstractDomain{FT}
    xlim::Tuple{FT, FT}
    ylim::Tuple{FT, FT}
    zlim::Tuple{FT, FT}
    dz_tuple::Union{Tuple{FT, FT}, Nothing}
    nelements::Tuple{Int, Int, Int}
    npolynomial::Int
    periodic::Tuple{Bool, Bool}
end

A struct holding the necessary information to construct a domain, a mesh, a 2d spectral element space (horizontal) x a 1d finite difference space (vertical), and the resulting coordinate field. This domain is not periodic along the z-axis. Note that only periodic domains are supported in the horizontal.

Fields

  • xlim: Domain interval limits along x axis, in meters

  • ylim: Domain interval limits along y axis, in meters

  • zlim: Domain interval limits along z axis, in meters

  • dz_tuple: Tuple for mesh stretching specifying target (dzbottom, dztop) (m). If nothing, no stretching is applied.

  • nelements: Number of elements to discretize interval, (nx, ny,nz)

  • npolynomial: Polynomial order for the horizontal directions

  • periodic: Flag indicating periodic boundaries in horizontal. only true is supported

  • space: The associated ClimaCore Space

source
ClimaLSM.Domains.ColumnType
Column{FT} <: AbstractDomain{FT}

A struct holding the necessary information to construct a domain, a mesh, a center and face space, etc. for use when a finite difference in 1D is suitable, as for a soil column model.

Fields

  • zlim: Domain interval limits, (zmin, zmax), in meters

  • nelements: Number of elements used to discretize the interval

  • dz_tuple: Tuple for mesh stretching specifying target (dzbottom, dztop) (m). If nothing, no stretching is applied.

  • boundary_tags: Boundary face identifiers

  • space: The associated ClimaCore Space

source
ClimaLSM.Domains.PlaneType
Plane{FT} <: AbstractDomain{FT}

A struct holding the necessary information to construct a domain, a mesh, a 2d spectral element space, and the resulting coordinate field. Note that only periodic domains are currently supported.

Fields

  • xlim: Domain interval limits along x axis, in meters

  • ylim: Domain interval limits along y axis, in meters

  • nelements: Number of elements to discretize interval, (nx, ny)

  • periodic: Flags for periodic boundaries; only true is supported

  • npolynomial: Polynomial order for both x and y

  • space: The associated ClimaCore Space

source
ClimaLSM.Domains.PointType
Point{FT} <: AbstractDomain{FT}

A domain for single column surface variables. For models such as ponds, snow, plant hydraulics, etc. Enables consistency in variable initialization across all domains.

Fields

  • z_sfc: Surface elevation relative to a reference (m)

  • space: The associated ClimaCore Space

source
ClimaLSM.Domains.LSMSingleColumnDomainType
LSMSingleColumnDomain{FT} <: AbstractLSMDomain{FT}

A mixed domain, consisting of a column domain with z-coordinates at the finite difference cell centers, and a point domain, with a single z coordinate at the top boundary of the column domain. For use in LSM modeling, where a subsurface finite difference space (for modeling soil hydrology and energy) and a surface space are both needed.

Fields

  • subsurface: The subsurface Column domain

  • surface: The surface Point domain

source
ClimaLSM.Domains.LSMMultiColumnDomainType
LSMMultiColumnDomain{FT} <: AbstractLSMDomain{FT}

A mixed domain, consisting of a hybdrid box domain with z-coordinates at the finite difference cell centers, and a plane domain, coinciding with the surface of the box.

For use in LSM modeling, where a subsurface finite difference space (for modeling soil hydrology and energy) and a surface space are both needed.

Fields

  • subsurface: The subsurface box domain

  • surface: The surface plane domain

source
ClimaLSM.Domains.LSMSphericalShellDomainType
LSMSphericalShellDomain{FT} <: AbstractLSMDomain{FT}

A mixed domain, consisting of a spherical shell domain with z-coordinates at the finite difference cell centers, and a spherical surface domain, coinciding with the surface of the shell.

For use in LSM modeling, where a subsurface finite difference space (for modeling soil hydrology and energy) and a surface space are both needed.

Fields

  • subsurface: The subsurface shell domain

  • surface: The surface domain

source
ClimaLSM.Domains.coordinatesFunction
coordinates(domain::AbstractDomain)

Returns the coordinate field for the domain.

source
coordinates(domain::AbstractLSMDomain{FT}) where {FT}

Returns the coordinates of the AbstractLSMDomain as a named tuple, with keys of subsurface and surface.

source
ClimaLSM.Domains.obtain_face_spaceFunction
obtain_face_space(cs::ClimaCore.Spaces.AbstractSpace)

Returns the face space, if applicable, for the center space cs.

source
obtain_face_space(cs::ClimaCore.Spaces.CenterExtrudedFiniteDifferenceSpace)

Returns the face space for the CenterExtrudedFiniteDifferenceSpace cs.

source
obtain_face_space(cs::ClimaCore.Spaces.CenterFiniteDifferenceSpace)

Returns the face space corresponding to the CenterFiniteDifferenceSpace cs.

source
ClimaLSM.Domains.obtain_surface_spaceFunction
obtain_surface_space(cs::ClimaCore.Spaces.AbstractSpace)

Returns the surface space, if applicable, for the center space cs.

source
obtain_surface_space(cs::ClimaCore.Spaces.CenterExtrudedFiniteDifferenceSpace)

Returns the horizontal space for the CenterExtrudedFiniteDifferenceSpace cs.

source
obtain_surface_space(cs::ClimaCore.Spaces.CenterFiniteDifferenceSpace)

Returns the top level of the face space corresponding to the CenterFiniteDifferenceSpace cs.

source
ClimaLSM.Domains.top_center_to_surfaceFunction
top_center_to_surface(center_field::ClimaCore.Fields.Field)

Creates and returns a ClimaCore.Fields.Field defined on the space corresponding to the surface of the space on which center_field is defined, with values equal to the those at the level of the top center.

For example, given a center_field defined on 1D center finite difference space, this would return a field defined on the Point space of the surface of the column. The value would be the value of the oroginal center_field at the topmost location. Given a center_field defined on a 3D extruded center finite difference space, this would return a 2D field corresponding to the surface, with values equal to the topmost level.

source

Models

ClimaLSM.AbstractImExModelType
AbstractImExModel{FT} <: AbstractModel{FT}

An abstract type for models which must be treated implicitly (and which may also have tendency terms that can be treated explicitly). This inherits all the default function definitions from AbstractModel, as well as make_imp_tendency and make_compute_imp_tendency defaults.

source
ClimaLSM.AbstractExpModelType
AbstractExpModel{FT} <: AbstractModel{FT}

An abstract type for models which must be treated explicitly. This inherits all the default function definitions from AbstractModel, as well as a make_imp_tendency default.

source
ClimaLSM.make_exp_tendencyFunction
make_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.

source
ClimaLSM.make_imp_tendencyFunction
make_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.

source
make_imp_tendency(model::AbstractModel)

Returns an imp_tendency that does nothing. This model type is not stepped explicity.

source
ClimaLSM.make_compute_exp_tendencyFunction
make_compute_exp_tendency(model::BucketModel{FT}) where {FT}

Creates the computeexptendency! function for the bucket model.

source
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))"

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

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

source
make_compute_exp_tendency(canopy::CanopyModel)

Creates and returns the computeexptendency! for the CanopyModel.

This allows for prognostic variables in each canopy component, and specifies that they will be stepped explicitly.

source
make_compute_exp_tendency(model::AbstractModel)

Return a compute_exp_tendency! function that updates state variables that we will be stepped explicitly.

compute_exp_tendency! should be compatible with SciMLBase.jl solvers.

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

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

source
ClimaLSM.make_compute_imp_tendencyFunction
make_compute_imp_tendency(model::AbstractModel)

Return a compute_imp_tendency! function that updates state variables that we will be stepped implicitly.

compute_imp_tendency! should be compatible with SciMLBase.jl solvers.

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

source
ClimaLSM.make_update_auxFunction
make_update_aux(model::BucketModel{FT}) where {FT}

Creates the update_aux! function for the BucketModel.

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

source
 ClimaLSM.make_update_aux(canopy::CanopyModel{FT, <:BeerLambertModel,
                                              <:FarquharModel,
                                              <:MedlynConductanceModel,
                                              <:PlantHydraulicsModel,},
                          ) where {FT}

Creates the update_aux! function for the CanopyModel; a specific method for update_aux! for the case where the canopy model components are of the type in the parametric type signature: AbstractRadiationModel, FarquharModel, MedlynConductanceModel, and PlantHydraulicsModel.

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.

source
make_update_aux(model::AbstractModel)

Return an update_aux! function that updates auxiliary parameters p.

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

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

source
ClimaLSM.make_set_initial_aux_stateFunction
ClimaLSM.make_set_initial_aux_state(model::BucketModel{FT}) where{FT}

Returns the setinitialaux_state! function, which updates the auxiliary state p in place with the initial values corresponding to Y(t=t0) = Y0.

In this case, we also use this method to update the initial values for the spatially varying parameter fields, read in from data files.

source
ClimaLSM.make_set_initial_aux_state(model::CanopyModel)

Returns the setinitialaux_state! function, which updates the auxiliary state p in place with the initial values corresponding to Y(t=t0) = Y0.

In this case, we also use this method to update the initial values for the spatially and temporally varying canopy parameter fields, read in from data files or otherwise prescribed.

source
make_set_initial_aux_state(model::AbstractModel)

Returns the setinitialaux_state! 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.

source
ClimaLSM.prognostic_varsFunction
prognostic_vars(model::PlantHydraulicsModel)

A function which returns the names of the prognostic variables of the PlantHydraulicsModel.

source
ClimaLSM.prognostic_vars(::AbstractCanopyComponent)

Returns the prognostic vars of the canopy component passed in as an argument.

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

source

prognostic_vars(m::AbstractModel)

Returns the prognostic variable symbols for the model in the form of a tuple.

source
prognostic_vars(soil::RichardsModel)

A function which returns the names of the prognostic variables of RichardsModel.

source
prognostic_vars(soil::EnergyHydrology)

A function which returns the names of the prognostic variables of EnergyHydrology.

source
ClimaLSM.prognostic_typesFunction
ClimaLSM.prognostic_types(model::PlantHydraulicsModel{FT}) where {FT}

Defines the prognostic types for the PlantHydraulicsModel.

source
ClimaLSM.prognostic_types(::AbstractCanopyComponent)

Returns the prognostic types of the canopy component passed in as an argument.

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

source

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).

source
prognostic_types(soil::EnergyHydrology{FT}) where {FT}

A function which returns the types of the prognostic variables of EnergyHydrology.

source
ClimaLSM.auxiliary_varsFunction
auxiliary_vars(model::PlantHydraulicsModel)

A function which returns the names of the auxiliary variables of the PlantHydraulicsModel, the transpiration stress factor β (unitless), the water potential ψ (m), the volume fluxcross section fa (1/s), and the volume fluxroot cross section in the roots fa_roots (1/s), where the cross section can be represented by an area index.

source
ClimaLSM.auxiliary_vars(::AbstractCanopyComponent)

Returns the auxiliary types of the canopy component passed in as an argument.

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

source

auxiliary_vars(m::AbstractModel)

Returns the auxiliary variable symbols for the model in the form of a tuple.

source
auxiliary_vars(soil::RichardsModel)

A function which returns the names of the auxiliary variables of RichardsModel.

Note that auxiliary variables are not needed for such a simple model. We could instead compute the conductivity and matric potential within the tendency function explicitly, rather than compute and store them in the auxiliary vector p. We did so in this case as a demonstration.

source
auxiliary_vars(soil::EnergyHydrology)

A function which returns the names of the auxiliary variables of EnergyHydrology.

source
ClimaLSM.auxiliary_typesFunction
ClimaLSM.auxiliary_types(model::PlantHydraulicsModel{FT}) where {FT}

Defines the auxiliary types for the PlantHydraulicsModel.

source
ClimaLSM.auxiliary_types(::AbstractCanopyComponent)

Returns the auxiliary types of the canopy component passed in as an argument.

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

source

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).

source
auxiliary_types(soil::EnergyHydrology{FT}) where {FT}

A function which returns the types of the auxiliary variables of EnergyHydrology.

source
ClimaLSM.initialize_prognosticFunction
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.

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

source
initialize_prognostic(model::AbstractModel, state::Union{ClimaCore.Fields.Field, Vector{FT}})

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.

source
ClimaLSM.initialize_auxiliaryFunction
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.

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

source
initialize_auxiliary(model::AbstractModel, state::Union{ClimaCore.Fields.Field, Vector{FT}})

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.

source
ClimaLSM.initializeFunction
ClimaLSM.initialize(model::BucketModel{FT}) where {FT}

Initializes the variables for the BucketModel.

Note that the BucketModel has prognostic variables that are defined on different subsets of the domain. Because of that, we have to treat them independently. In LSM models which are combinations of standalone component models, this is not needed, and we can use the default initialize. Here, however, we need to do some hardcoding specific to this model.

source
initialize(model::AbstractModel)

Creates the prognostic and auxiliary states structures, but with unset values; constructs and returns the coordinates for the model domain. We may need to consider this default more as we add diverse components and Simulations.

source
ClimaLSM.nameFunction
name(model::AbstractModel)

Returns a symbol of the model component name, e.g. :soil or :vegetation.

source
ClimaLSM.AbstractBCType
AbstractBC

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.

source
ClimaLSM.source!Function
ClimaLSM.source!(dY::ClimaCore.Fields.FieldVector,
                      src::MicrobeProduction,
                      Y::ClimaCore.Fields.FieldVector,
                      p::NamedTuple,
                      params)

A method which extends the ClimaLSM source! function for the case of microbe production of CO2 in soil.

source
 source!(dY::ClimaCore.Fields.FieldVector,
         src::AbstractSource,
         Y::ClimaCore.Fields.FieldVector,
         p::NamedTuple
         )::ClimaCore.Fields.Field

A stub function, which is extended by ClimaLSM.

source
ClimaLSM.source!(dY::ClimaCore.Fields.FieldVector,
                 src::RootExtraction,
                 Y::ClimaCore.Fields.FieldVector,
                 p::NamedTuple
                 model::RichardsModel)

An extension of the ClimaLSM.source! function, which computes source terms for the soil model; this method returns the water loss or gain due to roots when a plant hydraulic prognostic model is included.

source
ClimaLSM.source!(dY::ClimaCore.Fields.FieldVector,
                 src::RootExtraction,
                 Y::ClimaCore.Fields.FieldVector,
                 p::NamedTuple
                 model::EnergyHydrology)

An extension of the ClimaLSM.source! function, which computes source terms for the soil model; this method returns the water and energy loss/gain due to root extraction.

source
 source!(dY::ClimaCore.Fields.FieldVector,
         src::PhaseChange{FT},
         Y::ClimaCore.Fields.FieldVector,
         p::NamedTuple,
         model
         )

Computes the source terms for phase change.

source
ClimaLSM.AbstractBoundaryType
AbstractBoundary

An abstract type to indicate which boundary we are doing calculations for. Currently, we support the top boundary (TopBoundary) and bottom boundary (BottomBoundary).

source
ClimaLSM.TopBoundaryType
TopBoundary{} <: AbstractBoundary{}

A simple object which should be passed into a function to indicate that we are considering the top boundary of the soil.

source
ClimaLSM.BottomBoundaryType
BottomBoundary{} <: AbstractBoundary{}

A simple object which should be passed into a function to indicate that we are considering the bottom boundary of the soil.

source
ClimaLSM.boundary_fluxFunction
ClimaLSM.boundary_flux(
    bc::SoilCO2FluxBC,
    boundary::ClimaLSM.AbstractBoundary,
    Δz::ClimaCore.Fields.Field,
    Y::ClimaCore.Fields.FieldVector,
    p::NamedTuple,
    t::FT,
)::ClimaCore.Fields.Field where {FT}

A method of ClimaLSM.boundary_flux which returns 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.

source
ClimaLSM.boundary_flux(
bc::SoilCO2StateBC,
boundary::ClimaLSM.TopBoundary,
Δz::ClimaCore.Fields.Field,
Y::ClimaCore.Fields.FieldVector,
p::NamedTuple,
t::FT,
)::ClimaCore.Fields.Field where {FT}

A method of ClimaLSM.boundary_flux which returns the soilco2 flux in the case of a prescribed state BC at top of the domain.

source
ClimaLSM.boundary_flux(
    bc::SoilCO2StateBC,
    boundary::ClimaLSM.BottomBoundary,
    Δz::ClimaCore.Fields.Field,
    Y::ClimaCore.Fields.FieldVector,
    p::NamedTuple,
    t::FT,
)::ClimaCore.Fields.Field where {FT}

A method of ClimaLSM.boundary_flux which returns the soilco2 flux in the case of a prescribed state BC at bottom of the domain.

source
boundary_flux(bc::AbstractBC, bound_type::AbstractBoundary, Δz, _...)::ClimaCore.Fields.Field

A function which returns the correct boundary flux given any boundary condition (BC).

source
function ClimaLSM.boundary_flux(
    bc::RunoffBC,
    ::TopBoundary,
    Δz::FT,
    p::NamedTuple,
    t::FT,
    params,
)::ClimaCore.Fields.Field where {FT}

Extension of the ClimaLSM.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).

source
ClimaLSM.boundary_flux(bc::FluxBC,  _...)::ClimaCore.Fields.Field

A method of boundary fluxes which returns the desired flux.

We add a field of zeros in order to convert the bc (float) into a field.

source
ClimaLSM.boundary_flux(bc::RichardsAtmosDrivenFluxBC,
                       boundary::ClimaLSM.AbstractBoundary,
                       model::RichardsModel,
                       Δz::ClimaCore.Fields.Field,
                       Y::ClimaCore.Fields.FieldVector,
                       p::NamedTuple,
                       t,
                       )::ClimaCore.Fields.Field

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.

source
ClimaLSM.boundary_flux(rre_bc::MoistureStateBC,
                       ::ClimaLSM.TopBoundary,
                       model::AbstractSoilModel,
                       Δz::ClimaCore.Fields.Field,
                       Y::ClimaCore.Fields.FieldVector,
                       p::NamedTuple,
                       t,
                       )::ClimaCore.Fields.Field

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.

source
ClimaLSM.boundary_flux(rre_bc::MoistureStateBC,
                       ::ClimaLSM.BottomBoundary,
                       model::AbstractSoilModel,
                       Δz::ClimaCore.Fields.Field,
                       Y::ClimaCore.Fields.FieldVector,
                       p::NamedTuple,
                       t,
                       )::ClimaCore.Fields.Field

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.

source
ClimaLSM.boundary_flux(heat_bc::TemperatureStateBC,
                       ::ClimaLSM.TopBoundary,
                       model::EnergyHydrology,
                       Δz::ClimaCore.Fields.Field,
                       Y::ClimaCore.Fields.FieldVector,
                       p::NamedTuple,
                       t,
                       )::ClimaCore.Fields.Field

A method of boundary fluxes which converts a state boundary condition on temperature at the top of the domain into a flux of energy.

source
ClimaLSM.boundary_flux(heat_bc::TemperatureStateBC,
                       ::ClimaLSM.BottomBoundary,
                       model::EnergyHydrology,
                       Δz::ClimaCore.Fields.Field,
                       Y::ClimaCore.Fields.FieldVector,
                       p::NamedTuple,
                       t,
                       )::ClimaCore.Fields.Field

A method of boundary fluxes which converts a state boundary condition on temperature at the bottom of the domain into a flux of energy.

source
ClimaLSM.boundary_flux(bc::FreeDrainage,
                       boundary::ClimaLSM.BottomBoundary,
                       model::AbstractSoilModel,
                       Δz::ClimaCore.Fields.Field,
                       Y::ClimaCore.Fields.FieldVector,
                       p::NamedTuple,
                       t,
                       )::ClimaCore.Fields.Field

A method of boundary fluxes which enforces free drainage at the bottom of the domain.

source
ClimaLSM.diffusive_fluxFunction
diffusive_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.

source
ClimaLSM.get_ΔzFunction
get_Δz(z::ClimaCore.Fields.Field)

A function to return a tuple containing the distance between the top boundary and its closest center, and the bottom boundary and its closest center, both as Fields.

source
ClimaLSM.make_tendency_jacobianFunction

maketendencyjacobian(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 tendency_jacobian! should be used as Wfact! in ClimaTimeSteppers.jl and SciMLBase.jl.

source
ClimaLSM.make_update_jacobianFunction
make_update_jacobian(model::AbstractModel)

Creates and returns a function which updates 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, no implicit tendency is present, and hence the timestepping is entirely explicit.

source
ClimaLSM.make_update_jacobian(model::RichardsModel)

Creates and returns the update_jacobian! function for RichardsModel.

Using this Jacobian with a backwards Euler timestepper is equivalent to using the modified Picard scheme of Celia et al. (1990).

source
ClimaLSM.∂tendencyBC∂YFunction
∂tendencyBC∂Y(::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.

source
ClimaLSM.∂tendencyBC∂Y(
    model::RichardsModel,
    ::MoistureStateBC,
    boundary::ClimaLSM.TopBoundary,
    Δz,
    Y,
    p,
    t,

)

Computes and returns the derivative of the part of the implicit tendency in the top layer, due to the boundary condition, with respect to the state variable in the top layer.

For a diffusion equation like Richards equation with a single state variable, this is given by ∂T_N∂Y_N = [-∂/∂z(∂F_bc/∂Y_N)]_N, where N indicates the top layer cell index.

source
ClimaLSM.∂tendencyBC∂Y(
    ::AbstractSoilModel,
    ::AbstractSoilBC,
    boundary::ClimaLSM.TopBoundary,
    Δz,
    Y,
    p,
    t,

)

A default method which computes and returns the zero for the derivative of the part of the implicit tendency in the top layer, due to the boundary condition, with respect to the state variable in the top layer.

For a diffusion equation like Richards equation with a single state variable, this is given by ∂T_N∂Y_N = [-∂/∂z(∂F_bc/∂Y_N)]_N, where N indicates the top layer cell index.

If F_bc can be approximated as independent of Y_N, the derivative is zero.

source

Drivers

ClimaLSM.PrescribedAtmosphereType
PrescribedAtmosphere{FT, LP, SP, TA, UA, QA, RA, CA, DT} <: AbstractAtmosphericDrivers{FT}

Container for holding prescribed atmospheric drivers and other information needed for computing turbulent surface fluxes when driving land models in standalone mode.

Since not all models require co2 concentration, the default for that is nothing.

  • liquid_precip: Precipitation (m/s) function of time: positive by definition

  • snow_precip: Snow precipitation (m/s) function of time: positive by definition

  • T: Prescribed atmospheric temperature (function of time) at the reference height (K)

  • u: Prescribed wind speed (function of time) at the reference height (m/s)

  • q: Prescribed specific humidity (function of time) at the reference height (_)

  • P: Prescribed air pressure (function of time) at the reference height (Pa)

  • c_co2: CO2 concentration in atmosphere (mol/mol)

  • ref_time: Reference time - the datetime corresponding to t=0 for the simulation

  • h: Reference height (m), relative to surface elevation

  • gustiness: Minimum wind speed (gustiness; m/s)

source
ClimaLSM.PrescribedRadiativeFluxesType
PrescribedRadiativeFluxes{FT, SW, LW, DT, T, OD} <: AbstractRadiativeDrivers{FT}

Container for the prescribed radiation functions needed to drive land models in standalone mode.

  • SW_d: Downward shortwave radiation function of time (W/m^2): positive indicates towards surface

  • LW_d: Downward longwave radiation function of time (W/m^2): positive indicates towards surface

  • ref_time: Reference time - the datetime corresponding to t=0 for the simulation

  • θs: Sun zenith angle, in radians

  • orbital_data: Orbital Data for Insolation.jl

source
ClimaLSM.surface_fluxes_at_a_pointFunction
surface_fluxes_at_a_point(T_sfc::FT,
                          q_sfc::FT,
                          ρ_sfc::FT,
                          β_sfc::FT,
                          h_sfc::FT,
                          t::FT,
                          parameters,
                          atmos::PA,
                          ) where {FT <: AbstractFloat, PA <: PrescribedAtmosphere{FT}}

Computes turbulent surface fluxes at a point on a surface given (1) the surface temperature, specific humidity, and air density, (2) the time at which the fluxes are needed, (3) a factor βsfc which scales the evaporation from the potential rate, (4) the parameter set for the model, which must have fields `earthparamset, and roughness lengthsz0m, z_0b. (5) the prescribed atmospheric state, stored inatmos`.

This returns an energy flux and a liquid water volume flux, stored in a tuple with self explanatory keys.

source
Missing docstring.

Missing docstring for ClimaLSM.radiative_fluxes_at_a_point. Check Documenter's build log for details.

ClimaLSM.construct_atmos_tsFunction
construct_atmos_ts(
    atmos::PrescribedAtmosphere,
    t::FT,
    thermo_params,
) where {FT}

A helper function which constructs a Clima.Thermodynamics thermodynamic state given a PrescribedAtmosphere, a time at which the state is needed, and a set of Clima.Thermodynamics parameters thermo_params.

source
ClimaLSM.surface_air_densityFunction
ClimaLSM.surface_air_density(model::BucketModel, Y, p)

a helper function which computes and returns the surface air density for the bucket model.

source
ClimaLSM.surface_air_density(model::CanopyModel, Y, p)

A helper function which computes and returns the surface air density for the canopy model.

source
surface_air_density(
                    atmos::AbstractAtmosphericDrivers,
                    model::AbstractModel,
                    Y,
                    p,
                    t,
                    T_sfc,
                    )

A helper function which returns the surface air density for a given model, needed because different models compute and store surface air density in different ways and places.

We additionally include the atmos type as an argument because the surface air density computation may change between a coupled simulation and a prescibed atmos simulation.

Extending this function for your model is only necessary if you need to compute surface fluxes and radiative fluxes at the surface using the functions in this file.

source
ClimaLSM.surface_air_density(
    atmos::PrescribedAtmosphere{FT},
    model::EnergyHydrology{FT},
    Y,
    p,
    t,
    T_sfc
) where {FT}

Returns the surface air density field of the EnergyHydrology soil model for the PrescribedAtmosphere case.

This assumes the ideal gas law and hydrostatic balance to estimate the air density at the surface from the values of surface temperature and the atmospheric thermodynamic state, because the surface air density is not a prognostic variable of the soil model.

source
ClimaLSM.snow_precipitationFunction
snow_precipitation(atmos::PrescribedAtmosphere, p, t)

Returns the precipitation in snow (m of liquid water/s) at the surface.

source