evapotranspiration package

Submodules

evapotranspiration.penman_monteith_daily module

class evapotranspiration.penman_monteith_daily.PenmanMonteithDaily(elevation, latitude, **kwargs)

Bases: object

The class PenmanMonteithDaily calculates daily potential evapotranspiration according to the Penman-Monteith method as described in FAO 56 (Allen et al., 1998). Reference evapotranspiration for a hypothetical grass reference crop (\(h=12\) cm; \(albedo=0.23\), and \(LAI=2.88\)) is calculated by default. Wind and humidity observations at 2 meters height as well as soil heat flux density \(G=0.0\) MJ/m²day are also assumed by default. Default values can be changed in the keyword arguments (**kwargs) described below.

The class PenmanMonteithDaily solves equation 3 in FAO 56:

\[ET = \frac{\Delta (R_n - G) + \rho_a c_p \frac{e_s - e_a}{r_a}} {\lambda \left[ \Delta + \gamma \left( 1 + \frac{r_s}{r_a} \right) \right]} \tag{eq. 3, p. 19}\]
Parameters
Keyword Arguments
  • albedo (float) - albedo or canopy reflection coefficient (\(\alpha\)) [-]. Range: \(0.0 \leq \alpha \leq 1.0\). Default \(albedo=0.23\) for the hypothetical grass reference crop. Used in net_shortwave_radiation()

  • h (float) - crop height (h) [m]. Default \(h=0.12\) for the hypothetical grass reference crop. Required to calculate the zero plane displacement height (\(d\)) [m] and the roughness length governing momentum (\(z_{om}\)) [m], both necessary for the aerodynamic resistance (\(r_a\)) [s/m]. See aerodynamic_resistance_factor()

  • lai (float) - leaf area index (\(LAI\)) [-]. Default \(lai=2.88\) for the hypothetical grass reference crop. See BOX 5 in FAO 56 and bulk_surface_resistance()

  • rl (float) - bulk stomatal resistance of well-illuminated leaf (\(r_l\)) [s/m]. Default \(rl=100.0\) for any crop. See bulk_surface_resistance()

  • zm (float) - height of wind measurements (\(z_m\)) [m]. Default \(zm=2.0\). Required to calculate aerodynamic resistance (\(r_a\)) [s/m]. See aerodynamic_resistance_factor()

  • zh (float) - height of humidity measurements (\(z_h\)) [m]. Default \(zh=2.0\). Required to calculate aerodynamic resistance (\(r_a\)) [s/m]. See aerodynamic_resistance_factor()

  • g (float) - soil heat flux density (\(G\)) [MJ/m²day]. Default \(g=0.0\). This corresponds to \(G\) in eq. 3, p. 19 above. It can be also given with daily parameters in et0()

Note

Only elevation and latitude are mandatory parameters of PenmanMonteithDaily(). albedo, h, and lai are only necessary when calculating evapotranspiration for crops other than reference grass.

Variables
Object Constants:
  • e - ratio molecular weight of water vapour/dry air (\(\varepsilon\)) [-]. \(e = 0.622\)

  • r - specific gas constant [kJ/kg.K]. \(r = 0.287\)

  • k - von Karman constant (\(k\)) [-], see FAO 56 eq. 4. \(k=0.41\)

Object crop specific factors:
  • d_factor - factor of the zero plane displacement height (\(d\)) [-]. \(d\_factor = 2.0 / 3.0\)

  • zom_factor - factor of the roughness length governing momentum transfer (\(z_{om}\)) [-]. \(zom\_factor = 0.123\)

  • zoh_factor - factor of the roughness length governing transfer of heat and vapour (\(z_{oh}\)) [-]. \(zoh\_factor = 0.1\)

  • lai_active_factor - factor of the active (sunlit) leaf area index (\(LAI_{active}\)) [-] (it considers that generally only the upper half of dense clipped grass is actively contributing to the surface heat and vapour transfer). \(lai\_active\_factor = 0.5\)

Calculation with et0():

- pm = PenmanMonteithDaily(elevation, latitude, ...)
- et0 = pm.et0(...)

