Title: | Solve for Leaf Temperature Using Energy Balance |
---|---|
Description: | Implements models of leaf temperature using energy balance. It uses units to ensure that parameters are properly specified and transformed before calculations. It allows separate lower and upper surface conductances to heat and water vapour, so sensible and latent heat loss are calculated for each surface separately as in Foster and Smith (1986) <doi:10.1111/j.1365-3040.1986.tb02108.x>. It's straightforward to model leaf temperature over environmental gradients such as light, air temperature, humidity, and wind. It can also model leaf temperature over trait gradients such as leaf size or stomatal conductance. Other references are Monteith and Unsworth (2013, ISBN:9780123869104), Nobel (2009, ISBN:9780123741431), and Okajima et al. (2012) <doi:10.1007/s11284-011-0905-5>. |
Authors: | Chris Muir [aut, cre] |
Maintainer: | Chris Muir <[email protected]> |
License: | MIT + file LICENSE |
Version: | 1.0.6 |
Built: | 2025-02-19 03:44:27 UTC |
Source: | https://github.com/cdmuir/tealeaves |
d_wv: water vapour gradient (mol / m ^ 3)
.get_dwv(T_leaf, pars, unitless)
.get_dwv(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Water vapour gradient: The water vapour pressure differential from inside to outside of the leaf is the saturation water vapor pressure inside the leaf (p_leaf) minus the water vapor pressure of the air (p_air):
Note that water vapor pressure is converted from kPa to mol / m^3 using ideal gas law.
Symbol | R | Description | Units | Default |
|
p_air |
saturation water vapour pressure of air | kPa | calculated |
|
p_leaf |
saturation water vapour pressure inside the leaf | kPa | calculated |
|
R |
ideal gas constant | J / (mol K) | 8.3144598 |
|
RH |
relative humidity | % | 0.50 |
|
T_air |
air temperature | K | 298.15 |
|
T_leaf |
leaf temperature | K | input |
Value in mol / m^3 of class units
# Water vapour gradient: leaf_par <- make_leafpar() enviro_par <- make_enviropar() constants <- make_constants() pars <- c(leaf_par, enviro_par, constants) T_leaf <- set_units(300, K) T_air <- set_units(298.15, K) p_leaf <- set_units(35.31683, kPa) p_air <- set_units(31.65367, kPa) d_wv <- p_leaf / (pars$R * T_leaf) - pars$RH * p_air / (pars$R * T_air)
# Water vapour gradient: leaf_par <- make_leafpar() enviro_par <- make_enviropar() constants <- make_constants() pars <- c(leaf_par, enviro_par, constants) T_leaf <- set_units(300, K) T_air <- set_units(298.15, K) p_leaf <- set_units(35.31683, kPa) p_air <- set_units(31.65367, kPa) d_wv <- p_leaf / (pars$R * T_leaf) - pars$RH * p_air / (pars$R * T_air)
D_x: Calculate diffusion coefficient for a given temperature and pressure
.get_Dx(D_0, Temp, eT, P, unitless)
.get_Dx(D_0, Temp, eT, P, unitless)
D_0 |
Diffusion coefficient at 273.15 K (0 °C) and 101.3246 kPa |
Temp |
Temperature in Kelvin |
eT |
Exponent for temperature dependence of diffusion |
P |
Atmospheric pressure in kPa |
unitless |
Logical. Should function use parameters with |
According to Montieth & Unger (2013), eT is generally between 1.5 and 2. Their data in Appendix 3 indicate is reasonable for environmental physics.
Value in m/s of class
units
Monteith JL, Unsworth MH. 2013. Principles of Environmental Physics. 4th edition. Academic Press, London.
tealeaves:::.get_Dx( D_0 = set_units(2.12e-05, m^2/s), Temp = set_units(298.15, K), eT = set_units(1.75), P = set_units(101.3246, kPa), unitless = FALSE )
tealeaves:::.get_Dx( D_0 = set_units(2.12e-05, m^2/s), Temp = set_units(298.15, K), eT = set_units(1.75), P = set_units(101.3246, kPa), unitless = FALSE )
g_bw: Boundary layer conductance to water vapour (m / s)
.get_gbw(T_leaf, surface, pars, unitless)
.get_gbw(T_leaf, surface, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
surface |
Leaf surface (lower or upper) |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Symbol | R | Description | Units | Default |
|
leafsize |
Leaf characteristic dimension in meters | m | 0.1 |
|
D_w |
diffusion coefficient for water vapour | m / s |
calculated |
|
Sh |
Sherwood number | none | calculated |
Value in m / s of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_gbw(T_leaf, "lower", c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_gbw(T_leaf, "lower", c(cs, ep, lp), FALSE)
g_h: boundary layer conductance to heat (m / s)
.get_gh(T_leaf, surface, pars, unitless)
.get_gh(T_leaf, surface, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
surface |
Leaf surface (lower or upper) |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Symbol | R | Description | Units | Default |
|
leafsize |
Leaf characteristic dimension in meters | m | 0.1 |
|
D_h |
diffusion coefficient for heat in air | m / s |
calculated |
|
Nu |
Nusselt number | none | calculated |
Value in m/s of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_gh(T_leaf, "lower", c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_gh(T_leaf, "lower", c(cs, ep, lp), FALSE)
Gr: Grashof number
.get_gr(T_leaf, pars, unitless)
.get_gr(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Symbol | R | Description | Units | Default |
|
leafsize |
Leaf characteristic dimension in meters | m | 0.1 |
|
D_m |
diffusion coefficient of momentum in air | m / s |
calculated |
|
G |
gravitational acceleration | m / s |
9.8 |
|
t_air |
coefficient of thermal expansion of air | 1 / K | 1 / Temp |
|
Tv_air |
virtual air temperature | K | calculated |
|
Tv_leaf |
virtual leaf temperature | K | calculated |
A unitless number of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_gr(T_leaf, c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_gr(T_leaf, c(cs, ep, lp), FALSE)
g_tw: total conductance to water vapour (m/s)
.get_gtw(T_leaf, pars, unitless)
.get_gtw(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Total conductance to water vapor: The total conductance to water vapor () is the sum of the parallel lower (abaxial) and upper (adaxial) conductances:
The conductance to water vapor on each surface is a function of parallel stomatal () and cuticular (
) conductances in series with the boundary layer conductance (
). The stomatal, cuticular, and boundary layer conductance on the lower surface are:
See .get_gbw
for details on calculating boundary layer conductance. The equations for the upper surface are:
Note that the stomatal and cuticular conductances are given in units of (mol H2O) / (m
s Pa) (see
make_leafpar
) and converted to m/s using the ideal gas law. The total leaf stomatal () and cuticular (
) conductances are partitioned across lower and upper surfaces. The stomatal conductance on each surface depends on stomatal ratio (sr); the cuticular conductance is assumed identical on both surfaces.
Symbol | R | Description | Units | Default |
|
g_sw |
stomatal conductance to H2O | ( mol H2O) / (m s Pa) |
5 |
|
g_uw |
cuticular conductance to H2O | ( mol H2O) / (m s Pa) |
0.1 |
|
R |
ideal gas constant | J / (mol K) | 8.3144598 |
|
logit_sr |
stomatal ratio (logit transformed) | none | 0 = logit(0.5) |
|
T_air |
air temperature | K | 298.15 |
|
T_leaf |
leaf temperature | K | input |
Value in m/s of class units
# Total conductance to water vapor ## Hypostomatous leaf; default parameters leaf_par <- make_leafpar(replace = list(logit_sr = set_units(-Inf))) enviro_par <- make_enviropar() constants <- make_constants() pars <- c(leaf_par, enviro_par, constants) T_leaf <- set_units(300, K) ## Fixing boundary layer conductance rather than calculating gbw_lower <- set_units(0.1, m / s) gbw_upper <- set_units(0.1, m / s) # Lower surface ---- ## Note that pars$logit_sr is logit-transformed! Use stats::plogis() to convert to proportion. gsw_lower <- set_units(pars$g_sw * (set_units(1) - stats::plogis(pars$logit_sr)) * pars$R * ((T_leaf + pars$T_air) / 2), "m / s") guw_lower <- set_units(pars$g_uw * 0.5 * pars$R * ((T_leaf + pars$T_air) / 2), m / s) gtw_lower <- 1 / (1 / (gsw_lower + guw_lower) + 1 / gbw_lower) # Upper surface ---- gsw_upper <- set_units(pars$g_sw * stats::plogis(pars$logit_sr) * pars$R * ((T_leaf + pars$T_air) / 2), m / s) guw_upper <- set_units(pars$g_uw * 0.5 * pars$R * ((T_leaf + pars$T_air) / 2), m / s) gtw_upper <- 1 / (1 / (gsw_upper + guw_upper) + 1 / gbw_upper) ## Lower and upper surface are in parallel g_tw <- gtw_lower + gtw_upper
# Total conductance to water vapor ## Hypostomatous leaf; default parameters leaf_par <- make_leafpar(replace = list(logit_sr = set_units(-Inf))) enviro_par <- make_enviropar() constants <- make_constants() pars <- c(leaf_par, enviro_par, constants) T_leaf <- set_units(300, K) ## Fixing boundary layer conductance rather than calculating gbw_lower <- set_units(0.1, m / s) gbw_upper <- set_units(0.1, m / s) # Lower surface ---- ## Note that pars$logit_sr is logit-transformed! Use stats::plogis() to convert to proportion. gsw_lower <- set_units(pars$g_sw * (set_units(1) - stats::plogis(pars$logit_sr)) * pars$R * ((T_leaf + pars$T_air) / 2), "m / s") guw_lower <- set_units(pars$g_uw * 0.5 * pars$R * ((T_leaf + pars$T_air) / 2), m / s) gtw_lower <- 1 / (1 / (gsw_lower + guw_lower) + 1 / gbw_lower) # Upper surface ---- gsw_upper <- set_units(pars$g_sw * stats::plogis(pars$logit_sr) * pars$R * ((T_leaf + pars$T_air) / 2), m / s) guw_upper <- set_units(pars$g_uw * 0.5 * pars$R * ((T_leaf + pars$T_air) / 2), m / s) gtw_upper <- 1 / (1 / (gsw_upper + guw_upper) + 1 / gbw_upper) ## Lower and upper surface are in parallel g_tw <- gtw_lower + gtw_upper
H: sensible heat flux density (W / m^2)
.get_H(T_leaf, pars, unitless)
.get_H(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Symbol | R | Description | Units | Default |
|
c_p |
heat capacity of air | J / (g K) | 1.01 |
|
g_h |
boundary layer conductance to heat | m / s | calculated |
|
P_a |
density of dry air | g / m^3 | calculated |
|
T_air |
air temperature | K | 298.15 |
|
T_leaf |
leaf temperature | K | input |
Value in W / m of class
units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_H(T_leaf, c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_H(T_leaf, c(cs, ep, lp), FALSE)
h_vap: heat of vaporization (J / mol)
.get_hvap(T_leaf, unitless)
.get_hvap(T_leaf, unitless)
T_leaf |
Leaf temperature in Kelvin |
unitless |
Logical. Should function use parameters with |
Heat of vaporization: The heat of vaporization () is a function of temperature. I used data from on temperature and
from Nobel (2009, Appendix 1) to estimate a linear regression. See Examples.
Value in J/mol of class units
Nobel PS. 2009. Physicochemical and Environmental Plant Physiology. 4th Edition. Academic Press.
# Heat of vaporization and temperature ## data from Nobel (2009) T_K <- 273.15 + c(0, 10, 20, 25, 30, 40, 50, 60) h_vap <- 1e3 * c(45.06, 44.63, 44.21, 44.00, 43.78, 43.35, 42.91, 42.47) # (in J / mol) fit <- lm(h_vap ~ T_K) ## coefficients are 56847.68250 J / mol and 43.12514 J / (mol K) coef(fit) T_leaf <- 298.15 h_vap <- set_units(56847.68250, J / mol) - set_units(43.12514, J / mol / K) * set_units(T_leaf, K) ## h_vap at 298.15 K is 43989.92 [J/mol] set_units(h_vap, J / mol) tealeaves:::.get_hvap(set_units(298.15, K), FALSE)
# Heat of vaporization and temperature ## data from Nobel (2009) T_K <- 273.15 + c(0, 10, 20, 25, 30, 40, 50, 60) h_vap <- 1e3 * c(45.06, 44.63, 44.21, 44.00, 43.78, 43.35, 42.91, 42.47) # (in J / mol) fit <- lm(h_vap ~ T_K) ## coefficients are 56847.68250 J / mol and 43.12514 J / (mol K) coef(fit) T_leaf <- 298.15 h_vap <- set_units(56847.68250, J / mol) - set_units(43.12514, J / mol / K) * set_units(T_leaf, K) ## h_vap at 298.15 K is 43989.92 [J/mol] set_units(h_vap, J / mol) tealeaves:::.get_hvap(set_units(298.15, K), FALSE)
L: Latent heat flux density (W / m^2)
.get_L(T_leaf, pars, unitless)
.get_L(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Symbol | R | Description | Units | Default |
|
d_wv |
water vapour gradient | mol / m ^ 3 | calculated |
|
h_vap |
latent heat of vaporization | J / mol | calculated |
|
g_tw |
total conductance to H2O | ( mol H2O) / (m s Pa) |
calculated |
Value in W / m^2 of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_L(T_leaf, c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_L(T_leaf, c(cs, ep, lp), FALSE)
Nu: Nusselt number
.get_nu(T_leaf, surface, pars, unitless)
.get_nu(T_leaf, surface, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
surface |
Leaf surface (lower or upper) |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
The Nusselt number depends on a combination how much free or forced convection predominates. For mixed convection:
Symbol | R | Description | Units | Default |
|
a, b, c, d |
empirical coefficients | none | calculated |
|
Gr |
Grashof number | none | calculated |
|
Re |
Reynolds number | none | calculated |
A unitless number of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_nu(T_leaf, "lower", c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_nu(T_leaf, "lower", c(cs, ep, lp), FALSE)
P_a: density of dry air (g / m^3)
.get_Pa(T_leaf, pars, unitless)
.get_Pa(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Symbol | R | Description | Units | Default |
|
P |
atmospheric pressure | kPa | 101.3246 |
|
R_air |
specific gas constant for dry air | J / (kg K) | 287.058 |
|
T_air |
air temperature | K | 298.15 |
|
T_leaf |
leaf temperature | K | input |
Value in g / m of class
units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_Pa(T_leaf, c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_Pa(T_leaf, c(cs, ep, lp), FALSE)
Saturation water vapour pressure (kPa)
.get_ps(Temp, P, unitless)
.get_ps(Temp, P, unitless)
Temp |
Temperature in Kelvin |
P |
Atmospheric pressure in kPa |
unitless |
Logical. Should function use parameters with |
Goff-Gratch equation (see http://cires1.colorado.edu/~voemel/vp.html)
This equation assumes P = 1 atm = 101.3246 kPa, otherwise boiling temperature needs to change
Value in kPa of class units
http://cires1.colorado.edu/~voemel/vp.html
T_leaf <- set_units(298.15, K) P <- set_units(101.3246, kPa) tealeaves:::.get_ps(T_leaf, P, FALSE)
T_leaf <- set_units(298.15, K) P <- set_units(101.3246, kPa) tealeaves:::.get_ps(T_leaf, P, FALSE)
R_abs: total absorbed radiation (W / m^2)
.get_Rabs(pars, unitless)
.get_Rabs(pars, unitless)
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
The following treatment follows Okajima et al. (2012):
The incident longwave (aka thermal infrared) radiation is modeled from sky and air temperature where
is function of the air temperature and incoming solar shortwave radiation:
Symbol | R | Description | Units | Default |
|
abs_s |
absorbtivity of shortwave radiation (0.3 - 4 m) |
none | 0.80 |
|
abs_l |
absorbtivity of longwave radiation (4 - 80 m) |
none | 0.97 |
|
r |
reflectance for shortwave irradiance (albedo) | none | 0.2 |
|
s |
Stefan-Boltzmann constant | W / (m K ) |
5.67e-08 |
|
S_sw |
incident short-wave (solar) radiation flux density | W / m |
1000 |
|
S_lw |
incident long-wave radiation flux density | W / m |
calculated |
|
T_air |
air temperature | K | 298.15 |
|
T_sky |
sky temperature | K | calculated |
Value in W / m of class
units
Okajima Y, H Taneda, K Noguchi, I Terashima. 2012. Optimum leaf size predicted by a novel leaf energy balance model incorporating dependencies of photosynthesis on light and temperature. Ecological Research 27: 333-46.
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() ep$T_sky <- ep$T_sky(ep) tealeaves:::.get_Rabs(c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() ep$T_sky <- ep$T_sky(ep) tealeaves:::.get_Rabs(c(cs, ep, lp), FALSE)
Re: Reynolds number
.get_re(T_leaf, pars, unitless)
.get_re(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
Symbol | R | Description | Units | Default |
|
leafsize |
Leaf characteristic dimension in meters | m | 0.1 |
|
D_m |
diffusion coefficient of momentum in air | m / s |
calculated |
|
wind |
windspeed | m / s | 2 |
A unitless number of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_re(T_leaf, c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_re(T_leaf, c(cs, ep, lp), FALSE)
Sh: Sherwood number
.get_sh(T_leaf, surface, pars, unitless)
.get_sh(T_leaf, surface, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
surface |
Leaf surface (lower or upper) |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
The Sherwood number depends on a combination how much free or forced convection predominates. For mixed convection:
Symbol | R | Description | Units | Default |
|
a, b, c, d |
empirical coefficients | none | calculated |
|
Gr |
Grashof number | none | calculated |
|
Re |
Reynolds number | none | calculated |
A unitless number of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_sh(T_leaf, "lower", c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_sh(T_leaf, "lower", c(cs, ep, lp), FALSE)
S_r: longwave re-radiation (W / m^2)
.get_Sr(T_leaf, pars)
.get_Sr(T_leaf, pars)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
The factor of 2 accounts for re-radiation from both leaf surfaces (Foster and Smith 1986).
Symbol | R | Description | Units | Default |
|
abs_l |
absorbtivity of longwave radiation (4 - 80 m) |
none | 0.97 |
|
T_air |
air temperature | K | 298.15 |
|
s |
Stefan-Boltzmann constant | W / (m K ) |
5.67e-08 |
Note that leaf absorbtivity is the same value as leaf emissivity
Value in W / m of class
units
Foster JR, Smith WK. 1986. Influence of stomatal distribution on transpiration in low-wind environments. Plant, Cell & Environment 9: 751-9.
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_Sr(T_leaf, c(cs, ep, lp))
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) tealeaves:::.get_Sr(T_leaf, c(cs, ep, lp))
Calculate virtual temperature
.get_Tv(Temp, p, P, epsilon, unitless)
.get_Tv(Temp, p, P, epsilon, unitless)
Temp |
Temperature in Kelvin |
p |
water vapour pressure in kPa |
P |
Atmospheric pressure in kPa |
epsilon |
ratio of water to air molar masses (unitless) |
unitless |
Logical. Should function use parameters with |
Eq. 2.35 in Monteith & Unsworth (2013)
Symbol | R | Description | Units | Default |
|
epsilon |
ratio of water to air molar masses | unitless | 0.622 |
|
p |
water vapour pressure | kPa | calculated |
|
P |
atmospheric pressure | kPa | 101.3246 |
Value in K of class units
Monteith JL, Unsworth MH. 2013. Principles of Environmental Physics. 4th edition. Academic Press, London.
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) p <- ep$RH * tealeaves:::.get_ps(T_leaf, ep$P, FALSE) tealeaves:::.get_Tv(T_leaf, p, ep$P, cs$epsilon, FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) p <- ep$RH * tealeaves:::.get_ps(T_leaf, ep$P, FALSE) tealeaves:::.get_Tv(T_leaf, p, ep$P, cs$epsilon, FALSE)
Ar: Archimedes number
Ar(T_leaf, pars, unitless = FALSE)
Ar(T_leaf, pars, unitless = FALSE)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
The Archimedes number is a dimensionless number that describes when free or forced convection dominates.
Symbol | R | Description | Units | Default |
|
Gr |
Grashof number | none | calculated |
|
Re |
Reynolds number | none | calculated |
unitless = TRUE
: A unitless number of class numeric
unitless = FALSE
: A unitless number of class units
Also returns Reynolds and Grashof numbers
cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() pars <- c(cs, lp, ep) T_leaf <- set_units(298.15, "K") Ar(T_leaf, pars)
cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() pars <- c(cs, lp, ep) T_leaf <- set_units(298.15, "K") Ar(T_leaf, pars)
Constructor function for constants class. This function ensures that physical constant inputs are properly formatted.
constants(.x)
constants(.x)
.x |
A list to be constructed into constants. If units are not provided, they will be set without conversion. If units are provided, they will be checked and converted to units that tealeaves uses. |
Convert conductance units
convert_conductance(.g, Temp = NULL, P = NULL)
convert_conductance(.g, Temp = NULL, P = NULL)
.g |
Conductance in class units. Units must convertible to one of "m/s", "umol/m^2/s/Pa", or "mol/m^2/s" |
Temp |
A temperature value of class |
P |
A pressure value of class |
A list of three values of clas units
with units "m/s", "umol/m^2/s/Pa", and "mol/m^2/s".
g_sw <- set_units(10, "m/s") convert_conductance(g_sw, Temp = set_units(298.15, "K"), P = set_units(101.3246, "kPa")) g_sw <- set_units(4, "umol/m^2/s/Pa") convert_conductance(g_sw, Temp = set_units(298.15, "K"), P = set_units(101.3246, "kPa")) g_sw <- set_units(0.4, "mol/m^2/s") convert_conductance(g_sw, Temp = set_units(298.15, "K"), P = set_units(101.3246, "kPa"))
g_sw <- set_units(10, "m/s") convert_conductance(g_sw, Temp = set_units(298.15, "K"), P = set_units(101.3246, "kPa")) g_sw <- set_units(4, "umol/m^2/s/Pa") convert_conductance(g_sw, Temp = set_units(298.15, "K"), P = set_units(101.3246, "kPa")) g_sw <- set_units(0.4, "mol/m^2/s") convert_conductance(g_sw, Temp = set_units(298.15, "K"), P = set_units(101.3246, "kPa"))
Evaporation (mol / (m^2 s))
E(T_leaf, pars, unitless)
E(T_leaf, pars, unitless)
T_leaf |
Leaf temperature in Kelvin |
pars |
Concatenated parameters ( |
unitless |
Logical. Should function use parameters with |
The leaf evaporation rate is the product of the total conductance to water vapour (m / s) and the water vapour gradient (mol / m^3):
If unitless = TRUE
, T_leaf
is assumed in degrees K without checking.
unitless = TRUE
: A value in units of mol / (m ^ 2 / s) number of class numeric
unitless = FALSE
: A value in units of mol / (m ^ 2 / s) of class units
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) E(T_leaf, c(cs, ep, lp), FALSE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() T_leaf <- set_units(298.15, K) E(T_leaf, c(cs, ep, lp), FALSE)
Calculate leaf energy balance
energy_balance( tleaf, leaf_par, enviro_par, constants, quiet = FALSE, components = FALSE, set_units = FALSE )
energy_balance( tleaf, leaf_par, enviro_par, constants, quiet = FALSE, components = FALSE, set_units = FALSE )
tleaf |
Leaf temperature in Kelvin. If input is numeric, it will be automatically converted to |
leaf_par |
A list of leaf parameters. This can be generated using the |
enviro_par |
A list of environmental parameters. This can be generated using the |
constants |
A list of physical constants. This can be generated using the |
quiet |
Logical. Should a message appear about conversion from |
components |
Logical. Should leaf energy components be returned? Transpiration (in mol / (m^2 s)) also returned. |
set_units |
Logical. Should |
A numeric value in W / m^2. Optionally, a named list of energy balance components in W / m^2 and transpiration in mol / (m^2 s).
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() ep$T_sky <- ep$T_sky(ep) T_leaf <- set_units(298.15, K) energy_balance(T_leaf, lp, ep, cs, FALSE, TRUE, TRUE)
library(tealeaves) cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() ep$T_sky <- ep$T_sky(ep) T_leaf <- set_units(298.15, K) energy_balance(T_leaf, lp, ep, cs, FALSE, TRUE, TRUE)
Constructor function for enviro_par class. This function ensures that environmental parameter inputs are properly formatted.
enviro_par(.x)
enviro_par(.x)
.x |
A list to be constructed into enviro_par. If units are not provided, they will be set without conversion. If units are provided, they will be checked and converted to units that tealeaves uses. |
Constructor function for leaf_par class. This function ensures that leaf parameter inputs are properly formatted.
leaf_par(.x)
leaf_par(.x)
.x |
A list to be constructed into leaf_par. If units are not provided, they will be set without conversion. If units are provided, they will be checked and converted to units that tealeaves uses. |
Make lists of parameters of leaf, environmental, or constant parameters
make_leafpar
make_enviropar
make_constants
make_leafpar(replace = NULL) make_enviropar(replace = NULL) make_constants(replace = NULL)
make_leafpar(replace = NULL) make_enviropar(replace = NULL) make_constants(replace = NULL)
replace |
A named list of parameters to replace defaults. If |
Leaf parameters:
Symbol | R | Description | Units | Default |
|
leafsize |
Leaf characteristic dimension | m | 0.1 |
|
abs_l |
absorbtivity of longwave radiation (4 - 80 m) |
none | 0.97 |
|
abs_s |
absorbtivity of shortwave radiation (0.3 - 4 m) |
none | 0.50 |
|
g_sw |
stomatal conductance to H2O | ( mol H2O) / (m s Pa) |
5 |
|
g_uw |
cuticular conductance to H2O | ( mol H2O) / (m s Pa) |
0.1 |
|
logit_sr |
stomatal ratio (logit transformed) | none | 0 = logit(0.5) |
Environment parameters:
Symbol | R | Description | Units | Default |
|
P |
atmospheric pressure | kPa | 101.3246 |
|
r |
reflectance for shortwave irradiance (albedo) | none | 0.2 |
|
RH |
relative humidity | none | 0.50 |
|
S_sw |
incident short-wave (solar) radiation flux density | W / m |
1000 |
|
S_lw |
incident long-wave radiation flux density | W / m |
calculated |
|
T_air |
air temperature | K | 298.15 |
|
wind |
windspeed | m / s | 2 |
Constants:
Symbol | R | Description | Units | Default |
|
c_p |
heat capacity of air | J / (g K) | 1.01 |
|
D_h0 |
diffusion coefficient for heat in air at 0 °C | m / s |
19.0e-06 |
|
D_m0 |
diffusion coefficient for momentum in air at 0 °C | m / s |
13.3e-06 |
|
D_w0 |
diffusion coefficient for water vapour in air at 0 C | m / s |
21.2e-06 |
|
epsilon |
ratio of water to air molar masses | none | 0.622 |
|
eT |
exponent for temperature dependence of diffusion | none | 1.75 |
|
G |
gravitational acceleration | m / s |
9.8 |
|
Nu |
Nusselt number | none | calculated |
|
R |
ideal gas constant | J / (mol K) | 8.3144598 |
|
R_air |
specific gas constant for dry air | J / (kg K) | 287.058 |
|
s |
Stefan-Boltzmann constant | W / (m K ) |
5.67e-08 |
|
Sh |
Sherwood number | none | calculated |
make_leafpar
: An object inheriting from class leaf_par
make_enviropar
: An object inheriting from class enviro_par
make_constants
: An object inheriting from class constants
library(tealeaves) # Use defaults cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() # Replace defaults ep <- make_enviropar( replace = list( T_air = set_units(300, K) ) ) lp <- make_leafpar( replace = list( leafsize = set_units(c(0.1, 0.2), m) ) )
library(tealeaves) # Use defaults cs <- make_constants() ep <- make_enviropar() lp <- make_leafpar() # Replace defaults ep <- make_enviropar( replace = list( T_air = set_units(300, K) ) ) lp <- make_leafpar( replace = list( leafsize = set_units(c(0.1, 0.2), m) ) )
Get vector of parameter names
parameter_names(which)
parameter_names(which)
which |
A character string indicating which parameter names to retrieve, "constants", "enviro", or "leaf". Partial matching allowed. |
parameter_names("leaf")
parameter_names("leaf")
tealeaves
packageSolve for Leaf Temperature Using Energy Balance
See the README on GitHub
An example output from the tealeaves
function.
tl_example1
tl_example1
A data frame with 150 rows and 20 variables:
tleaves
: find leaf temperatures for multiple parameter setstleaves
: find leaf temperatures for multiple parameter sets
tleaf
: find leaf temperatures for a single parameter set
tleaves( leaf_par, enviro_par, constants, progress = TRUE, quiet = FALSE, set_units = TRUE, parallel = FALSE ) tleaf(leaf_par, enviro_par, constants, quiet = FALSE, set_units = TRUE)
tleaves( leaf_par, enviro_par, constants, progress = TRUE, quiet = FALSE, set_units = TRUE, parallel = FALSE ) tleaf(leaf_par, enviro_par, constants, quiet = FALSE, set_units = TRUE)
leaf_par |
A list of leaf parameters. This can be generated using the |
enviro_par |
A list of environmental parameters. This can be generated using the |
constants |
A list of physical constants. This can be generated using the |
progress |
Logical. Should a progress bar be displayed? |
quiet |
Logical. Should messages be displayed? |
set_units |
Logical. Should |
parallel |
Logical. Should parallel processing be used via |
tleaves
:
A tibble with the following units
columns
Input: | |
abs_l |
Absorbtivity of longwave radiation (unitless) |
abs_s |
Absorbtivity of shortwave radiation (unitless) |
g_sw |
Stomatal conductance to H2O ( mol H2O / (m^2 s Pa)) |
g_uw |
Cuticular conductance to H2O ( mol H2O / (m^2 s Pa)) |
leafsize Leaf characteristic dimension |
(m) |
logit_sr |
Stomatal ratio (logit transformed; unitless) |
P |
Atmospheric pressure (kPa) |
RH |
Relative humidity (unitless) |
S_lw |
incident long-wave radiation flux density (W / m^2) |
S_sw |
incident short-wave (solar) radiation flux density (W / m^2) |
T_air |
Air temperature (K) |
wind |
Wind speed (m / s) |
Output: | |
T_leaf |
Equilibrium leaf tempearture (K) |
value |
Leaf energy balance (W / m^2) at tleaf |
convergence |
Convergence code (0 = converged) |
R_abs |
Total absorbed radiation (W / m^2; see .get_Rabs ) |
S_r |
Thermal infrared radiation loss (W / m^2; see .get_Sr ) |
H |
Sensible heat flux density (W / m^2; see .get_H ) |
L |
Latent heat flux density (W / m^2; see .get_L ) |
E |
Evapotranspiration (mol H2O/ (m^2 s)) |
tleaf
:
A data.frame with the following numeric columns:
T_leaf |
Equilibrium leaf temperature (K) |
value |
Leaf energy balance (W / m^2) at tleaf |
convergence |
Convergence code (0 = converged) |
R_abs |
Total absorbed radiation (W / m^2; see .get_Rabs ) |
S_r |
Longwave re-radiation (W / m^2; see .get_Sr ) |
H |
Sensible heat flux density (W / m^2; see .get_H ) |
L
|
Latent heat flux density (W / m^2; see .get_L ) |
E |
Evapotranspiration (mol H2O/ (m^2 s)) |
# tleaf for single parameter set: leaf_par <- make_leafpar() enviro_par <- make_enviropar() constants <- make_constants() tleaf(leaf_par, enviro_par, constants) # tleaves for multiple parameter set: enviro_par <- make_enviropar( replace = list( T_air = set_units(c(293.15, 298.15), K) ) ) tleaves(leaf_par, enviro_par, constants)
# tleaf for single parameter set: leaf_par <- make_leafpar() enviro_par <- make_enviropar() constants <- make_constants() tleaf(leaf_par, enviro_par, constants) # tleaves for multiple parameter set: enviro_par <- make_enviropar( replace = list( T_air = set_units(c(293.15, 298.15), K) ) ) tleaves(leaf_par, enviro_par, constants)