Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
interception__bucket.F90
Go to the documentation of this file.
2
3 use iso_c_binding
5 use datetime, only : mmdd2doy
6 use exceptions
7 use parameters, only : params
9 use fstring
11 implicit none
12
13 private
14
18! public :: IS_GROWING_SEASON
19
20! logical (c_bool) :: GROWING_SEASON = .true._c_bool
21
22 integer (c_int), allocatable :: ilandusecodes(:)
23 real (c_float), allocatable :: interception_a_value_growing_season(:)
24 real (c_float), allocatable :: interception_b_value_growing_season(:)
25 real (c_float), allocatable :: interception_n_value_growing_season(:)
26 real (c_float), allocatable :: interception_a_value_nongrowing_season(:)
27 real (c_float), allocatable :: interception_b_value_nongrowing_season(:)
28 real (c_float), allocatable :: interception_n_value_nongrowing_season(:)
29 real (c_float), allocatable :: bucket_interception_storage_max_growing_season(:)
30 real (c_float), allocatable :: bucket_interception_storage_max_nongrowing_season(:)
31
32 !> Form of the bucket interception: I = A + P*B^n
33
34! logical (c_bool), allocatable :: IS_GROWING_SEASON(:)
35
36 character( len=2 ), parameter :: date_delims = "/-"
37 real (c_float), parameter :: nodata_value = -9999._c_float
38
39contains
40
41 subroutine interception_bucket_initialize( active_cells )
42
43 logical (c_bool), intent(in) :: active_cells(:,:)
44
45 ! [ LOCALS ]
46 integer (c_int) :: inumberoflanduses
47 logical (c_bool) :: larelengthsequal
48 character (len=:), allocatable :: stemp
49 type (fstring_list_t) :: sl_temp_list
50 type (fstring_list_t) :: sl_growing_season_begin
51 type (fstring_list_t) :: sl_growing_season_end
52 character (len=32) :: str_buffer
53 real (c_float), allocatable :: temp_values(:)
54 integer (c_int) :: indx
55 integer (c_int) :: status
56
57 !> Determine how many landuse codes are present
58 stemp = "LU_Code"
59 call params%get_parameters( skey=stemp, &
60 ivalues=ilandusecodes )
61
62 inumberoflanduses = count( ilandusecodes > 0 )
63
64 !> retrieve growing season interception amount: 'a' term
65 call sl_temp_list%clear()
66 call sl_temp_list%append("Growing_season_interception")
67 call sl_temp_list%append("Growing_season_interception_a")
68 call sl_temp_list%append("Interception_growing")
69 call sl_temp_list%append("Interception_a_term_growing_season")
70
71 call params%get_parameters( slkeys=sl_temp_list, &
72 fvalues=temp_values, &
73 lfatal=true )
74
75 if (all( temp_values > ftinyval) ) then
76
77 larelengthsequal = ( ubound(temp_values,1) == ubound(ilandusecodes,1) )
78
79 if ( .not. larelengthsequal ) &
80 call warn( smessage="The number of landuses does not match the number of interception values " &
81 //"specified for the 'a' term for use during the growing season.", &
82 smodule=__file__, iline=__line__, lfatal=true )
83
84 call move_alloc( temp_values, interception_a_value_growing_season )
85
86 endif
87
88 !> retrieve growing season interception amount: 'b' term
89 call sl_temp_list%clear()
90 call sl_temp_list%append("Growing_season_interception_b")
91 call sl_temp_list%append("Interception_growing_b_term")
92 call sl_temp_list%append("Interception_b_term_growing_season")
93
94 call params%get_parameters( slkeys=sl_temp_list, &
95 fvalues=temp_values, &
96 lfatal=false )
97
98 if (all( temp_values <= ftinyval) ) then
99
100 allocate(interception_b_value_growing_season( ubound( ilandusecodes, 1) ), stat=status )
102
103 else
104
105 larelengthsequal = ( ubound(temp_values,1) == ubound(ilandusecodes,1) )
106
107 if ( .not. larelengthsequal ) &
108 call warn( smessage="The number of landuses does not match the number of interception values " &
109 //"specified for the 'b' term for use during the growing season.", &
110 smodule=__file__, iline=__line__, lfatal=false )
111
112 call move_alloc( temp_values, interception_b_value_growing_season )
113
114 endif
115
116 !> retrieve growing season interception amount: 'n' term
117 call sl_temp_list%clear()
118 call sl_temp_list%append("Growing_season_interception_n")
119 call sl_temp_list%append("Interception_growing_n_term")
120 call sl_temp_list%append("Interception_n_term_growing_season")
121
122 call params%get_parameters( slkeys=sl_temp_list, &
123 fvalues=temp_values, &
124 lfatal=false )
125
126 if (all( temp_values <= ftinyval) ) then
127
128 allocate(interception_n_value_growing_season( ubound( ilandusecodes, 1) ), stat=status )
130
131 else
132
133 larelengthsequal = ( ubound(temp_values,1) == ubound(ilandusecodes,1) )
134
135 if ( .not. larelengthsequal ) &
136 call warn( smessage="The number of landuses does not match the number of interception values " &
137 //"specified for the 'n' term for use during the growing season.", &
138 smodule=__file__, iline=__line__, lfatal=true )
139
140 call move_alloc( temp_values, interception_n_value_growing_season )
141
142 endif
143
144 !> retrieve nongrowing season interception amount: 'a' term
145 call sl_temp_list%clear()
146 call sl_temp_list%append("Nongrowing_season_interception")
147 call sl_temp_list%append("Nongrowing_season_interception_a")
148 call sl_temp_list%append("Interception_nongrowing")
149 call sl_temp_list%append("Interception_a_term_nongrowing_season")
150
151 call params%get_parameters( slkeys=sl_temp_list, &
152 fvalues=temp_values, &
153 lfatal=true )
154
155 if (all( temp_values > ftinyval) ) then
156
157 larelengthsequal = ( ubound(temp_values,1) == ubound(ilandusecodes,1) )
158
159 if ( .not. larelengthsequal ) &
160 call warn( smessage="The number of landuses does not match the number of interception values " &
161 //"specified for the 'a' term for use during the nongrowing season.", &
162 smodule=__file__, iline=__line__, lfatal=false )
163
164 call move_alloc( temp_values, interception_a_value_nongrowing_season )
165
166 endif
167
168 !> retrieve nongrowing season interception amount: 'b' term
169 call sl_temp_list%clear()
170 call sl_temp_list%append("Nongrowing_season_interception_b")
171 call sl_temp_list%append("Interception_nongrowing_b_term")
172 call sl_temp_list%append("Interception_b_term_nongrowing_season")
173
174 call params%get_parameters( slkeys=sl_temp_list, &
175 fvalues=temp_values, &
176 lfatal=false )
177
178 if (all( temp_values <= ftinyval) ) then
179
180 allocate(interception_b_value_nongrowing_season( ubound( ilandusecodes, 1) ), stat=status )
182
183 else
184
185 larelengthsequal = ( ubound(temp_values,1) == ubound(ilandusecodes,1) )
186
187 if ( .not. larelengthsequal ) &
188 call warn( smessage="The number of landuses does not match the number of interception values " &
189 //"specified for the 'b' term for use during the nongrowing season.", &
190 smodule=__file__, iline=__line__, lfatal=false )
191
192 call move_alloc( temp_values, interception_b_value_nongrowing_season )
193
194 endif
195
196 !> retrieve nongrowing season interception amount: 'n' term
197 call sl_temp_list%clear()
198 call sl_temp_list%append("Nongrowing_season_interception_n")
199 call sl_temp_list%append("Interception_nongrowing_n_term")
200 call sl_temp_list%append("Interception_n_term_nongrowing_season")
201
202 call params%get_parameters( slkeys=sl_temp_list, &
203 fvalues=temp_values, &
204 lfatal=false )
205
206 if (all( temp_values <= ftinyval) ) then
207
208 allocate(interception_n_value_nongrowing_season( ubound( ilandusecodes, 1) ), stat=status )
210
211 else
212
213 larelengthsequal = ( ubound(temp_values,1) == ubound(ilandusecodes,1) )
214
215 if ( .not. larelengthsequal ) &
216 call warn( smessage="The number of landuses does not match the number of interception values " &
217 //"specified for the 'n' term for use during the nongrowing season.", &
218 smodule=__file__, iline=__line__, lfatal=false )
219
220 call move_alloc( temp_values, interception_n_value_nongrowing_season )
221
222 endif
223
224 !> retrieve interception storage max (NONGROWING)
225 call sl_temp_list%clear()
226 call sl_temp_list%append("Interception_storage_max_nongrowing")
227 call sl_temp_list%append("Interception_storage_max_nongrowing_season")
228 call sl_temp_list%append("Interception_Storage_Maximum_nongrowing")
229
230 call params%get_parameters( slkeys=sl_temp_list, &
232 lfatal=false )
233
234 larelengthsequal = ( ubound(bucket_interception_storage_max_nongrowing_season,1) == ubound(ilandusecodes,1) )
235
236 if ( .not. larelengthsequal ) then
237 call warn( smessage="The number of landuses does not match the number of interception storage " &
238 //"maximum values for the NONGROWING season " &
239 //"('interception_storage_max_nongrowing').", &
240 shints="A default value equal to the 'Growing_season_interception_a' was assigned for the" &
241 //" maximum interception storage", &
242 smodule=__file__, iline=__line__, lfatal=false )
243 allocate(temp_values(inumberoflanduses), stat=status)
245 call move_alloc(temp_values, bucket_interception_storage_max_nongrowing_season)
246 endif
247
248 !> retrieve interception storage max (GROWING)
249 call sl_temp_list%clear()
250 call sl_temp_list%append("Interception_storage_max_growing")
251 call sl_temp_list%append("Interception_storage_max_growing_season")
252 call sl_temp_list%append("Interception_Storage_Maximum_growing")
253
254 call params%get_parameters( slkeys=sl_temp_list, &
256 lfatal=false )
257
258 larelengthsequal = ( ubound(bucket_interception_storage_max_growing_season,1) == ubound(ilandusecodes,1) )
259
260 if ( .not. larelengthsequal ) then
261 call warn( smessage="The number of landuses does not match the number of interception storage " &
262 //"maximum values for the GROWING season " &
263 //"('interception_storage_max_growing').", &
264 shints="A default value equal to the 'Nongrowing_season_interception_a' was" &
265 //" assigned for the maximum interception storage", &
266 smodule=__file__, iline=__line__, lfatal=false )
267 allocate(temp_values(inumberoflanduses), stat=status)
269 call move_alloc(temp_values, bucket_interception_storage_max_growing_season)
270 endif
271
272 end subroutine interception_bucket_initialize
273
274!--------------------------------------------------------------------------------------------------
275
276 elemental subroutine interception_bucket_calculate( iLanduseIndex, &
277 fPrecip, &
278 fFog, &
279 fCanopy_Cover_Fraction, &
280 it_is_growing_season, &
281 fInterception )
282
283 integer (c_int), intent(in) :: ilanduseindex
284 real (c_float), intent(in) :: fprecip
285 real (c_float), intent(in) :: ffog
286 real (c_float), intent(in) :: fcanopy_cover_fraction
287 logical (c_bool), intent(in) :: it_is_growing_season
288 real (c_float), intent(out) :: finterception
289
290 ! [ LOCALS ]
291 real (c_float) :: fpotentialinterception
292 real (c_float) :: precip_plus_fog
293
294 precip_plus_fog = fprecip + ffog
295
296 if ( it_is_growing_season ) then
297
298 fpotentialinterception = interception_a_value_growing_season( ilanduseindex ) ! &
299! + INTERCEPTION_B_VALUE_GROWING_SEASON( iLanduseIndex ) &
300! * precip_plus_fog &
301! ** INTERCEPTION_N_VALUE_GROWING_SEASON( iLanduseIndex )
302
303 else
304
305 fpotentialinterception = interception_a_value_nongrowing_season( ilanduseindex ) ! &
306 ! + INTERCEPTION_B_VALUE_NONGROWING_SEASON( iLanduseIndex ) &
307 ! * precip_plus_fog &
308 ! ** INTERCEPTION_N_VALUE_NONGROWING_SEASON( iLanduseIndex )
309
310 endif
311
312 finterception = min( fpotentialinterception, fprecip + ffog ) * fcanopy_cover_fraction
313
314 end subroutine interception_bucket_calculate
315
316end module interception__bucket
This module contains physical constants and convenience functions aimed at performing unit conversion...
logical(c_bool), parameter, public true
logical(c_bool), parameter, public false
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)
subroutine, public warn(smessage, smodule, iline, shints, lfatal, iloglevel, lecho)
real(c_float), dimension(:), allocatable, public bucket_interception_storage_max_growing_season
real(c_float), parameter nodata_value
real(c_float), dimension(:), allocatable interception_n_value_growing_season
character(len=2), parameter date_delims
Form of the bucket interception: I = A + P*B^n.
real(c_float), dimension(:), allocatable, public bucket_interception_storage_max_nongrowing_season
integer(c_int), dimension(:), allocatable ilandusecodes
elemental subroutine, public interception_bucket_calculate(ilanduseindex, fprecip, ffog, fcanopy_cover_fraction, it_is_growing_season, finterception)
real(c_float), dimension(:), allocatable interception_b_value_growing_season
real(c_float), dimension(:), allocatable interception_b_value_nongrowing_season
subroutine, public interception_bucket_initialize(active_cells)
real(c_float), dimension(:), allocatable interception_a_value_growing_season
real(c_float), dimension(:), allocatable interception_a_value_nongrowing_season
real(c_float), dimension(:), allocatable interception_n_value_nongrowing_season
type(parameters_t), public params
type(date_range_t), public sim_dt