Calculation with et0_frame() given a pandas.DataFrame() as input parameter:

- pm = PenmanMonteithDaily(elevation, latitude, ...)
- df = pm.et0_frame(df, ...)
f1 = None

f1 = (specific heat at constant pressure) * (mean air density at constant pressure) / (1.01 * r * aerodynamic_resistance_factor()). FAO 56 Box 6

f2 = None

\(f_1 = \frac{rs}{f_{ra}}\) with \(f_{ra}\) = aerodynamic_resistance_factor()

reset()

Reset the following output attributes before calculating \(ETo\): \(doy\), \(u2\), \(ld\), \(s\), \(pc\), \(mn\), \(es\), \(ea\), \(ra\), \(rs\), \(rs0\), \(rns\), \(rnl\), \(rn\), \(etr\), \(etw\), and \(et\)

static atmospheric_pressure(z)

Return the atmospheric pressure (\(P\)) [kPa] as a function of the elevation above sea level as defined in FAO 56 (eq. 7, p. 31):

\[P = 101.3\left(\frac{293-0.0065z}{293}\right)^{5.26}\]

The atmospheric pressure (\(P\)) is the pressure exerted by the weight of the earth’s atmosphere. Evaporation at high altitudes is promoted due to low atmospheric pressure as expressed in the psychrometric constant. The effect is, however, small and in the calculation procedures, the average value for a location is sufficient. A simplification of the ideal gas law, assuming \(20\) °C for a standard atmosphere, can be employed to calculate \(P\) (FAO 56).

Parameters

z (float or np.array) – elevation above sea level [m]

Returns

(float or np.array) atmospheric pressure (\(P\)) [kPa]

static latent_heat_of_vaporization(temperature=20)

Return the latent heat of vaporization (\(\lambda\)) [MJ/kg] as described in FAO 56 (Annex 3, eq. 3-1, p. 223):

\[\lambda = 2.501-(2.361 * 10^{-3})T\]
Parameters

temperature (float or np.array) – air temperature (\(T\)) [°C]. Default \(temperature=20\)

Returns

(float or np.array) latent heat of vaporization (\(\lambda\)) [MJ/kg]. Default \(\lambda=2.45378\)

static psychrometric_constant(p, **kwargs)

Return the psychrometric constant (\(\gamma\)) [kPa/°C] according to FAO 56 eq. 8, p. 32:

\[\gamma = \frac{c_p P}{\varepsilon \lambda}\]

or, using default values:

\[\gamma = a_{psy} \cdot P\]
Parameters

p (float or np.array) – atmospheric pressure (\(P\)) [kPa]

Keyword Arguments
  • lamda (float) - latent heat of vaporization (\(\lambda\)) [MJ/kg]. Default \(lamda=2.45\). See Used in latent_heat_of_vaporization()

  • cp (float) - specific heat at constant pressure (\(c_p\)) [MJ/kg]. Default \(cp=1.013e^{-3}\)

  • epsilon (float) - ratio molecular weight of water vapour/dry air (\(\epsilon\)) [-]. Default \(epsilon=0.622\)

  • a_psy (float) - coefficient depending on the type of the ventilation of the bulb [1/°C]. Examples:

    • \(a_{psy} = 0.000665\) (default)

    • \(a_{psy} = 0.000662\) for ventilated (Asmann type) psychrometers, with an air movement of some 5 m/s

    • \(a_{psy} = 0.000800\) for natural ventilated psychrometers (about 1 m/s)

    • \(a_{psy} = 0.001200\) for non-ventilated psychrometers installed indoors

The method uses \(a_{psy}\) if given, otherwise eq. 8 (see above) with given or default values. Default values correspond to \(a_{psy} = 0.000665\) as argument.

Returns

(float or np.array) psychrometric constant (\(\gamma\)) [kPa/°C]

static saturation_vapour_pressure(*temperature)

Return the saturation vapour pressure (\(e_s\)) [kPa] according to FAO 56 (eq. 11, p. 36):

\[e^{°}(T) = 0.6108 exp \left[\frac{17.27 T}{T + 237.3}\right]\]
Parameters

