Soil Models

Soil Models

ClimaLand.Soil.AbstractSoilModelType
AbstractSoilModel{FT} <: ClimaLand.AbstractImExModel{FT}

The abstract type for all soil models.

Currently, we only have plans to support a RichardsModel, simulating the flow of liquid water through soil via the Richardson-Richards equation, and a fully integrated soil heat and water model, with phase change.

source
ClimaLand.Soil.RichardsModelType
RichardsModel

A model for simulating the flow of water in a porous medium by solving the Richardson-Richards Equation.

A variety of boundary condition types are supported, including FluxBC, RichardsAtmosDrivenFluxBC, MoistureStateBC, and FreeDrainage (only for the bottom of the domain).

If you wish to simulate soil hydrology under the context of a prescribed precipitation volume flux (m/s) as a function of time, the RichardsAtmosDrivenFluxBC type should be chosen. Please see the documentation for more details.

  • parameters: the parameter set

  • domain: the soil domain, using ClimaCore.Domains

  • boundary_conditions: the boundary conditions, of type AbstractSoilBoundaryConditions

  • sources: A tuple of sources, each of type AbstractSoilSource

  • lateral_flow: A boolean flag which, when false, turns off the horizontal flow of water

source
ClimaLand.Soil.EnergyHydrologyType
EnergyHydrology <: AbstractSoilModel

A model for simulating the flow of water and heat in a porous medium by solving the Richardson-Richards equation and the heat equation, including terms for phase change.

A variety of boundary condition types are supported, including FluxBC, MoistureStateBC/TemperatureStateBC, FreeDrainage (only for the bottom of the domain), and an AtmosDrivenFluxBC (under which radiative fluxes and turbulent surface fluxes are computed and used as boundary conditions). Please see the documentation for this boundary condition type for more details.

  • parameters: The parameter sets

  • domain: the soil domain, using ClimaCore.Domains

  • boundary_conditions: the boundary conditions for RRE and heat, of type AbstractSoilBoundaryConditions

  • sources: A tuple of sources, each of type AbstractSoilSource

  • lateral_flow: A boolean flag which, when false, turns off the horizontal flow of water and heat

source

Soil Parameter Structs

ClimaLand.Soil.RichardsParametersType
RichardsParameters{F <: Union{<: AbstractFloat, ClimaCore.Fields.Field}, C <: AbstractSoilHydrologyClosure}

