This function calculates daily radiation, condensation, and evaporation fluxes.
calc_daily_evap(
lat,
n,
elv = 0,
y = 0,
sf = 1,
tc = 23,
sw = 1,
ke = 0.0167,
keps = 23.44,
komega = 283,
kw = 0.26
)
double, decimal degrees.
double, day of year.
double, elevation, m A.S.L. Default: \(0\).
double, year. Default: \(0\).
double, fraction of sunshine hours. Default: \(1\).
double, mean daily air temperature, degrees C. Default: \(23.0\).
double, evaporative supply rate, mm/hr. Default: \(1.0\).
double, eccentricity of earth's orbit. Default: \(0.01670\), 2000CE (Berger, 1978).
double, obliquity of earth's elliptic. Default: \(23.44\), 2000CE (Berger, 1978).
double, lon. of perihelion, degrees Default: \(283\), 2000CE (Berger, 1978).
double, PET entrainment, \((1 + kw) * EET\) Default: \(0.26\) (Priestley-Taylor, 1972)
Returns a list
object with the following variables:
nu_deg ............ true anomaly, degrees
lambda_deg ........ true longitude, degrees
dr ................ distance factor, unitless
delta_deg ......... declination angle, degrees
hs_deg ............ sunset angle, degrees
ra_j.m2 ........... daily extraterrestrial radiation, J/m^2
tau ............... atmospheric transmittivity, unitless
ppfd_mol.m2 ....... daily photosyn photon flux density, mol/m^2
hn_deg ............ net radiation hour angle, degrees
rn_j.m2 ........... daily net radiation, J/m^2
rnn_j.m2 .......... daily nighttime net radiation, J/m^2
econ_m3.j ......... water to energy conversion, m^3/J
cond_mm ........... daily condensation, mm
eet_mm ............ daily equilibrium evapotranspiration, mm
pet_mm ............ daily potential evapotranspiration, mm
hi_deg ............ intersection hour angle, degrees
aet_mm ............ daily actual evapotranspiration, mm
Berger, A.L., 1978. Long-term variations of daily insolation and Quaternary climatic changes. Journal of Atmospheric Sciences, 35(12), pp.2362-2367. doi:10.1175/1520-0469(1978)035<2362:LTVODI>2.0.CO;2
Priestley, C.H.B. and Taylor, R.J., 1972. On the assessment of surface heat flux and evaporation using large-scale parameters. Monthly weather review, 100(2), pp.81-92. doi:10.1175/1520-0493(1972)100<0081:OTAOSH>2.3.CO;2
evap <- splash::calc_daily_evap(lat = 37.7,
n = 172,
elv = 142,
y = 2000,
sf = 1,
tc = 23.0,
sw = 0.9)
cat(sprintf("Evaporation values:\n"))
#> Evaporation values:
cat(sprintf(" s: %0.6f Pa/K\n", evap$s_pa.k))
#> s: 169.896093 Pa/K
cat(sprintf(" Lv: %0.6f MJ/kg\n", (1e-6) * evap$lv_j.kg))
#> Lv: 2.446687 MJ/kg
cat(sprintf(" Patm: %0.6f bar\n", (1e-5) * evap$patm_pa))
#> Patm: 0.996308 bar
cat(sprintf(" pw: %0.6f kg/m^3\n", evap$pw_kg.m3))
#> pw: 997.583620 kg/m^3
cat(sprintf(" gamma: %0.6f Pa/K\n", evap$gam_pa.k))
#> gamma: 66.729719 Pa/K
cat(sprintf(" Econ: %0.6f mm^3/J\n", (1e9) * evap$econ_m3.j))
#> Econ: 0.294167 mm^3/J
cat(sprintf(" Cn: %0.6f mm\n", evap$cond_mm))
#> Cn: 0.885192 mm
cat(sprintf(" rx: %0.6f\n", evap$rx))
#> rx: 0.001334
cat(sprintf(" hi: %0.6f degrees\n", evap$hi_deg))
#> hi: 20.959320 degrees
cat(sprintf(" EET: %0.6f mm\n", evap$eet_mm))
#> EET: 6.405468 mm
cat(sprintf(" PET: %0.6f mm\n", evap$pet_mm))
#> PET: 8.070889 mm
cat(sprintf(" AET: %0.6f mm\n", evap$aet_mm))
#> AET: 7.972788 mm