temperature (float or np.array) – air temperature (\(T\)) [°C]

Returns

(float or np.array) saturation vapour pressure (\(e_s\)) [kPa]

static slope_of_saturation_vapour_pressure_curve(*temperature)

Return the slope of saturation vapour pressure curve (\(\Delta\)) [kPa/°C] according to FAO 56 (eq. 13, p. 37):

\[\Delta = 4098\left[\frac{0.6108exp\left(\frac{17.27 T}{T + 237.3}\right)}{(T + 237.3)^{2}}\right]\]
Parameters

temperature (float or np.array) – air temperature (\(T\)) [°C]

Returns

(float or np.array) slope of saturation vapour pressure curve (\(\Delta\)) [kPa/°C]

static actual_vapour_pressure(**kwargs)

Return the actual vapour pressure (\(e_a\)) [kPa] as defined in FAO 56 (p. 37 , 38 , and 39):

Keyword Arguments
  • rh_min (float) - 0.0 to 100.0 [%]

  • rh_max (float) - 0.0 to 100.0 [%]

  • es_min (float) - saturation vapour pressure for \(t\_min\) [kPa]

  • es_max (float) - saturation vapour pressure for \(t\_max\) [kPa]

  • t_min (float) - minimum air temperature [°C]

  • t_max (float) - maximum air temperature [°C]

  • t_dew (float) - dew point temperature [°C]

  • t_wet (float) - wet bulb temperature [°C]

  • t_dry (float) - dry bulb temperature [°C]

  • apsy (float) - coefficient depending on the type of ventilation of the wet bulb [-]

Returns

(float or np.array) actual vapour pressure (\(e_a\)) [kPa]

aerodynamic_resistance_factor()

Return the aerodynamic resistance (\(r_a\)) [s/m] as defined in FAO 56 (eq. 4, p. 20):

\[r_a = \frac{ \ln \left( \frac{z_m - d}{z_{om}} \right) \ln \left( \frac{z_h - d}{z_{oh}} \right) } { k^2 u_z }\]

where (see PenmanMonteithDaily()):

\(u_z\) — the wind speed [m/s] at height \(z\) (see et0())

\(k\) — von Karman’s constant [-]

\(zm\) — height of wind measurements [m]

\(zh\) — height of air humidity measurements [m]

The aerodynamic resistance factor \(f_{r_a}\) is constant for a given crop:

\[f_{r_a} = \frac{ \ln \left( \frac{z_m - d}{z_{om}} \right) \ln \left( \frac{z_h - d}{z_{oh}} \right) } { k^2}\]

with the zero plane displacement height (\(d\)):

\[d = f_d \cdot h\]

and roughness length governing momentum transfer (\(z_{om}\)):

\[z_{om} = f_{zom} \cdot h\]

where:

\(f_d\) — defined in d_factor

\(f_{zom}\) — defined in in zom_factor

Returns

(float) aerodynamic resistance factor \(f_{r_a}\)

bulk_surface_resistance()

Return (bulk) surface resistance (\(r_s\)) [s/m] as defined in FAO 56 (eq. 5, p. 21):

\[r_s = \frac{ r_l } { LAI_{active} }\]

where:

\(r_l\) — the bulk stomatal resistance of the well-illuminated leaf [s/m]

\(LAI_{active}\) — the active (sunlit) leaf area index [m² (leaf area) / m² (soil surface)]

A general equation for \(LAI_{active}\) is:

\[LAI_{active} = 0.5 LAI\]

with:

\[LAI = 24 h\]

where \(h\) is an optional input parameter in PenmanMonteithDaily.

Returns

(float) (bulk) surface resistance \(r_s\) [s/m]

static to_u2(uz, z)

Return the calculated wind speed at 2 meters above ground surface (\(u_2\)) [m/s]:

\[u_2 = \frac{ 4.87 u_z}{ \ln{(67.8 z - 5.42)}}\]
Parameters
  • uz (float or np.array) – measured wind speed at \(z\) meters above ground surface [m/s]

  • z (float) – height of measurement above ground surface [m]

Returns

(float or np.array) wind speed at 2 meters above ground surface [m/s]

