5.13.2-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_pet Module Reference

Module for calculating reference/potential evapotranspiration [mm d-1]. More...

Functions/Subroutines

elemental pure real(dp) function, public pet_hargreaves (harsamcoeff, harsamconst, tavg, tmax, tmin, latitude, doy)
 Reference Evapotranspiration after Hargreaves.
 
elemental pure real(dp) function, public pet_priestly (prietayparam, rn, tavg)
 Reference Evapotranspiration after Priestly-Taylor.
 
elemental pure real(dp) function, public pet_penman (net_rad, tavg, act_vap_pressure, aerodyn_resistance, bulksurface_resistance, a_s, a_sh)
 Reference Evapotranspiration after Penman-Monteith.
 
elemental pure real(dp) function, public extraterr_rad_approx (doy, latitude)
 Approximation of extraterrestrial radiation.
 
elemental pure real(dp) function, public slope_satpressure (tavg)
 slope of saturation vapour pressure curve
 
elemental pure real(dp) function, public sat_vap_pressure (tavg)
 calculation of the saturation vapour pressure
 

Detailed Description

Module for calculating reference/potential evapotranspiration [mm d-1].

This module calculates PET [mm/d] based on one of the methods

  • Hargreaves-Samani (1982)
  • Priestly-Taylor (1972)
  • Penman-Monteith FAO (1998)
    Authors
    Matthias Zink, Christoph Schneider, Matthias Cuntz
    Date
    Apr 2014

Function/Subroutine Documentation

◆ extraterr_rad_approx()

elemental pure real(dp) function, public mo_pet::extraterr_rad_approx ( integer(i4), intent(in)  doy,
real(dp), intent(in)  latitude 
)

Approximation of extraterrestrial radiation.

Approximation of extraterrestrial radiation at the top of the atmosphere \( R_a \) after Duffie and Beckman (1980). \( R_a \) is converted from \( [J\;m^{-2}\;d^{-1}] \) in \( [mm\;d^{-1}]\) .