A struct for storing parameters of the RichardModel.

  • ν: The porosity of the soil (m^3/m^3)

  • hydrology_cm: The hydrology closure model: vanGenuchten or BrooksCorey

  • K_sat: The saturated hydraulic conductivity (m/s)

  • S_s: The specific storativity (1/m)

  • θ_r: The residual water fraction (m^3/m^3

source
ClimaLand.Soil.EnergyHydrologyParametersType
EnergyHydrologyParameters{FT <: AbstractFloat, 
                          F <: Union{<:AbstractFloat,
                                     ClimaCore.Fields.Field},
                          C,
                          PSE,}

A parameter structure for the integrated soil water and energy equation system.

  • κ_dry: The dry soil thermal conductivity, W/m/K

  • κ_sat_frozen: The saturated thermal conductivity of frozen soil, W/m/K

  • κ_sat_unfrozen: The saturated thermal conductivity of unfrozen soil, W/m/K

  • ρc_ds: The volumetric heat capacity of dry soil, J/m^3/K (per volume dry soil, not per volume soil solids)

  • ν: The porosity of the soil (m^3/m^3)

  • ν_ss_om: The volumetric fraction of the soil solids in organic matter (m^3/m^3)

  • ν_ss_quartz: The volumetric fraction of the soil solids in quartz (m^3/m^3)

  • ν_ss_gravel: The volumetric fraction of the soil solids in gravel (m^3/m^3)

  • α: The parameter α used in computing Kersten number, unitless

  • β: The parameter β used in computing Kersten number, unitless

  • hydrology_cm: The soil hydrology closure model: van Genuchten or Brooks and Corey

  • K_sat: The saturated hydraulic conductivity (m/s)

  • S_s: The specific storativity (1/m)

  • θ_r: The residual water fraction (m^3/m^3

  • Ω: Ice impedance factor for the hydraulic conductivity

  • γ: Coefficient of viscosity factor for the hydraulic conductivity

  • γT_ref: Reference temperature for the viscosity factor

  • PAR_albedo: Soil PAR Albedo

  • NIR_albedo: Soil NIR Albedo

  • emissivity: Soil Emissivity

  • z_0m: Roughness length for momentum

  • z_0b: Roughness length for scalars

  • d_ds: Maximum dry soil layer thickness under evaporation (m)

  • earth_param_set: Physical constants and clima-wide parameters

source

Soil Hydrology Parameterizations

ClimaLand.Soil.volumetric_liquid_fractionFunction
volumetric_liquid_fraction(ϑ_l::FT, ν_eff::FT, θ_r::FT) where {FT}

A pointwise function returning the volumetric liquid fraction given the augmented liquid fraction and the effective porosity.

source
ClimaLand.Soil.pressure_headFunction
pressure_head(
    cm::vanGenuchten{FT},
    θ_r::FT,
    ϑ_l::FT,
    ν_eff::FT,
    S_s::FT,
) where {FT}

A point-wise function returning the pressure head in variably saturated soil, using the van Genuchten matric potential if the soil is not saturated, and an approximation of the positive pressure in the soil if the soil is saturated.

source
pressure_head(
    cm::BrooksCorey{FT},
    θ_r::FT,
    ϑ_l::FT,
    ν_eff::FT,
    S_s::FT,
) where {FT}

A point-wise function returning the pressure head in variably saturated soil, using the Brooks and Corey matric potential if the soil is not saturated, and an approximation of the positive pressure in the soil if the soil is saturated.

source
ClimaLand.Soil.hydraulic_conductivityFunction
 hydraulic_conductivity(cm::vanGenuchten{FT}, K_sat::FT, S::FT) where {FT}

A point-wise function returning the hydraulic conductivity, using the van Genuchten formulation.

source
 hydraulic_conductivity(cm::BrooksCorey{FT}, K_sat::FT, S::FT) where {FT}

A point-wise function returning the hydraulic conductivity, using the Brooks and Corey formulation.

source
ClimaLand.Soil.impedance_factorFunction
impedance_factor(
    f_i::FT,
    Ω::FT
) where {FT}

Returns the multiplicative factor reducing conductivity when a fraction of ice f_i is present.

Only for use with the EnergyHydrology model.

source
ClimaLand.Soil.viscosity_factorFunction
viscosity_factor(
    T::FT,
    γ::FT,
    γT_ref::FT,
) where {FT}

Returns the multiplicative factor which accounts for the temperature dependence of the conductivity.

Only for use with the EnergyHydrology model.

source
ClimaLand.Soil.matric_potentialFunction
 matric_potential(cm::vanGenuchten{FT}, S::FT) where {FT}

A point-wise function returning the matric potential, using the van Genuchten formulation.

source
 matric_potential(cm::BrooksCorey{FT}, S::FT) where {FT}

A point-wise function returning the matric potential, using the Brooks and Corey formulation.

source
ClimaLand.Soil.dψdϑFunction

dψdϑ(cm::vanGenuchten{FT}, ϑ, ν, θr, Ss)

Computes and returns the derivative of the pressure head with respect to ϑ for the van Genuchten formulation.

source

dψdϑ(cm::BrooksCorey{FT}, ϑ, ν, θr, Ss)

Computes and returns the derivative of the pressure head with respect to ϑ for the Brooks and Corey formulation.

source
ClimaLand.Soil.inverse_matric_potentialFunction
 inverse_matric_potential(cm::vanGenuchten{FT}, ψ::FT) where {FT}

A point-wise function returning the effective saturation, given the matric potential, using the van Genuchten formulation.

source
 inverse_matric_potential(cm::BrooksCorey{FT}, ψ::FT) where {FT}

A point-wise function returning the effective saturation, given the matric potential, using the Brooks and Corey formulation.

source
ClimaLand.Soil.AbstractSoilHydrologyClosureType
AbstractSoilHydrologyClosure{FT <: AbstractFloat}

The abstract type of soil hydrology closure, of which vanGenuchten{FT} and BrooksCorey{FT} are the two supported concrete types.

To add a new parameterization, methods are required for:

  • matric_potential,
  • inversematricpotential,
  • pressure_head,
  • dψdϑ,
  • hydraulic_conductivity.
source
ClimaLand.Soil.vanGenuchtenType
vanGenuchten{FT} <: AbstractSoilHydrologyClosure{FT}

The van Genuchten soil hydrology closure, chosen when the hydraulic conductivity and matric potential are modeled using the van Genuchten parameterization (van Genuchten 1980; see also Table 8.2 of G. Bonan 2019).

  • α: The inverse of the air entry potential (1/m)

  • n: The van Genuchten pore-size distribution index (unitless)

  • m: The van Genuchten parameter m = 1 - 1/n (unitless)

  • S_c: A derived parameter: the critical saturation at which capillary flow no longer replenishes the surface

source
ClimaLand.Soil.BrooksCoreyType

BrooksCorey{FT} <: AbstractSoilHydrologyClosure{FT}

The Brooks and Corey soil hydrology closure, chosen when the hydraulic conductivity and matric potential are modeled using the Brooks and Corey parameterization (Brooks and Corey, 1964, 1966; see also Table 8.2 of G. Bonan 2019).

  • c: The pore-size distribution index (unitless)

  • ψb: The air entry matric potential, when S=1 (m)

  • S_c: A derived parameter: the critical saturation at which capillary flow no longer replenishes the surface

source

Soil Heat Parameterizations

ClimaLand.Soil.κ_solidFunction
κ_solid(ν_ss_om::FT,
        ν_ss_quartz::FT,
        κ_om::FT,
        κ_quartz::FT,
        κ_minerals::FT) where {FT}

Computes the thermal conductivity of the solid material in soil. The _ss_ subscript denotes that the volumetric fractions of the soil components are referred to the soil solid components, not including the pore space.

source
ClimaLand.Soil.κ_sat_frozenFunction
function κ_sat_frozen(
    κ_solid::FT,
    ν::FT,
    κ_ice::FT
) where {FT}

Computes the thermal conductivity for saturated frozen soil.

source
ClimaLand.Soil.κ_sat_unfrozenFunction
function κ_sat_unfrozen(
    κ_solid::FT,
    ν::FT,
    κ_l::FT
) where {FT}

Computes the thermal conductivity for saturated unfrozen soil.

source
ClimaLand.Soil.κ_satFunction
κ_sat(
    θ_l::FT,
    θ_i::FT,
    κ_sat_unfrozen::FT,
    κ_sat_frozen::FT
) where {FT}

Compute the expression for saturated thermal conductivity of soil matrix.

source
ClimaLand.Soil.κ_dryFunction
function κ_dry(ρp::FT,
               ν::FT,
               κ_solid::FT,
               κ_air::FT;
               a = FT(0.053)) where {FT}

Computes the thermal conductivity of dry soil according to the model of Balland and Arp.

source
ClimaLand.Soil.kersten_numberFunction
kersten_number(
    θ_i::FT,
    S_r::FT,
    parameters::EnergyHydrologyParameters{FT},
) where {FT}

Compute the expression for the Kersten number, using the Balland and Arp model.

source
ClimaLand.Soil.relative_saturationFunction
relative_saturation(
        θ_l::FT,
        θ_i::FT,
        ν::FT
) where {FT}

Compute the expression for relative saturation. This is referred to as θ_sat in Balland and Arp's paper.

source
ClimaLand.Soil.volumetric_internal_energyFunction
volumetric_internal_energy(θ_i::FT, ρc_s::FT, T::FT,
                             parameters::EnergyHydrologyParameters{FT}) where {FT}

A pointwise function for computing the volumetric internal energy of the soil, given the volumetric ice content, volumetric heat capacity, and temperature.

source
ClimaLand.Soil.volumetric_internal_energy_liqFunction
volumetric_internal_energy_liq(T::FT, parameters::EnergyHydrologyParameters{FT}) where {FT}

A pointwise function for computing the volumetric internal energy of the liquid water in the soil, given the temperature T.

source
ClimaLand.Soil.temperature_from_ρe_intFunction
temperature_from_ρe_int(ρe_int::FT, θ_i::FT, ρc_s::FT
                        parameters::EnergyHydrologyParameters{FT}) where {FT}

A pointwise function for computing the temperature from the volumetric internal energy, volumetric ice content, and volumetric heat capacity of the soil.

source
ClimaLand.Soil.phase_change_sourceFunction
phase_change_source(
    θ_l::FT,
    θ_i::FT,
    T::FT,
    τ::FT,
    params::EnergyHydrologyParameters{FT},
) where {FT}

Returns the source term (1/s) used for converting liquid water and ice into each other during phase changes. Note that there are unitless prefactors multiplying this term in the equations.

Note that these equations match what is in Dall'Amico (for θstar, ψ(T), ψw0). We should double check them in the case where we have ϑl > θl, but they should be very close to the form we want regardless.

source
ClimaLand.Soil.thermal_timeFunction
thermal_time(ρc::FT, Δz::FT, κ::FT) where {FT}

Returns the thermal timescale for temperature differences across a typical thickness Δz to equilibrate.

source

Soil Surface Parameterizations

Missing docstring.

Missing docstring for ClimaLand.soil.soil_resistance. Check Documenter's build log for details.

ClimaLand.Soil.dry_soil_layer_thicknessFunction
dry_soil_layer_thickness(S_l_sfc::FT, S_c::FT, d_ds::FT) where {FT}

Returns the maximum dry soil layer thickness that can develop under evaporation; this is used when computing the soil resistance to evaporation according to Swenson et al (2012).

source
ClimaLand.Soil.soil_tortuosityFunction
soil_tortuosity(θ_l::FT, θ_i::FT, ν::FT) where {FT}

Computes the tortuosity of water vapor in a porous medium, as a function of porosity ν and the volumetric liquid water and ice contents, θ_l and θ_i.

See Equation (1) of : Shokri, N., P. Lehmann, and D. Or (2008), Effects of hydrophobic layers on evaporation from porous media, Geophys. Res. Lett., 35, L19407, doi:10.1029/ 2008GL035230.

source

Soil Runoff Types and Methods

ClimaLand.Soil.Runoff.TOPMODELRunoffType
TOPMODELRunoff{FT <: AbstractFloat, F <: ClimaCore.Fields.Field} <: AbstractRunoffModel

The TOPMODEL surface runoff parameterization, which is affects the surface boundary condition of the soil model.

The runoff flux is given by Equation 8 of with fsat given by Equation (11), of Niu et al. (2005), "A simple TOPMODEL-based runoff parameterization (SIMTOP) for use in global climate models".

  • f_over: A calibrated parameter defining how subsurface runoff decays with depth to water table (1/m ; calibrated)

  • f_max: The maximum saturated fraction of a grid cell, computed from the topographic index CDF per grid cell.

  • subsurface_source: The subsurface source term corresponding to this implementation of TOPMODEL.

source
ClimaLand.Soil.Runoff.TOPMODELSubsurfaceRunoffType
TOPMODELSubsurfaceRunoff{FT} <: AbstractSoilSource{FT}

The TOPMODEL subsurface runoff parameterization, which is implemented as a sink term in the soil equations.

The runoff flux is given by Equation 12 of Niu et al. (2005), "A simple TOPMODEL-based runoff parameterization (SIMTOP) for use in global climate models".

  • R_sb: The subsurface runoff flux (m/s) when the depth to the water table = 1/f_over; calibrated

  • f_over: A calibrated parameter defining how subsurface runoff decays with depth to water table (1/m ; calibrated)

source
ClimaLand.Soil.Runoff.update_runoff!Function
update_runoff!(p, runoff::NoRunoff, _...)

Updates the runoff variables in the cache p.soil in place in the case of NoRunoff: sets infiltration = precipitation.

source
update_runoff!(p, runoff::TOPMODELRunoff, Y,t, model::AbstractSoilModel)

Updates the runoff model variables in place in p.soil for the TOPMODELRunoff parameterization: p.soil.Rs p.soil.Rss p.soil.h∇ p.soil.infiltration

source

Soil BC Methods and Types

ClimaLand.Soil.MoistureStateBCType

MoistureStateBC <: AbstractWaterBC

A simple concrete type of boundary condition, which enforces a state boundary condition ϑ_l = f(p,t) at either the top or bottom of the domain.

source
ClimaLand.Soil.HeatFluxBCType

HeatFluxBC <: AbstractHeatBC

A simple concrete type of boundary condition, which enforces a normal flux boundary condition f(p,t) at either the top or bottom of the domain.

source
ClimaLand.Soil.WaterFluxBCType

WaterFluxBC <: AbstractWaterBC

A simple concrete type of boundary condition, which enforces a normal flux boundary condition f(p,t) at either the top or bottom of the domain.

source
ClimaLand.Soil.TemperatureStateBCType

TemperatureStateBC <: AbstractHeatBC

A simple concrete type of boundary condition, which enforces a state boundary condition T = f(p,t) at either the top or bottom of the domain.

source
ClimaLand.Soil.FreeDrainageType
FreeDrainage <: AbstractWaterBC

A concrete type of soil boundary condition, for use at the BottomBoundary only, where the flux is set to be F = -K∇h = -K.

source
ClimaLand.Soil.RichardsAtmosDrivenFluxBCType

RichardsAtmosDrivenFluxBC{F <: PrescribedPrecipitation, R <: AbstractRunoffModel} <: AbstractWaterBC

A concrete type of boundary condition intended only for use with the RichardsModel, which uses a prescribed precipitation rate (m/s) to compute the infiltration into the soil.

A runoff model is used to simulate surface and subsurface runoff and this is accounted for when setting boundary conditions. In order to run the simulation without runoff, choose runoff = NoRunoff() - this is also the default.

If you wish to simulate precipitation and runoff in the full EnergyHydrology model, you must use the AtmosDrivenFluxBC type.

  • precip: The prescribed liquid water precipitation rate f(t) (m/s); Negative by convention.

  • runoff: The runoff model. The default is no runoff.

source
ClimaLand.Soil.AtmosDrivenFluxBCType
AtmosDrivenFluxBC{
    A <: AbstractAtmosphericDrivers,
    B <: AbstractRadiativeDrivers,
    R <: AbstractRunoffModel
} <: AbstractEnergyHydrologyBC

A concrete type of soil boundary condition for use at the top of the domain. This holds the conditions for the atmosphere AbstractAtmosphericDrivers, for the radiation state AbstractRadiativeDrivers. This is only supported for the EnergyHydrology model.

This choice indicates the Monin-Obukhov Surface Theory will be used to compute the sensible and latent heat fluxes, as well as evaporation, and that the net radiation and precipitation will also be computed. The net energy and water fluxes are used as boundary conditions.

A runoff model is used to simulate surface and subsurface runoff and this is accounted for when setting boundary conditions. The default is to have no runoff accounted for.

  • atmos: The atmospheric conditions driving the model

  • radiation: The radiative fluxes driving the model

  • runoff: The runoff model. The default is no runoff.

source
ClimaLand.Soil.WaterHeatBCType
WaterHeatBC{W <: AbstractWaterBC, H <: AbstractHeatBC} <:
   AbstractEnergyHydrologyBC

A general struct used to store the boundary conditions for Richards and the soil heat equations separately; useful when the boundary conditions for each component are independent of each other.

source
ClimaLand.Soil.soil_boundary_fluxes!Function
soil_boundary_fluxes!(bc::WaterHeatBC, boundary::TopBoundary, model, Δz, Y, p, t)

updates the boundary fluxes for ϑl and ρeint.

source
soil_boundary_fluxes!(
    bc::AtmosDrivenFluxBC{
        <:PrescribedAtmosphere,
        <:PrescribedRadiativeFluxes,
    },
    boundary::ClimaLand.TopBoundary,
    model::EnergyHydrology,
    Δz,
    Y,
    p,
    t,
)

Returns the net volumetric water flux (m/s) and net energy flux (W/m^2) for the soil EnergyHydrology model at the top of the soil domain.

If you wish to compute surface fluxes taking into account the presence of a canopy, snow, etc, as in a land surface model, this is not the correct method to be using.

This function calls the turbulent_fluxes and net_radiation functions, which use the soil surface conditions as well as the atmos and radiation conditions in order to compute the surface fluxes using Monin Obukhov Surface Theory.

source
soil_boundary_fluxes!(
    bc::AtmosDrivenFluxBC{<:PrescribedAtmosphere, <:CanopyRadiativeFluxes},
    boundary::ClimaLand.TopBoundary,
    soil::EnergyHydrology{FT},
    Δz,
    Y,
    p,
    t,
) where {FT}

A method of ClimaLand.Soil.soil_boundary_fluxes! which is used for integrated land surface models; this computes and returns the net energy and water flux at the surface of the soil for use as boundary conditions.

source

Soil Source Types

ClimaLand.Soil.AbstractSoilSourceType
AbstractSoilSource{FT} <:  ClimaLand.AbstractSource{FT}

An abstract type for types of source terms for the soil equations.

In standalone mode, the only supported source type is freezing and thawing. ClimaLand.jl creates additional sources to include as necessary e.g. root extraction (not available in stand alone mode).

source
Missing docstring.

Missing docstring for ClimaLand.Soil.RootExtraction. Check Documenter's build log for details.

Soil Jacobian Structures

ClimaLand.Soil.RichardsTridiagonalWType
RichardsTridiagonalW{R, J, W, T} <: ClimaLand.AbstractTridiagonalW

A struct containing the necessary information for constructing a tridiagonal Jacobian matrix (W) for solving Richards equation, treating only the vertical diffusion term implicitly.

Note that the diagonal, upper diagonal, and lower diagonal entry values are stored in this struct and updated in place.

  • dtγ_ref: Reference to dtγ, which is specified by the ODE solver

  • ∂ϑₜ∂ϑ: Diagonal entries of the Jacobian stored as a ClimaCore.Fields.Field

  • W_column_arrays: Array of tridiagonal matrices containing W for each column

  • temp1: An allocated cache used to evaluate ldiv!

  • temp2: An allocated cache used to evaluate ldiv!

  • transform: A flag indicating whether this struct is used to compute Wfact_t or Wfact

  • ones_face_space: A pre-allocated cache storing ones on the face space

source