static extraterrestrial_radiation(dr, ws, lat, sd)

Return the extraterrestrial radiation (\(R_a\)) [MJ/m²day] as defined in FAO 56 (eq. 21, p. 46):

\[R_a = \frac{24(60)}{\pi} G_{sc} d_r [ \omega_s \sin(\varphi) \sin(\delta) + \cos(\varphi) \cos(\delta) \sin(\omega_s)]\]
Parameters
Returns

(float or np.array) daily extraterrestrial radiation (\(R_a\)) [MJ/m²day]

static inverse_relative_distance_earth_sun(day)

Return the inverse relative distance Earth-Sun (\(d_r\)) [-] as defined in FAO 56 (eq. 23, p. 46):

\[d_r = 1 + 0.033 \cos{ \left( \frac{2 \pi}{365} J \right)}\]
Parameters

day (int or np.array) – day of the year (\(J\)) [-]. Range: \(1 \leq J \leq 366\)

Returns

(float or np.array) inverse relative distance Earth-Sun (\(d_r\)) [-]

static solar_declination(day)

Return the solar declination (\(\delta\)) [rad] as defined in FAO 56 (eq. 24, p. 46):

\[\delta = 0.409 \sin{ \left( \frac{2 \pi}{365} J - 1.39\right)}\]
Parameters

day (int) – day of the year (\(J\)) [-]. Range: \(1 \leq J \leq 366\)

Returns

(float or np.array) solar declination (\(\delta\)) [rad]

static sunset_hour_angle(lat, sd)

Return the sunset hour angle (\(\omega_s\)) [rad] as defined in FAO 56 (eq. 25, p. 46):

\[\omega_s = \arccos{ \left[-tan(\varphi)tan(\delta)\right]}\]
Parameters
  • lat (float or np.array) – latitude (\(\varphi\)) [rad]

  • sd (float or np.array) – solar declination (\(\delta\)) [rad]. See solar_declination()

Returns

(float or np.array) sunset hour angle (\(\omega_s\)) [rad]

static daylight_hours(ws)

Return the daylight hours (\(N\)) [hour] as defined in FAO 56 (eq. 34, p. 49):

\[N = \frac{24}{\pi} \omega_s\]
Parameters

ws (float or np.numpy) – sunset hour angle (\(\omega_s\)) [rad]. See sunset_hour_angle()

Returns

(float or np.numpy) daylight hours (\(N\)) [hour]

static clear_sky_shortwave_radiation(ra, elevation=0.0, a_s=0.25, b_s=0.5)

Return the clear-sky shortwave radiation (\(R_{so}\)) [MJ/m²day]. It is required for computing net_longwave_radiation().

For near sea level or when calibrated values for \(a_s\) and \(b_s\) are available (FAO 56, eq. 36, p. 51):

\[R_{so} = (a_s + b_s ) R_a\]

When calibrated values for \(a_s\) and \(b_s\) are not available (FAO 56, eq. 37, p. 51):

\[R_{so} = (0.75 + 2 * 10^{−5} z) R_a\]

where \(z\) is the station elevation above sea level [m].

Parameters
  • ra (float or np.numpy) – extraterrestrial radiation (\(R_a\)) [MJ/m²day]. See extraterrestrial_radiation()

  • elevation (float or np.numpy) – meters above sea level see (\(z\)) [m]. See elevation

  • a_s (float or np.numpy) – regression constant (\(a_s\)) [-]. Default \(a_s=0.25\). It expresses the fraction of extraterrestrial radiation reaching the earth on overcast days (\(n = 0\))

  • b_s (float or np.numpy) – regression constant (\(b_s\)) [-]. Default \(b_s=0.50\). The expression \(a_s+b_s\) indicates the fraction of extraterrestrial radiation reaching the earth on clear days (\(n = N\))

Returns

(float or np.numpy) daily clear-sky shortwave radiation (\(R_{so}\)) [MJ/m²day]

static shortwave_radiation(ra, n, mn, a_s=0.25, b_s=0.5)

Return the daily shortwave radiation (\(R_s\)) [MJ/m²day] according to the Angstrom formula as described in FAO 56 (eq. 35, p. 50):