\[ R_a = \frac{86400}{ \pi \cdot \lambda} \cdot E_0 \cdot d_r \cdot (\omega \cdot \sin(latitude) \cdot \sin(\delta) + \cos(latitude) \cdot \cos(\delta) \cdot \sin(\omega) \]

where \( E_0=1367\;J\;m^{-2}\;s^{-1} \) is the solar constant and It is dependent on the following sub equations: The relative distance Earth-Sun:

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

in which doy is the day of the year. The solar declination [radians] defined by

\[ \delta = 0.4093 \cdot \sin\left( \frac{2 \cdot \pi \cdot doy}{365} - 1.405 \right) \]

The sunset hour angle [radians]:

\[ \omega = \arccos( - \tan(latitude) * \tan(\delta) ) \]

< \( \lambda = 2.45 \cdot 10^6\;J\;m^{-2}\;mm^{-1} \) is the latent heat of vaporization.

Parameters
[in]integer(i4) :: doyday of year [-]
[in]real(dp) :: latitudelatitude [rad]
Returns
real(dp) :: extraterr_rad_approx — extraterrestrial radiation approximation \([W\;m^{-2}]\)
Authors
Matthias Zink
Date
Apr 2014

Definition at line 314 of file mo_pet.f90.

References mo_mhm_constants::duffiedelta1, mo_mhm_constants::duffiedelta2, mo_mhm_constants::duffiedr, and extraterr_rad_approx().

Referenced by extraterr_rad_approx(), and pet_hargreaves().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ pet_hargreaves()

elemental pure real(dp) function, public mo_pet::pet_hargreaves ( real(dp), intent(in)  harsamcoeff,
real(dp), intent(in)  harsamconst,
real(dp), intent(in)  tavg,
real(dp), intent(in)  tmax,
real(dp), intent(in)  tmin,
real(dp), intent(in)  latitude,
integer(i4), intent(in)  doy 
)

Reference Evapotranspiration after Hargreaves.

Calculates the Reference Evapotranspiration \( [mm\;d^{-1}] \) based on the Hargreaves-Samani (1982) model for a given cell by applying the equation

\[ PET = HarSamCoeff * R_a * (T_{avg} + HarSamConst) * \sqrt{ T_{max} - T_{min}} \]

where \( R_a\;[W\;m^{-2}]\) is the incoming solar radiation and \( T_{avg}, T_{max} \) and \( T_{min}\) \( [ ^{\circ}C]\) are the mean, maximum, and minimum daily temperatures at the given day, respectively.

Note
Hargreaves, G.H., and Samani, Z.A. (1982). "Estimating potential evapotranspiration." Tech. Note, J. Irrig. and drain. Engrg., ASCE, 108(3):225-230.
Parameters
[in]real(dp) :: HarSamCoeffcoefficient of Hargreaves-Samani equation [-]
[in]real(dp) :: HarSamConstconstant of Hargreaves-Samani equation [-]
[in]real(dp) :: tavgdaily men temperature \([^{\circ}C]\)
[in]real(dp) :: tmaxdaily maximum of temp \([^{\circ}C]\)
[in]real(dp) :: tmindaily minimum of temp \([^{\circ}C]\)
[in]real(dp) :: latitudelatitude of the cell for Ra estimation \([radians]\)
[in]integer(i4) :: doyday of year for Ra estimation
Returns
real(dp) :: pet_hargreaves — Hargreaves-Samani pot. evapotranspiration [mm d-1]
Authors
Matthias Zink
Date
Dec 2012

Definition at line 73 of file mo_pet.f90.

References extraterr_rad_approx(), and pet_hargreaves().

Referenced by mo_meteo_handler::get_corrected_pet(), and pet_hargreaves().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ pet_penman()

elemental pure real(dp) function, public mo_pet::pet_penman ( real(dp), intent(in)  net_rad,
real(dp), intent(in)  tavg,
real(dp), intent(in)  act_vap_pressure,
real(dp), intent(in)  aerodyn_resistance,
real(dp), intent(in)  bulksurface_resistance,
real(dp), intent(in)  a_s,
real(dp), intent(in)  a_sh 
)

Reference Evapotranspiration after Penman-Monteith.

Calculates the reference evapotranspiration \( [mm\;d^{-1}] \) based on the Penman-Monteith model for every given cell by applying the equation

\[ PET = \frac{1}{\lambda} \cdot \frac{\Delta \cdot R_n + \rho \cdot c_p \cdot (e_s-e) \cdot \frac{a_sh}{r_a}} {\Delta + \gamma \cdot \frac{a_sh}{a_s} \cdot \left( 1 + \frac{r_s}{r_a} \right) } \]

where \(R_n\;[W\;m^{-2}]\) is the net solar radiation, \(\Delta\;[kPa\;K^{-1}]\) is the slope of the saturation-vapour pressure curve, \( \lambda\;[MJ\;kg^{-1}] \) is the latent heat of vaporization, \( (e_s-e)\;[kPa] \) is the vapour pressure deficit of the air, \( \rho\;[kg\;m^{-3}] \) is the mean atmospheric density, \( c_p=1005.0\;J\;kg^{-1}\;K^{-1} \) is the specific heat of the air, \( \gamma [kPa\;K^{-1}] \) is the psychrometric constant, \( r_s [s m^{-1}] \) is the bulk canopy resistance, \( r_a [s m^{-1}] \) is the aerodynamic resistance, \( a_s [1] \) is the fraction of one-sided leaf area covered by stomata (1 if stomata are on one side only, 2 if they are on both sides) and \( a_{sh} [-] \) is the fraction of projected area exchanging sensible heat with the air (2) Implementation refers to the so-called Penman-Montheith equation for transpiration. Adjusting the arguments \( a_{sh} \) and \( a_s \) we obtain the corrected MU equation (for details see Schymanski and Or, 2017). If \( a_{sh} = 1 = a_s \) Penman-Montheith equation for transpiration is preserved. For reproducing characteristics of symmetrical amphistomatous leaves use \( a_{sh} = 2 = a_s \), in which case the classic PM equation is only missing a factor of 2 in the nominator, as pointed out by Jarvis and McNaughton (1986, Eq. A9). These analytical solutions eliminated the non-linearity problem of the saturation vapour pressure curve, but they do not consider the dependency of the long-wave component of the soil surface or leaf energy balance ( \( R_l \)) on soil or leaf temperature ( \( T_l \)). We assume that net radiation equals the absorbed short-wave radiation, i.e. \( R_N = R_s \) (p.79 in Monteith and Unsworth, 2013).

Parameters
[in]real(dp) :: net_radnet radiation \([W m^{-2}]\)
[in]real(dp) :: tavgaverage daily temperature \([^{\circ}C]\)
[in]real(dp) :: act_vap_pressureactual vapur pressure \([kPa]\)
[in]real(dp) :: aerodyn_resistanceaerodynmaical resistance \(s\;m^{-1}\)
[in]real(dp) :: bulksurface_resistancebulk surface resistance \(s\;m^{-1}\)
[in]real(dp) :: a_sfraction of one-sided leaf area covered by stomata \(1\)
[in]real(dp) :: a_shfraction of projected area exchanging sensible heat with the air \(1\)
Returns
real(dp) :: pet_penman — Reference Evapotranspiration [mm d-1]
Authors
Matthias Zink
Date
Apr 2014

Definition at line 235 of file mo_pet.f90.

References pet_penman(), sat_vap_pressure(), and slope_satpressure().

Referenced by mo_meteo_handler::get_corrected_pet(), and pet_penman().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ pet_priestly()

elemental pure real(dp) function, public mo_pet::pet_priestly ( real(dp), intent(in)  prietayparam,
real(dp), intent(in)  rn,
real(dp), intent(in)  tavg 
)

Reference Evapotranspiration after Priestly-Taylor.

Calculates the Reference Evapotranspiration \( [mm\;d^{-1}] \) based on the Priestly-Taylor (1972) model for every given cell by applying the equation

\[ PET = \alpha * \frac{\Delta}{(\gamma + \Delta)} * R_n \]

where \(R_n\;[W\;m^{-2}]\) is the net solar radiation \(\Delta = f(T_{avg})\) is the slope of the saturation-vapour pressure curve and \(\alpha\) is a emperical coefficient.

Parameters
[in]real(dp) :: PrieTayParamPriestley-Taylor coefficient \( \alpha [-] \)
[in]real(dp) :: Rnnet solar radiation \( [W\;m^{-2}] \)
[in]real(dp) :: Tavg
Returns
real(dp) :: pet_priestly — Priestley-Taylor pot. evapotranspiration [mm d-1]
Authors
Matthias Zink
Date
Apr 2014

Definition at line 149 of file mo_pet.f90.

References pet_priestly(), and slope_satpressure().

Referenced by mo_meteo_handler::get_corrected_pet(), and pet_priestly().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ sat_vap_pressure()

elemental pure real(dp) function, public mo_pet::sat_vap_pressure ( real(dp), intent(in)  tavg)

calculation of the saturation vapour pressure

Calculation of the saturation vapour pressure

\[ e_s(T_a) = 0.6108 \cdot \exp \left( \frac{17.27 \cdot T_a}{T_a + 237.3} \right) \]

Parameters
[in]real(dp) :: tavgtemperature [degC]
Returns
real(dp) :: sat_vap_pressure — saturation vapour pressure [kPa]
Authors
Matthias Zink
Date
Apr 2014

Definition at line 426 of file mo_pet.f90.

References sat_vap_pressure(), mo_mhm_constants::tetens_c1, mo_mhm_constants::tetens_c2, and mo_mhm_constants::tetens_c3.

Referenced by pet_penman(), sat_vap_pressure(), and slope_satpressure().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ slope_satpressure()

elemental pure real(dp) function, public mo_pet::slope_satpressure ( real(dp), intent(in)  tavg)

slope of saturation vapour pressure curve

slope of saturation vapour pressure curve after Tetens

\[ \Delta = \frac{0.6108 * e_s(T_a)}{e^(2 \cdot \log(T_a + 237.3))} \]

Parameters
[in]real(dp) :: tavgaverage daily temperature \([^{\circ}C]\)
Returns
real(dp) :: slope_satpressure — slope of saturation vapour pressure curve \([kPa\;K{-1}]\)
Authors
Matthias Zink
Date
Apr 2014

Definition at line 384 of file mo_pet.f90.

References sat_vap_pressure(), mo_mhm_constants::satpressureslope1, slope_satpressure(), and mo_mhm_constants::tetens_c3.

Referenced by mo_mrm_riv_temp_class::get_lat_heat(), pet_penman(), pet_priestly(), and slope_satpressure().

Here is the call graph for this function:
Here is the caller graph for this function: