Soil Models
Soil Models
ClimaLSM.Soil.AbstractSoilModel
— TypeAbstractSoilModel{FT} <: ClimaLSM.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.
ClimaLSM.Soil.RichardsModel
— TypeRichardsModel
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 setdomain
: the soil domain, using ClimaCore.Domainsboundary_conditions
: the boundary conditions, of type AbstractSoilBoundaryConditionssources
: A tuple of sources, each of type AbstractSoilSourcelateral_flow
: A boolean flag which, when false, turns off the horizontal flow of water
ClimaLSM.Soil.EnergyHydrology
— TypeEnergyHydrology <: 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 setsdomain
: the soil domain, using ClimaCore.Domainsboundary_conditions
: the boundary conditions for RRE and heat, of type AbstractSoilBoundaryConditionssources
: A tuple of sources, each of type AbstractSoilSourcelateral_flow
: A boolean flag which, when false, turns off the horizontal flow of water and heat
Soil Parameter Structs
ClimaLSM.Soil.RichardsParameters
— TypeRichardsParameters{FT <: AbstractFloat}
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 BrooksCoreyK_sat
: The saturated hydraulic conductivity (m/s)S_s
: The specific storativity (1/m)θ_r
: The residual water fraction (m^3/m^3
ClimaLSM.Soil.EnergyHydrologyParameters
— TypeEnergyHydrologyParameters{FT <: AbstractFloat}
A parameter structure for the integrated soil water and energy equation system. In this simplest form, we assume the conductivity and volumetric heat capacity of the soil are constant.
κ_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ν
: 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, unitlesshydrology_cm
: The soil hydrology closure model: van Genucthen or Brooks and CoreyK_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 factorPAR_albedo
: Soil PAR AlbedoNIR_albedo
: Soil NIR Albedoemissivity
: Soil Emissivityz_0m
: Roughness length for momentumz_0b
: Roughness length for scalarsd_ds
: Maximum dry soil layer thickness under evaporation (m)earth_param_set
: Physical constants and clima-wide parameters
Soil Hydrology Parameterizations
ClimaLSM.Soil.volumetric_liquid_fraction
— Functionvolumetric_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.
ClimaLSM.Soil.pressure_head
— Functionpressure_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.
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.
ClimaLSM.Soil.hydraulic_conductivity
— Function 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.
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.
ClimaLSM.Soil.impedance_factor
— Functionimpedance_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.
ClimaLSM.Soil.viscosity_factor
— Functionviscosity_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.
ClimaLSM.Soil.effective_saturation
— Functioneffective_saturation(porosity::FT, ϑ_l::FT, θr::FT) where {FT}
A point-wise function computing the effective saturation.
ClimaLSM.Soil.matric_potential
— Function matric_potential(cm::vanGenuchten{FT}, S::FT) where {FT}
A point-wise function returning the matric potential, using the van Genuchten formulation.
matric_potential(cm::BrooksCorey{FT}, S::FT) where {FT}
A point-wise function returning the matric potential, using the Brooks and Corey formulation.
ClimaLSM.Soil.dψdϑ
— Functiondψdϑ(cm::vanGenuchten{FT}, ϑ, ν, θr, Ss)
Computes and returns the derivative of the pressure head with respect to ϑ for the van Genuchten formulation.
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.
ClimaLSM.Soil.inverse_matric_potential
— Function 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.
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.
ClimaLSM.Soil.AbstractSoilHydrologyClosure
— TypeAbstractSoilHydrologyClosure{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.
ClimaLSM.Soil.vanGenuchten
— TypevanGenuchten{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
ClimaLSM.Soil.BrooksCorey
— TypeBrooksCorey{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
Soil Heat Parameterizations
ClimaLSM.Soil.volumetric_heat_capacity
— Functionvolumetric_heat_capacity(
θ_l::FT,
θ_i::FT,
parameters::EnergyHydrologyParameters{FT},
) where {FT}
Compute the expression for volumetric heat capacity.
ClimaLSM.Soil.κ_solid
— Functionκ_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.
ClimaLSM.Soil.κ_sat_frozen
— Functionfunction κ_sat_frozen(
κ_solid::FT,
ν::FT,
κ_ice::FT
) where {FT}
Computes the thermal conductivity for saturated frozen soil.
ClimaLSM.Soil.κ_sat_unfrozen
— Functionfunction κ_sat_unfrozen(
κ_solid::FT,
ν::FT,
κ_l::FT
) where {FT}
Computes the thermal conductivity for saturated unfrozen soil.
ClimaLSM.Soil.κ_sat
— Functionκ_sat(
θ_l::FT,
θ_i::FT,
κ_sat_unfrozen::FT,
κ_sat_frozen::FT
) where {FT}
Compute the expression for saturated thermal conductivity of soil matrix.
ClimaLSM.Soil.κ_dry
— Functionfunction κ_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.
ClimaLSM.Soil.kersten_number
— Functionkersten_number(
θ_i::FT,
S_r::FT,
parameters::EnergyHydrologyParameters{FT},
) where {FT}
Compute the expression for the Kersten number, using the Balland and Arp model.
ClimaLSM.Soil.relative_saturation
— Functionrelative_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.
ClimaLSM.Soil.volumetric_internal_energy
— Functionvolumetric_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.
ClimaLSM.Soil.volumetric_internal_energy_liq
— Functionvolumetric_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.
ClimaLSM.Soil.temperature_from_ρe_int
— Functiontemperature_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.
ClimaLSM.Soil.thermal_conductivity
— Functionthermal_conductivity(
κ_dry::FT,
K_e::FT,
κ_sat::FT
) where {FT}
Compute the expression for thermal conductivity of soil matrix.
ClimaLSM.Soil.phase_change_source
— Functionphase_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.
ClimaLSM.Soil.thermal_time
— Functionthermal_time(ρc::FT, Δz::FT, κ::FT) where {FT}
Returns the thermal timescale for temperature differences across a typical thickness Δz to equilibrate.
Soil Surface Parameterizations
Missing docstring for ClimaLSM.soil.soil_resistance
. Check Documenter's build log for details.
ClimaLSM.Soil.dry_soil_layer_thickness
— Functiondry_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).
ClimaLSM.Soil.soil_tortuosity
— Functionsoil_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.
Soil Runoff Types and Methods
ClimaLSM.Soil.NoRunoff
— TypeNoRunoff <: AbstractRunoffModel
A concrete type of soil runoff model; the default choice which does not include the effects of runoff.
ClimaLSM.Soil.subsurface_runoff_source
— Functionsubsurface_runoff_source(runoff::AbstractRunoffModel)::Union{Nothing, AbstractSoilSource}
A function which returns the soil source for the runoff model runoff
; the default returns nothing in which case no source is added.
ClimaLSM.Soil.soil_surface_infiltration
— Functionsoil_surface_infiltration(::NoRunoff, net_water_flux, _...)
A function which computes the infiltration into the soil for the default of NoRunoff
.
If net_water_flux = P+E
, where P
is the precipitation and E
is the evaporation (both negative if towards the soil), this returns P+E
as the water boundary flux for the soil.
Soil BC Methods and Types
ClimaLSM.Soil.AbstractSoilBC
— TypeAbstractSoilBC <: ClimaLSM. AbstractBC
An abstract type for soil-specific types of boundary conditions, like free drainage.
ClimaLSM.Soil.MoistureStateBC
— TypeMoistureStateBC <: AbstractSoilBC
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.
ClimaLSM.Soil.FluxBC
— TypeFluxBC <: AbstractSoilBC
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.
ClimaLSM.Soil.TemperatureStateBC
— TypeTemperatureStateBC <: AbstractSoilBC
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.
ClimaLSM.Soil.FreeDrainage
— TypeFreeDrainage <: AbstractSoilBC
A concrete type of soil boundary condition, for use at the BottomBoundary only, where the flux is set to be F = -K∇h = -K
.
ClimaLSM.Soil.RichardsAtmosDrivenFluxBC
— TypeRichardsAtmosDrivenFluxBC{R <: AbstractRunoffModel} <: AbstractSoilBC
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 preciptation 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.
ClimaLSM.Soil.AtmosDrivenFluxBC
— TypeAtmosDrivenFluxBC{
A <: AbstractAtmosphericDrivers,
B <: AbstractRadiativeDrivers,
R <: AbstractRunoffModel
} <: AbstractSoilBC
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 modelradiation
: The radiative fluxes driving the modelrunoff
: The runoff model. The default is no runoff.
Soil Source Types
ClimaLSM.Soil.AbstractSoilSource
— TypeAbstractSoilSource{FT} <: ClimaLSM.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. ClimaLSM.jl creates additional sources to include as necessary e.g. root extraction (not available in stand alone mode).
ClimaLSM.Soil.PhaseChange
— TypePhaseChange{FT} <: AbstractSoilSource{FT}
PhaseChange source type.
Missing docstring for ClimaLSM.Soil.RootExtraction
. Check Documenter's build log for details.
Soil Jacobian Structures
ClimaLSM.Soil.RichardsTridiagonalW
— TypeRichardsTridiagonalW{R, J, W, T} <: ClimaLSM.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.FieldW_column_arrays
: Array of tridiagonal matrices containing W for each columntemp1
: 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 Wfactones_face_space
: A pre-allocated cache storing ones on the face space