Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
et__jensen_haise.F90
Go to the documentation of this file.
1!> @file
2!> Contains a single module, et_jensen_haise, which
3!> calculates potential evapotranspiration by means of the Jensen-Haise (1963) method.
4
5
6!> Calculates potential evapotranspiration by means of the
7!> Jensen-Haise (1963) method.
9!!****h* SWB/et_jensen_haise
10! NAME
11! et_jensen_haise.f95 - Evapotranspiration calculation using the
12! Jensen-Haise method.
13!
14! SYNOPSIS
15! This module calculates evapotranspiration using the Jensen-Haise
16! method.
17!
18! NOTES
19! Original method is documented in:
20!
21! Jensen, M.E., H.R. Haise. 1963. Estimating evapotranspiration from solar
22! radiation. Journal of Irrigation and Drainage Engineering 89(IR4):15-41.
23!
24!!***
25
26 use iso_c_binding, only : c_short, c_int, c_float, c_double
35
36 implicit none
37
38 private
39
40 public :: et_jh_calculate
41
42contains
43
44elemental function et_jh_calculate( iDayOfYear, iNumDaysInYear, fLatitude, fTMin, fTMax, &
45 fAs, fBs, fSunPct ) result(fReferenceET0)
46
47 integer (c_int), intent(in) :: idayofyear
48 integer (c_int), intent(in) :: inumdaysinyear
49 real (c_float), intent(in) :: flatitude
50 real (c_float), intent(in) :: ftmin
51 real (c_float), intent(in) :: ftmax
52 real (c_float), intent(in), optional :: fas
53 real (c_float), intent(in), optional :: fbs
54 real (c_float), intent(in), optional :: fsunpct
55 real (c_double) :: freferenceet0
56
57 ! [ LOCALS ]
58 real (c_double) :: dra
59 real (c_double) :: ddelta
60 real (c_double) :: domega_s
61 real (c_double) :: dd_r
62 real (c_double) :: drs
63 real (c_float) :: ftavg
64
65 real (c_double) :: das_l
66 real (c_double) :: dbs_l
67 real (c_double) :: dsunpct_l
68
69 if (present( fas ) ) then
70 das_l = fas
71 else
72 das_l = 0.25_c_double
73 endif
74
75 if (present( fbs ) ) then
76 dbs_l = fbs
77 else
78 dbs_l = 0.5_c_double
79 endif
80
81 if (present( fsunpct ) ) then
82 dsunpct_l = fsunpct
83 else
84 dsunpct_l = estimate_percent_of_possible_sunshine__psun(ftmax=f_to_k( ftmax ), ftmin=f_to_k( ftmin ) )
85 endif
86
87
88 ftavg = ( ftmin + ftmax ) / 2_c_float
89
90! fD_r = 1.0_c_float + 0.033_c_float * cos( TWOPI * iDayOfYear / iNumDaysInYear )
91
92 dd_r = relative_earth_sun_distance__d_r( idayofyear=idayofyear, inumdaysinyear=inumdaysinyear )
93
94! fDelta = 0.4093_c_float * sin( (TWOPI * iDayOfYear / iNumDaysInYear ) - 1.405_c_float )
95
96 ddelta = solar_declination_simple__delta( idayofyear=idayofyear, inumdaysinyear=inumdaysinyear )
97
98! fOmega_s = acos( -tan(fLatitude) * tan(fDelta) )
99
100 domega_s = sunrise_sunset_angle__omega_s( real( flatitude, c_double), ddelta )
101
102! dSo = 2.44722_c_float * 15.392_c_float * dD_r * ( dOmega_s * sin(dLatitude) * sin(dDelta) + &
103! sin(dOmega_s) * cos(dLatitude) * cos(dDelta) )
104
105 dra = extraterrestrial_radiation__ra(dlatitude=real( flatitude, c_double), &
106 ddelta=ddelta, &
107 domega_s=domega_s, &
108 ddsubr=dd_r )
109
110! dSn = dSo * ( 1.0_c_float - fAlbedo_l ) * ( fAs_l + fBs_l * fSunPct_l / 100.0_c_float )
111
112 drs = solar_radiation__rs(dra=dra, das=das_l, dbs=dbs_l, fpctsun=dsunpct_l )
113
114 if ( ftavg <= rfreezing ) then
115 freferenceet0 = 0.0_c_float
116 else
117 ! 'equivalent_evaporation' yields the depth of water that could be evaporated for the
118 ! given solar radiation in mm
119 freferenceet0 = ( 0.014_c_float * ftavg - 0.37_c_float ) * equivalent_evaporation( drs ) / mm_per_in
120! fReferenceET0 = UNIT_CONV * ( 0.025_c_float * fTAvg + 0.078_c_float ) * fSn
121 end if
122
123!> @todo Check these equations against original paper.
124
125!! original paper gives ET0 as:
126
127!! Etp = ( 0.014 * Tavg - 0.37 ) * Rs
128!
129! where Rs is the solar radiation in inches per day
130
131end function et_jh_calculate
132
133end module et__jensen_haise
This module contains physical constants and convenience functions aimed at performing unit conversion...
real(c_double), parameter, public mm_per_in
real(c_float), parameter, public rfreezing
Calculates potential evapotranspiration by means of the Jensen-Haise (1963) method.
elemental real(c_double) function, public et_jh_calculate(idayofyear, inumdaysinyear, flatitude, ftmin, ftmax, fas, fbs, fsunpct)
elemental real(c_double) function equivalent_evaporation(rr)
elemental real(c_double) function solar_declination_simple__delta(idayofyear, inumdaysinyear)
Calculate the solar declination for a given day of the year.
elemental real(c_float) function estimate_percent_of_possible_sunshine__psun(ftmax, ftmin)
Estimate percent of possible sunshine.
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.
elemental real(c_double) function solar_radiation__rs(dra, das, dbs, fpctsun)
Calculate solar radiation by means of the Angstrom formula.