Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
et__hargreaves_samani.F90
Go to the documentation of this file.
1!> @file
2!> Contains a single module, et_hargreaves, which
3!> calculates potential evapotranspiration by means of the Hargreaves-Samani (1985) method.
4
5
6!> Calculates potential evapotranspiration by means of the
7!> Hargreaves-Samani (1985) method.
9!!****h* SWB/et_hargreaves
10! NAME
11! et_hargreaves.f95 - Evapotranspiration calculation using the
12! Hargreaves method.
13!
14! SYNOPSIS
15! This module calculates evapotranspiration using the Hargreaves
16! method.
17!
18! NOTES
19!
20! Reference:
21!
22! Allen, R.G., and others, 2006, FAO Irrigation and Drainage Paper No. 56,
23! "Crop Evapotranspiration (Guidelines for computing crop water
24! requirements)", Food and Agriculture Organization, Rome, Italy.
25!
26!!***
27
28 use iso_c_binding, only : c_short, c_int, c_float, c_double
31 use parameters, only : params
33
34 implicit none
35
36
37 ! ET parameters -- default values are from Hargreaves and Samani (1985)
38 real (c_double) :: et_slope ! = 0.0023
39 real (c_double) :: et_exponent ! = 0.5
40 real (c_double) :: et_constant ! = 17.8
41
42contains
43
44subroutine et_hargreaves_initialize( ) !pConfig, sRecord )
45 !! Configures the module, using the command line in 'sRecord'
46 ! [ ARGUMENTS ]
47! type (T_MODEL_CONFIGURATION), pointer :: pConfig ! pointer to data structure that contains
48 ! model options, flags, and other settings
49! character (len=*),intent(inout) :: sRecord
50
51 ! [ LOCALS ]
52 character (len=256) :: sOption
53 integer (c_int) :: iStat
54 real (c_float) :: rValue
55
56 real (kind=c_float), allocatable :: fet_slope(:)
57 real (kind=c_float), allocatable :: fet_exponent(:)
58 real (kind=c_float), allocatable :: fet_constant(:)
59
60 call params%get_parameters( skey="Hargreaves_ET_slope", fvalues=fet_slope)
61 call params%get_parameters( skey="Hargreaves_ET_exponent", fvalues=fet_exponent)
62 call params%get_parameters( skey="Hargreaves_ET_constant", fvalues=fet_constant)
63
64 if ( fet_slope(1) > ftinyval) then
65 et_slope = real(fet_slope(1), c_double)
66 else
67 et_slope = 0.0023_c_double
68 endif
69
70 if ( fet_exponent(1) > ftinyval) then
71 et_exponent = real(fet_exponent(1), c_double)
72 else
73 et_exponent = 0.5_c_double
74 endif
75
76 if ( fet_constant(1) > ftinyval) then
77 et_constant = real(fet_constant(1), c_double)
78 else
79 et_constant = 17.8_c_double
80 endif
81
82 print *, "****************************************************************************"
83 print *, "initializing Hargreaves ET parameters"
84 print *, "****************************************************************************"
85 print *, " ET_SLOPE = ", et_slope
86 print *, " ET_EXPONENT = ", et_exponent
87 print *, " ET_CONSTANT = ", et_constant
88
89 ! write(UNIT=LU_LOG,FMT=*) "Configuring Hargreaves PET model"
90
91! if (pConfig%rSouthernLatitude <= rNO_DATA_NCDC &
92! .or. pConfig%rNorthernLatitude <= rNO_DATA_NCDC) then
93
94! call Chomp( sRecord,sOption )
95! read ( unit=sOption, fmt=*, iostat=iStat ) rValue
96! call Assert( iStat == 0, "Could not read the southerly latitude" )
97! pConfig%rSouthernLatitude = dpTWOPI * rValue / 360.0_c_float
98
99! call Chomp( sRecord,sOption )
100! read ( unit=sOption, fmt=*, iostat=iStat ) rValue
101! call Assert( iStat == 0, "Could not read the northerly latitude" )
102! pConfig%rNorthernLatitude = dpTWOPI * rValue / 360.0_c_float
103
104! else
105
106! call echolog("Southern and northern latitude values have been determined" &
107! //"~from project grid bounds and projection parameters. The values supplied" &
108! //"~with the Hargreaves PET option will be ignored...")
109
110! endif
111
112end subroutine et_hargreaves_initialize
113
114!------------------------------------------------------------------------------
115
116impure elemental function et_hargreaves_calculate( iDayOfYear, iNumDaysInYear, fLatitude, fTMin, fTMax ) result(fReferenceET0)
117 !! Computes the potential ET for each cell, based on TMIN and TMAX.
118 !! Stores cell-by-cell PET values in the model grid.
119
120 integer (c_int),intent(in) :: idayofyear
121 integer (c_int),intent(in) :: inumdaysinyear
122 real (c_float), intent(in) :: flatitude
123 real (c_float), intent(in) :: ftmin
124 real (c_float), intent(in) :: ftmax
125 real (c_double) :: freferenceet0
126
127 ! [ LOCALS ]
128 real (c_double) :: fdelta, fomega_s, fd_r, fra
129 real (c_double) :: dlatitude_radians
130
131 dlatitude_radians = flatitude * degrees_to_radians
132
133 fd_r =relative_earth_sun_distance__d_r(idayofyear,inumdaysinyear)
134 fdelta = solar_declination__delta(idayofyear, inumdaysinyear)
135
136 fomega_s = sunrise_sunset_angle__omega_s(dlatitude_radians, fdelta)
137
138 ! NOTE that the following equation returns extraterrestrial radiation in
139 ! MJ / m**2 / day. The Hargreaves equation requires extraterrestrial
140 ! radiation to be expressed in units of mm / day.
141 fra = extraterrestrial_radiation__ra(dlatitude_radians, fdelta, fomega_s, fd_r)
142
143 freferenceet0 = et0_hargreaves( equivalent_evaporation(fra), ftmin, ftmax)
144
145end function et_hargreaves_calculate
146
147
148elemental function et0_hargreaves( rRa, rTMinF, rTMaxF ) result(rET_0)
149
150 ! [ ARGUMENTS ]
151 real (c_double),intent(in) :: rra
152 real (c_float),intent(in) :: rtminf
153 real (c_float),intent(in) :: rtmaxf
154
155 ! [ RETURN VALUE ]
156 real (c_double) :: ret_0
157
158 ! [ LOCALS ]
159 real (c_double) :: rtdelta
160 real (c_double) :: rtavg
161
162 rtavg = (rtminf + rtmaxf) / 2.0_c_double
163
164 rtdelta = abs(f_to_k(real(rtmaxf, c_double)) - f_to_k(real(rtminf, c_double)))
165
166 ret_0 = max(rzero, &
167 mm_to_in( et_slope * rra * (f_to_c(rtavg) + et_constant) * (rtdelta**et_exponent) ) )
168! mm_to_in( 0.0023_c_float * rRa * (F_to_C(rTavg) + 17.8_c_float) * sqrt(rTDelta)) )
169
170end function et0_hargreaves
171
172
173end module et__hargreaves_samani
This module contains physical constants and convenience functions aimed at performing unit conversion...
real(c_float), parameter, public rzero
real(c_double), parameter, public degrees_to_radians
Calculates potential evapotranspiration by means of the Hargreaves-Samani (1985) method.
impure elemental real(c_double) function et_hargreaves_calculate(idayofyear, inumdaysinyear, flatitude, ftmin, ftmax)
elemental real(c_double) function et0_hargreaves(rra, rtminf, rtmaxf)
elemental real(c_double) function equivalent_evaporation(rr)
type(parameters_t), public params
elemental real(c_double) function solar_declination__delta(idayofyear, inumdaysinyear)
Calculate the solar declination for a given day of the year.
elemental real(c_double) function relative_earth_sun_distance__d_r(idayofyear, inumdaysinyear)
Calculate the inverse relative Earth-Sun distance for a given day of the year.
elemental real(c_double) function sunrise_sunset_angle__omega_s(dlatitude, ddelta)
Calculate sunrise/sunset angle, in RADIANS.
elemental real(c_double) function extraterrestrial_radiation__ra(dlatitude, ddelta, domega_s, ddsubr)
Calculate extraterrestrial radiation given latitude and time of year.