\[R_s = \left( a_s + b_s \frac{n}{N} \right) R_a\]

Depending on atmospheric conditions (humidity, dust) and solar declination (latitude and month), the Angstrom values \(a_s\) and \(b_s\) will vary. Where no actual solar radiation data are available and no calibration has been carried out for improved \(a_s\) and \(b_s\) parameters, the values \(a_s = 0.25\) and \(b_s = 0.50\) are recommended.

Parameters
  • ra (float or np.array) – extraterrestrial radiation (\(R_a\)) [MJ/m²day]. See extraterrestrial_radiation()

  • n (float or np.array) – actual duration of sunshine or cloudless hours (\(n\)) [hour]

  • mn (float, np.array) – maximum possible duration of sunshine or daylight hours (\(N\)) [hour] See daylight_hours()

  • a_s (float or np.numpy) – regression constant (\(as\)) [-]. Default \(a_s=0.25\). It expresses the fraction of extraterrestrial radiation reaching the earth on overcast days (\(n = 0\))

  • b_s (float or np.numpy) – regression constant (\(bs\)) [-]. Default \(b_s=0.50\). The expression \(a_s+b_s\) indicates the fraction of extraterrestrial radiation reaching the earth on clear days (\(n = N\))

Returns

(float, np.array) daily total shortwave radiation (\(R_s\)) [MJ/m²day] reaching the earth

Note

If shortwave radiation (i.e., solar radiation) measurements are available, shortwave_radiation() function is no needed. Measurements of shortwave radiation may be directly used as input data in et0().

static net_shortwave_radiation(rs, albedo)

The net shortwave radiation (\(R_{ns}\)) [MJ/m²day] resulting from the balance between incoming and reflected solar radiation as defined in FAO 56 (eq. 38, p. 51):

\[R_{ns} = (1 − \alpha) R_s\]
Parameters
  • rs (float or np.array) – daily shortwave radiation (\(R_s\)) [MJ/m²day]. See shortwave_radiation()

  • albedo (float or np.array) – albedo or reflection coefficient (\(\alpha\) [-]). Range: \(0.0 \leq \alpha \leq 1.0\) (\(\alpha=0.23\) for the hypothetical grass reference crop). See PenmanMonteithDaily and et0()

Returns

(float or np.array) daily net shortwave radiation (\(R_{ns}\)) [MJ/m²day] reaching the earth

static net_longwave_radiation(t_min, t_max, rs, rs0, ea=None)

Return the net outgoing longwave radiation (\(R_{nl}\)) [MJ/m²day] as defined in FAO 56 (eq. 39, p. 52):

\[R_{nl} = \sigma\left[\frac{T_{max,K}^4 + T_{min,K}^4}{2}\right](0.34-0.14\sqrt{e_a})\left(1.35 \frac{R_s}{R_{so}}-0.35\right)\]
Parameters
  • t_min (float or np.array) – minimum daily air temperature (\(T_{max}\)) [°C]

  • t_max (float or np.array) – maximum daily air temperature (\(T_{min}\)) [°C]

  • rs (float or np.array) – shortwave radiation (\(R_s\)) [MJ/m²day]. See shortwave_radiation()

  • rs0 (float or np.array) – clear-sky shortwave radiation (\(R_{so}\)) [MJ/m²day]. See clear_sky_shortwave_radiation()

  • ea (float or np.array) – actual vapour pressure (\(e_a\)) [kPa]

Returns

(float or np.array) daily net outgoing longwave radiation (\(R_{nl}\)) [MJ/m²day]

Note

The \(R_s/R_{so}\) term in the equation above must be limited so that \(R_s/R_{so} \leq 1.0\).

et0(**kwargs)

Returns potential evapotranspiration (\(ETo\)) [mm/day] as described in FAO 56. Reference (grass) potencial evapotranspiration is returned for default constructor values. If values in **kwargs are arrays, their lengths must be the same.

