EvaporationModel documentation
Config
EvaporationModel.VegetationParameters
— TypeVegetationParameters(;...)
Default vegetation parameters for the Jarvis-Stewart model. The fields are:
vegtype
: The type of vegetation (e.g.,Crops
,ShortGrass
, etc.)r_smin
: Minimum stomatal resistance [s/m]g_d
: Coefficient relating vapour pressure deficit to stomatal resistance [Pa⁻¹]
Based on get_default_value
function, values of r_smin
and g_d
are set based on vegtype
. Default values can be overridden by passing them as keyword arguments.
Examples
using EvaporationModel
# Get default parameters
params_default = VegetationParameters(vegtype=Crops())
# Adapt a default value
params = VegetationParameters(vegtype=Crops(), r_smin=300.0)
params.r_smin == 300.0
# output
true
EvaporationModel.VegetationType
— TypeVegetationType
Abstract type defined to represent different types of vegetation.
Below example provided on how to check all the subtypes (which are structs):
subtypes(VegetationType)
EvaporationModel.get_default_value
— Methodget_default_value(vegtype::VegetationType)
Assign default values for minimal stomatal resistance (r_smin
) and coefficient relating vapour pressure deficit to stomatal resistance (g_d
) based on vegetation type.
The values come from Table 8.1 of the IFS Cy49r1 documentation Part IV: Phyiscal processes. Not part of public API, only used in VegetationParameters
.
Evaporation
EvaporationModel.penman_monteith
— Methodpenman_monteith(Temp, p, VPD, A, r_a, r_s; kwargs...)
Compute evaporation (ET) and latent heat flux (LE) using the Penman-Monteith equation.
Arguments
Temp
: Temperature [K]p
: Pressure [Pa]VPD
: Vapor pressure deficit [Pa]A
: Available energy ($R_n -G$) [W/m²]r_a
: Aerodynamic resistance [s/m]r_s
: Surface resistance [s/m]kwargs
: Additional keyword arguments passed to theBigleaf.potential_ET
function.
Returns
ET
: Potential evapotranspiration [kg/(m² * s)]λE
: Latent heat flux [W/m²]
See also
Bigleaf.potential_ET
for more details on the Penman-Monteith equation.
Examples
using EvaporationModel
Temp = 273.15 + 30 # K
p = 101325.0 # Pa
VPD = 30.0 * 100 # Pa
A = 500.0 # W/m²
r_a = 100.0 # s/m
r_s = 150.0 # s/m
ET, λE = penman_monteith(Temp, p, VPD, A, r_a, r_s)
round(λE; digits = 2) ≈ 380.68
# output
true
EvaporationModel.total_evaporation
— Methodtotal_evaporation(T_a, p_a, VPD, A, A_c, A_s, r_aa, r_ac, r_as, r_sc, r_ss, f_wet)
Compute total evaporation (ET) / latent heat flux (λE) using a multi-source model accounting for bare soil evaporation, transpiration and interception
Arguments
T_a
: Air temperature (at $z_a$) [K]p_a
: Air pressure (at $z_a$)[Pa]VPD_a
: Vapor pressure deficit (at $z_a$) [Pa]A
: Total available energy [W/m²]A_c
: Available energy for canopy [W/m²]A_s
: Available energy for soil [W/m²]r_aa
: Aerodynamic resistance between canopy source height and observation height [s/m]r_ac
: Boundary layer resistance i.e. excess resistance to heat transfer [s/m]r_as
: Aerodynamic resistance between soil and canopy source height [s/m]r_sc
: Surface resistance for canopy [s/m]r_ss
: Surface resistance for soil [s/m]f_wet
: Fraction of canopy that is wet [-]
Returns
λE
: Total latent heat flux [W/m²]λE_p
: Potential latent heat flux [W/m²]
Details
For a model description, see TO DO ADD FULL MODEL DESCRIPTION IN DOCS.
The calculation is an extension of the model from Shuttleworth & Wallace (1985) to include interception. Equations are written in the notation of Lhomme et al. (2012), see e.g. equations (16) and (33)
Ground heat flux
EvaporationModel.compute_g_from_r_n
— Methodcompute_g_from_r_n(R_n, lai)
Compute the ground heat flux [W/m²] from net radiation (and LAI)
Arguments
R_n
: The net radiation [W/m²].lai
: The leaf area index [m²/m²].
Details
This implementation follows the approach of the METRIC model. See Equation 27 of Allen et al., 2007
EvaporationModel.compute_harmonic_sum
— Methodcompute_harmonic_sum(t::AbstractVector, a_bn::AbstractVector, ϕ::AbstractVector,
ω::Real, Δt::Real)
Applies broadcasting of the function for when t::AbstractVector
EvaporationModel.compute_harmonic_sum
— Methodcompute_harmonic_sum(t::Real, a_bn::AbstractVector, ϕ::AbstractVector,
ω::AbstractVector, Δt::Int)
Compute the sum of harmonic terms \Gamma_s
as defined in equation 1 of Murray and Verhoef (2007)
Soil
EvaporationModel.c_1
— Methodc_1(w_1, w_sat, b, c_1sat)
Compute force coefficient c_1
[-] of force restore framework for soil mositure. See equation 20 of Noilhan & Mahfouf, 1996.
Arguments
EvaporationModel.c_1sat
— Methodc_1sat(x_clay)
Compute the value for force-restore coefficient c_1
when w_g = w_sat
. See equation 32 of Noilhan & Mahfouf, 1996.
Arguments
x_clay
: The percentage of clay in the soil [%]
Returns
c_1sat
: The computed value ofc_1sat
[-]
EvaporationModel.c_2
— Methodc_2(w_2, w_sat, c2_ref)
Compute restore coefficient c_2
[-] of force restore framework for soil moisture See equation 21 of Noilhan & Mahfouf, 1996.
Arguments
w_2
: The second layer soil mositure [m³ m⁻³]w_sat
: The saturated soil moisture [m³ m⁻³]c_2ref
: Seec_2ref
EvaporationModel.c_2ref
— Methodc_2ref(x_clay)
Compute the value for force-restore coefficient c_2
[-] when w_2 = 0.5 w_sat
, c_2ref
based on the percentage of clay in the soil. See equation 33 of Noilhan & Mahfouf, 1996.
Arguments
x_clay
: The percentage of clay in the soil [%].
EvaporationModel.c_3
— Methodc_3(x_clay)
Compute the coefficient for graviational drainage c_3
[m] based on the percentage of clay in the soil. See equation 34 of Noilhan & Mahfouf, 1996.
Arguments
x_clay
: The percentage of clay in the soil [%].
EvaporationModel.compute_a
— Methodcompute_a(x_clay)
Compute a
[-], a parameter for for w_geq
calculation, based on percentage of clay in the soil. See equation 35 of Noilhan & Mahfouf, 1996.
Arguments
x_clay
: The percentage of clay in the soil [%].
EvaporationModel.compute_b
— Methodcompute_b(approach::Clay, x_clay)
compute_b(approach::VanGenuchten, n)
Compute b
[-], the Brooks-Corey/Clapp-Hornberger parameter (see equation 1 of Clapp & Hornberger for its definition), based on percentage clay or the van Genuchten paramter n
.
Arguments
approach
: calculation approach, subtype ofbMethod
.
With approach = Clay()
:
x_clay
: The percentage of clay in the soil [%]
With approach = VanGenuchten()
:
n
: The Van Genuchten parametern
[-]
Details
For approach = Clay()
, Equation (30) of Noilhan & Mahfouf, 1996 is used.
Forapproach = VanGenuchten()
, the parameter equivalence between the Brooks-Corey and van Genuchten, is based on Morel-Seytoux et al., 1996. Note that in this paper, M
is equivalent to b
.
EvaporationModel.compute_p
— Methodcompute_p(x_clay)
Compute p
[-], a parameter for for w_geq
calculation, based on the percentage of clay in the soil. See equation 36 of Noilhan & Mahfouf, 1996.
Arguments
x_clay
: The percentage of clay in the soil [%].
EvaporationModel.w_geq
— Methodw_geq(w_2, w_sat, a, p)
Compute w_geq
[m³ m⁻³], the equilibrium surface soil moisture (i.e. when capillary and gravitational forces are in equilibrium). See equation 19 of Noilhan & Mahfouf, 1996.
Arguments
Resistances
EvaporationModel.beta_to_r_ss
— Methodbeta_to_r_ss(beta, r_as)
Calculate soil surface resistance $r_{ss}$ from soil evaporation efficiency β
Arguments
beta
: Soil evaporation efficiency [-]r_as
: Aerodynamic resistance between soil and canopy source height [s/m]
Details
Equation used derived from equivalency between equation (6) and (7) of Merlin et al. (2016)
EvaporationModel.r_ss_to_beta
— Methodr_ss_to_beta(r_ss, r_as)
Calculate soil evaporation efficiency β soil surface resistance $r_{ss}$
Arguments
r_ss
: Soil surface resistance [s/m]r_as
: Aerodynamic resistance between soil and canopy source height [s/m]
Details
Equation used derived from equivalency between equation (6) and (7) of Merlin et al. (2016)
EvaporationModel.soil_aerodynamic_resistance
— Functionsoil_aerodynamic_resistance(::Choudhury1988soil, ustar, h, d_c, z_0mc, z_0ms, η=3)
Calculate soil aerodynamic resistance, which is between the soil surface and the canopy source height ($z_m = z_{0mc} + d_{c}$)
Arguments
- approach: The approach to use for calculating soil aerodynamic resistance
- ustar: Friction velocity [m/s]
- h: Canopy height [m]
- d_c: Canopy displacement height [m]
- z_0mc: Roughness length for momentum transfer for canopy [m]
- z_0ms: Roughness length for momentum transfer for soil [m]
- η: Extinction coefficient of $K$ in canopy, default = 3 [-]
Returns
- r_as: Aerodynamic resistance between soil and canopy source height [s/m]
Details
With approach = Choudhury1988soil()
, equation (25) of Choudhury and Monteith (1988) is used.
EvaporationModel.soil_evaporation_efficiency
— Methodsoil_evaporation_efficiency(approach::Pielke92, w_1, w_fc)
soil_evaporation_efficiency(approach::Martens17, w_1, w_res, w_c)
Calculate soil evaporation efficiency β, a factor between 0 and 1 which scales the potential soil evaporation
Arguments
approach
: The approach to use for calculating soil evaporation efficiency, a subtype ofSoilEvaporationEfficiencyMethod
w_1
: Soil moisture from first layer [m³/m³]
With approach = Pielke92()
, the following inputs are
w_fc
: Soil moisture at field capacity [m³/m³]
With approach = Martens17()
, the following inputs are
w_res
: Residual soil moisture [m³/m³]w_c
: Critical soil moisture [m³/m³]
Details
With approach = Pielke92()
, β is calculated using euqation (7) of Lee & Pielke (1992)
With approach = Martens17()
, β is calculated using equation (6) of Martens et al. (2017).
Examples
using EvaporationModel
w_fc = 0.35
w_1 = w_fc / 2
β_p = soil_evaporation_efficiency(Pielke92(), w_1, w_fc)
β_p ≈ 0.25
# output
true
EvaporationModel.surface_resistance
— Methodsurface_resistance(::JarvisStewart, SW_in, VPD, T_a, w_2, w_fc, w_wilt, LAI, g_d, r_smin)
Calculate surface resistance
Arguments
approach
: The approach to use for calculating surface resistanceSW_in
: Incoming solar radiation [W/m²]VPD
: Vapor pressure deficit [Pa]T_a
: Air temperature [K]w_2
: Root zone soil moisture [m³/m³]w_fc
: Soil moisture at field capacity [m³/m³]w_wilt
: Soil moisture at wilting point [m³/m³]LAI
: Leaf area index [m²/m²]g_d
: Parameter relatingVPD
to surface conductance [Pa⁻¹]r_smin
: Minimum surface resistance [s/m]T_opt
: Optimum temperature for stomatal conductance [K], default = 298.0 Kr_smax
: Maximum surface resistance [s/m], default = 500_000 s/m
Returns
r_s
: Surface resistance [s/m]
Details
With approach = JarvisStewart()
, the Jarvis-Stewart model is used as given in equations (8.9), (8.10), (8.11) and (8.14) of the IFS Cy47r1 documentation Part IV: Physical processes. The temperature constraint (f_4
) comes from Noilhan and Planton (1989) equation (37), with the default value of $T_{opt}$ set to 298.0 K
EvaporationModel.ustar_from_u
— Functionustar_from_u(u, z_obs, d, z_0m, ψ_m=0)
Calculate the friction velocity from the wind speed at a given height assuming a logarithmic wind profile.
Arguments
u
: Wind speed at the measurement height [m/s]z_obs
: Height of the wind speed measurement [m]d
: Displacement height [m]z_0m
: Roughness length for momentum transfer [m]ψ_m
: Stability correction for momentum transfer [m], default = 0
Returns
u_star
: Friction velocity [m/s]
Details
The friction velocity is calculated using the logarithmic wind profile equation, given by:
$u(z) = \frac{u^*}{k} \ln \left( \frac{z - d}{z_{0m}} \right) - \psi_m$
For more info on how to calculate $\psi_m$, see the Bigleaf package documentation
Utils
EvaporationModel.compute_amplitude_and_phase
— Methodcompute_amplitude_and_phase(an::AbstractVector, bn::AbstractVector)
Translate the fourier coefficients from sin-cos form (see fourier_series) to amplitude phase form:
$f(t) = a_0 + \sum_{n=1}^M \left( a_{bn} \sin(n \omega t + \phi) \right)$
Specifically, it returns both a_{bn}
and \phi
vectors. FYI: proof of conversion found here. Crucial to use atan2 function
EvaporationModel.fourier_series
— Methodfourier_series(t, coeffs, ω)
Function that gives the Fourier series, as presented in equation 5 of Murray & Verhoef
$f(t) = a_0 + \sum_{n=1}^M \left( a_n \cos(n \omega t) + b_n \sin(n \omega t) \right)$
Arguments
t
: Timestamp in secondscoeffs::ComponentArray
: array holding the coefficientsa0
: Mean componentan
: Cosine coefficientsbn
: Sine coefficients
ω
: Radial frequency [rad/s]
EvaporationModel.local_to_solar_time
— Methodlocal_to_solar_time(time, timezone, lon)
Arguments
time
: DateTime object given the current timetimezone
: Timezone in hours ahead of UTC (e.g. +1 for Brussels)lon
: Longitude in degrees
Returns
t_sol
: DateTime object of the local solar time
Details
Caluculations based on equations from NOAA and pveducation. For the equation ot time, the equation from NOAA is used.
Examples
Below an example to check if calculation matches calculator provided here
using Dates
using EvaporationModel
lon = 150 # °
timezone = 10 # UTC+10
time = DateTime(2003, 1, 5, 12, 30)
t_sol = local_to_solar_time(time, timezone, lon)
true_t_sol_hour = 12
true_t_sol_minute = 24
# Check if correct within minute
[hour(t_sol) - true_t_sol_hour, minute(t_sol) - true_t_sol_minute] ≤ [0, 1]
# output
true
EvaporationModel.seconds_since_solar_noon
— Methodseconds_since_solar_noon(t_sol)
Given the local solar time t_sol
, this function returns t_diff
, which is the number of seconds since solar noon, which takes places at 12:00:00 of that day (in local solar time).