Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
growing_degree_day_baskerville_emin.F90
Go to the documentation of this file.
2
3 use iso_c_binding
4
6 use datetime , only : mmdd2doy
7 use exceptions
9 use parameters
11 use fstring_list
12 implicit none
13
14 private
15
18
19 real (c_float), allocatable :: gdd_base(:)
20 real (c_float), allocatable :: gdd_max(:)
21 integer (c_int), allocatable :: gdd_reset_date(:)
22 type (t_netcdf4_file), pointer :: pncfile
23
24contains
25
26 subroutine growing_degree_day_be_initialize( is_cell_active, landuse_index )
27
28 logical (c_bool), intent(in) :: is_cell_active(:,:)
29 integer (c_int), intent(in) :: landuse_index(:)
30
31 ! [ LOCALS ]
32 integer (c_int) :: status
33 integer (c_int) :: indx
34 type (fstring_list_t) :: parameter_list
35 type (fstring_list_t) :: gdd_reset_val_list
36 character (len=32) :: sbuf
37 real (c_float), allocatable :: gdd_base_l(:)
38 real (c_float), allocatable :: gdd_max_l(:)
39 integer (c_int) :: number_of_landuse_codes
40 integer (c_int), allocatable :: landuse_code(:)
41
42 allocate( gdd_base( count( is_cell_active ) ), stat=status )
43 call assert( status == 0, "Problem allocating memory", __file__, __line__ )
44
45 allocate( gdd_max( count( is_cell_active ) ), stat=status )
46 call assert( status == 0, "Problem allocating memory", __file__, __line__ )
47
48 !> create string list that allows for alternate heading identifiers for the landuse code
49 call parameter_list%append("LU_Code")
50 call parameter_list%append("Landuse_Code")
51 call parameter_list%append("Landuse_Lookup_Code")
52
53 !> Determine how many landuse codes are present
54 call params%get_parameters( slkeys=parameter_list, ivalues=landuse_code )
55 number_of_landuse_codes = count( landuse_code >= 0 )
56 call parameter_list%clear()
57
58 call parameter_list%append("GDD_Base_Temp")
59 call parameter_list%append("GDD_Base_Temperature")
60 call parameter_list%append("GDD_Base")
61
62 call params%get_parameters( slkeys=parameter_list, fvalues=gdd_base_l )
63 call parameter_list%clear()
64
65 call parameter_list%append("GDD_Max_Temp")
66 call parameter_list%append("GDD_Maximum_Temperature")
67 call parameter_list%append("GDD_Maximum_Temp")
68 call parameter_list%append("GDD_Max")
69
70 call params%get_parameters( slkeys=parameter_list, fvalues=gdd_max_l )
71 call parameter_list%clear()
72
73 call parameter_list%append("GDD_Reset_Date")
74 call parameter_list%append("GDD_Reset")
75
76 call params%get_parameters( slkeys=parameter_list, slvalues=gdd_reset_val_list )
77 call parameter_list%clear()
78
79 allocate( gdd_reset_date( count( is_cell_active ) ), stat=status )
80 call assert( status==0, "Problem allocating memory.", __file__, __line__ )
81
82 if ( gdd_reset_val_list%count == number_of_landuse_codes &
83 .and. gdd_reset_val_list%count_matching("<NA>") == 0 ) then
84
85 ! retrieve gdd reset values; convert mm/dd to DOY
86 do indx=1, gdd_reset_val_list%count
87 sbuf = gdd_reset_val_list%get( indx )
88
89 where ( landuse_index == indx )
90 gdd_reset_date = mmdd2doy( sbuf, "GDD_RESET_DATE" )
91 end where
92
93 enddo
94
95 else
96
97 ! if no GDD_RESET_DATE found in parameter tables, assign the default: reset GDD at end of calendar year
98 gdd_reset_date = 365_c_int
99
100 endif
101
102
103 if ( ubound( gdd_max_l, 1 ) == number_of_landuse_codes &
104 .and. gdd_max_l(1) > rtinyval ) then
105
106 do indx=1, ubound( landuse_index, 1)
107 if( landuse_index( indx ) >= lbound( gdd_max, 1) .and. landuse_index( indx ) <= ubound( gdd_max, 1) ) then
108 gdd_max( indx ) = gdd_max_l( landuse_index( indx ) )
109 endif
110 enddo
111
112 else
113
114 ! if no GDD_MAX found in parameter tables, assign the default FAO-56 value
115 gdd_max = 86.0_c_float
116
117 endif
118
119 if ( ubound( gdd_base_l, 1 ) == number_of_landuse_codes &
120 .and. gdd_base_l(1) > rtinyval ) then
121
122 do indx=1, ubound( landuse_index, 1)
123 if( landuse_index( indx ) >= lbound( gdd_base, 1) .and. landuse_index( indx ) <= ubound( gdd_base, 1) ) then
124 gdd_base( indx ) = gdd_base_l( landuse_index( indx ) )
125 endif
126 enddo
127
128 else
129
130 ! if no GDD_Base found in parameter tables, assign the default FAO-56 value
131 gdd_base = 50.0_c_float
132
133 endif
134
136
137!--------------------------------------------------------------------------------------------------
138
139 subroutine growing_degree_day_be_calculate( gdd, tmean, tmin, tmax, order )
140
141 ! [ ARGUMENTS ]
142 real (c_float), intent(inout) :: gdd(:)
143 real (c_float), intent(in) :: tmean(:)
144 real (c_float), intent(in) :: tmin(:)
145 real (c_float), intent(in) :: tmax(:)
146 integer (c_int), intent(in) :: order(:)
147
148 ! [ LOCALS ]
149 real (c_float) :: delta
150 real (c_float) :: tmax_l
151 real (c_float) :: dd
152 real (c_float) :: w
153 real (c_float) :: at
154 real (c_float) :: a
155 integer (c_int) :: indx
156
157 do indx=lbound(order,1),ubound(order,1)
158
159 if ( sim_dt%iDOY == gdd_reset_date( order(indx) ) ) gdd(indx) = 0.0_c_float
160
161 tmax_l = min( tmax(indx), gdd_max( order(indx) ) )
162
163 if ( tmax_l <= gdd_base( order(indx) ) ) then
164
165 dd = 0.0_c_float
166
167 elseif ( tmin(indx) >= gdd_base( order(indx) ) ) then
168
169 dd = min(( tmean(indx) - gdd_base( order(indx) )), (gdd_max( order(indx) ) - gdd_base( order(indx) ) ) )
170
171 else
172
173 w = ( tmax_l - tmin(indx)) / 2.0_c_float
174
175 at = ( gdd_base( order(indx) ) - tmean(indx) ) / w
176
177 if ( at > 1 ) at = 1.0_c_float
178 if ( at < -1 ) at = -1._c_float
179
180 a = asin( at )
181
182 dd = ( ( w * cos( a ) ) - ( ( gdd_base( order(indx) ) - tmean(indx) ) &
183 * ( real( pi / 2._c_double, c_float ) - a ) ) ) / pi
184
185 endif
186
187 gdd(indx) = gdd(indx) + dd
188
189 enddo
190
192
193!--------------------------------------------------------------------------------------------------
194
This module contains physical constants and convenience functions aimed at performing unit conversion...
real(c_float), parameter, public rtinyval
real(c_double), parameter, public pi
This module contains the DATETIME_T class and associated time and date-related routines,...
Definition datetime.F90:9
integer(c_int) function, public mmdd2doy(smmdd, sinputitemname)
integer(c_int), dimension(:), allocatable, public gdd_reset_date
subroutine, public growing_degree_day_be_calculate(gdd, tmean, tmin, tmax, order)
real(c_float), dimension(:), allocatable, public gdd_base
subroutine, public growing_degree_day_be_initialize(is_cell_active, landuse_index)
real(c_float), dimension(:), allocatable, public gdd_max
Provide support for use of netCDF files as input for time-varying, gridded meteorlogic data,...
type(parameters_t), public params
type(date_range_t), public sim_dt