Keyword Arguments
  • date (str, datetime.date, datetime.datetime, pandas.TimeStamp, or np.array)

  • doy (int or np.array) - day of the year (\(J\)) [-]. Range: \(1 \leq J \leq 366\). It is not used if date is given

  • u2 (float or np.array) - wind speed at 2 meters above ground surface [m/s]

  • uz (float or np.array) - measured wind speed at \(z\) meters above ground surface [m/s]

  • z (float or np.array) - height of measurement above ground surface [m]

  • t_mean (float or np.array) - daily mean air temperature [°C]

  • t_min (float or np.array) - daily minimum air temperature [°C]

  • t_max (float or np.array) - daily maximum air temperature [°C]

  • rh_mean (float or np.array) - daily mean relative humidity [%]

  • rh_min (float or np.array) - daily minimum relative humidity [%]

  • rh_max (float or np.array) - daily maximum relative humidity [%]

  • rs (float or np.array) - solar or shortwave radiation [MJ/m²day]

  • n (float or np.array) - daily actual duration of sunshine or cloudless hours [hour]

  • g (float or np.array) - soil heat flux density [MJ/m²day]. If not given, g defined in PenmanMonteithDaily() will be used

  • a_s (float or np.array) - see shortwave_radiation(). Default \(a_s = 0.25\)

  • b_s (float or np.array) - see shortwave_radiation(). Default \(b_s = 0.50\)

  • negative_rnl (bool) - allow negative net longwave radiation. Default \(negative\_rnl=True\)

  • negative_et0 (bool) - allow negative reference evapotranspiration. Default \(negative\_et0=True\)

Returns

(float or np.array) potential evapotranspiration (\(ETo\)) [mm/day]

Cases:

  • If date and doy are given, \(doy\) is disregarded

  • if \(uz\) is given, \(z\) must also be given

  • if \(u2\) and (\(uz\), \(z\)) are given, both \(uz\) and \(z\) are disregarded

  • if \(rs\) and \(n\) are given, \(n\) will be disregarded

  • The best options for air temperature are, in this order: 1) t_min, t_max, and t_mean, 2) t_min, t_max, and 3) tmean

  • The best options for relative air humidity are, in this order: 1) rh_max and rh_min, 2) rh_max, and 3) rh_mean

Example 1:

>>> from evapotranspiration.penman_monteith_daily import PenmanMonteithDaily
>>> pm = PenmanMonteithDaily(elevation=100, latitude=50.80)
>>> et0 = pm.et0(doy=187, u2=2.078, t_min=12.3, t_max=21.5, rh_min=63, rh_max=84, n=9.25)
>>> print(et0)
3.872968723753793

Example 2:

>>> from evapotranspiration.penman_monteith_daily import PenmanMonteithDaily
>>> pm = PenmanMonteithDaily(elevation=100, latitude=50.80)
>>> et0 = pm.et0(date='2001-07-06', u2=2.078, t_min=12.3, t_max=21.5, rh_min=63, rh_max=84, n=9.25)
>>> print(et0)
3.872968723753793

Example 3:

>>> from evapotranspiration.penman_monteith_daily import PenmanMonteithDaily
>>> pm = PenmanMonteithDaily(elevation=100, latitude=50.80)
>>> date=np.array(['2001-07-06', '2001-07-06'])
>>> u2=np.array([2.078, 2.078])
>>> t_min=np.array([12.3, 12.3])
>>> t_max=np.array([21.5, 21.5])
>>> rh_min=np.array([63, 63])
>>> rh_max=np.array([84, 84])
>>> n=np.array([9.25, 9.25])
>>> et0 = pm.et0(date=date, u2=u2, t_min=t_min, t_max=t_max, rh_min=rh_min, rh_max=rh_max, n=n)
>>> print(et0)
[3.87296872 3.87296872]
et0_frame(df, **kwargs)

Return the input DataFrame extended by et0() and further calculation parameters.

Parameters

df (pandas.DataFrame) – pandas DataFrame with columns corresponding to the inputs described in et0()

Keyword Arguments
  • show_all (bool) - show all results if \(True\), otherwise set parameter=True to show individual parameters. For example \(doy=True\), \(ld=True\), etc. See PenmanMonteithDaily()

Returns

(pandas.DataFrame) DataFrame

Module contents