Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
crop_coefficients__fao56.F90
Go to the documentation of this file.
1!> @file
2!> Contains a single module, \ref crop_coefficients__fao56, which
3!> provides support for modifying reference ET through the use of
4!> crop coefficients
5
6!> Update crop coefficients for crop types in simulation.
7
9
10 use iso_c_binding, only : c_bool, c_short, c_int, c_float, c_double
14 use data_catalog, only : dat
16 use datetime
17 use logfiles, only : logs, log_all
18 use exceptions, only : assert, warn, die
19 use parameters, only : params
20 use simulation_datetime, only : sim_dt
21 use fstring, only : ascharacter, squote, operator(.contains.)
22 use fstring_list
23 implicit none
24
25 private
26
33 public :: kcb_min, kcb_ini, kcb_mid, kcb_end
36
37 enum, bind(c)
39 end enum
40
41 enum, bind(c)
43 end enum
44
45 enum, bind(c)
48 end enum
49
50 enum, bind(c)
51 enumerator :: kcb_ini=13, kcb_mid, kcb_end, kcb_min
52 end enum
53
54 enum, bind(c)
55 enumerator :: jan = 1, feb, mar, apr, may, jun, jul, aug, sep, oct, nov, dec
56 end enum
57
58 enum, bind(c)
60 end enum
61
62 ! Private, module level variables
63 ! kept at a landuse code level (i.e. same value applies to all cells with same LU codes)
64 integer (c_int), allocatable :: landuse_code(:)
65! real (c_float), allocatable :: REW(:,:)
66! real (c_float), allocatable :: TEW(:,:)
67 real (c_float), allocatable :: kcb_l(:,:)
68 integer (c_int), allocatable :: kcb_method(:)
69 real (c_float), allocatable :: growth_stage_shift_days(:)
70 real (c_float), allocatable :: growth_stage_length_in_days(:,:)
71 real (c_float), allocatable :: growth_stage_gdd(:,:)
72 type (datetime_t), allocatable :: growth_stage_date(:,:)
73
74 integer (c_int) :: lu_soils_csv
75
76contains
77
79
80 ! [ LOCALS ]
81 ! type (FSTRING_LIST_T) :: slREW, slTEW
82 type (fstring_list_t) :: sllist
83 type (datetime_t) :: dt
84 type (datetime_t) :: temp_date
85 ! integer (c_int), allocatable :: iTEWSeqNums(:)
86 ! integer (c_int), allocatable :: iREWSeqNums(:)
87 integer (c_int) :: inumberoftew, inumberofrew
88 integer (c_int) :: inumberoflanduses
89 integer (c_int) :: iindex, iindex2
90 integer (c_int) :: istat
91 real (c_float) :: growing_cycle_length
92
93 character (len=10) :: smmddyyyy
94 character (len=:), allocatable :: stext
95
96 type (fstring_list_t) :: slplantingdate
97 type (datetime_t) :: dtplantingdate
98 character (len=:), allocatable :: plantingdate_str
99
100 real (c_float), allocatable :: l_ini_l(:)
101 real (c_float), allocatable :: l_dev_l(:)
102 real (c_float), allocatable :: l_mid_l(:)
103 real (c_float), allocatable :: l_late_l(:)
104 real (c_float), allocatable :: l_fallow_l(:)
105
106 real (c_float), allocatable :: gdd_plant_l(:)
107 real (c_float), allocatable :: gdd_ini_l(:)
108 real (c_float), allocatable :: gdd_dev_l(:)
109 real (c_float), allocatable :: gdd_mid_l(:)
110 real (c_float), allocatable :: gdd_late_l(:)
111
112 real (c_float), allocatable :: kcb_max(:)
113
114 real (c_float), allocatable :: kcb_ini_l(:)
115 real (c_float), allocatable :: kcb_mid_l(:)
116 real (c_float), allocatable :: kcb_end_l(:)
117 real (c_float), allocatable :: kcb_min_l(:)
118
119 real (c_float), allocatable :: kcb_jan(:)
120 real (c_float), allocatable :: kcb_feb(:)
121 real (c_float), allocatable :: kcb_mar(:)
122 real (c_float), allocatable :: kcb_apr(:)
123 real (c_float), allocatable :: kcb_may(:)
124 real (c_float), allocatable :: kcb_jun(:)
125 real (c_float), allocatable :: kcb_jul(:)
126 real (c_float), allocatable :: kcb_aug(:)
127 real (c_float), allocatable :: kcb_sep(:)
128 real (c_float), allocatable :: kcb_oct(:)
129 real (c_float), allocatable :: kcb_nov(:)
130 real (c_float), allocatable :: kcb_dec(:)
131
132 real (c_float) :: fkcb_initial
133 real (c_float) :: frz_initial
134
135 real (c_float), parameter :: near_zero = 1.0e-9_c_float
136
137 type (data_catalog_entry_t), pointer :: pinitial_percent_soil_moisture
138
139 !> create string list that allows for alternate heading identifiers for the landuse code
140 sllist = create_list("LU_Code, Landuse_Code, Landuse_Lookup_Code")
141
142 !> Determine how many landuse codes are present
143 call params%get_parameters( slkeys=sllist, ivalues=landuse_code )
144 inumberoflanduses = count( landuse_code >= 0 )
145 !> @todo Implement thorough input error checking:
146 !! are all soils in grid included in table values?
147 !> is soil suffix vector continuous?
148
149 ! Retrieve and populate the Readily Evaporable Water (REW) table values
150 ! CALL PARAMS%get_parameters( fValues=REW, sPrefix="REW_", iNumRows=iNumberOfLanduses )
151
152 ! Retrieve and populate the Total Evaporable Water (TEW) table values
153 ! CALL PARAMS%get_parameters( fValues=TEW, sPrefix="TEW_", iNumRows=iNumberOfLanduses )
154
155 !> @TODO What should happen if the TEW / REW header entries do *not* fall in a
156 !! logical sequence of values? In other words, if the user has columns named
157 !! REW_1, REW_3, REW_5, only the values associated with "REW_1" would be retrieved.
158 !! Needless to say, this would be catastrophic.
159
160 call params%get_parameters( skey="Planting_date", slvalues=slplantingdate )
161
162 call params%get_parameters( skey="L_ini", fvalues=l_ini_l)
163 call params%get_parameters( skey="L_dev", fvalues=l_dev_l)
164 call params%get_parameters( skey="L_mid", fvalues=l_mid_l)
165 call params%get_parameters( skey="L_late", fvalues=l_late_l)
166 call params%get_parameters( skey="L_fallow", fvalues=l_fallow_l)
167
168 call params%get_parameters( skey="GDD_plant", fvalues=gdd_plant_l)
169 call params%get_parameters( skey="GDD_ini", fvalues=gdd_ini_l)
170 call params%get_parameters( skey="GDD_dev", fvalues=gdd_dev_l)
171 call params%get_parameters( skey="GDD_mid", fvalues=gdd_mid_l)
172 call params%get_parameters( skey="GDD_late", fvalues=gdd_late_l)
173
174 call params%get_parameters( skey="Kcb_ini", fvalues=kcb_ini_l)
175 call params%get_parameters( skey="Kcb_mid", fvalues=kcb_mid_l)
176 call params%get_parameters( skey="Kcb_end", fvalues=kcb_end_l)
177 call params%get_parameters( skey="Kcb_min", fvalues=kcb_min_l)
178
179 call params%get_parameters( skey="Kcb_Jan", fvalues=kcb_jan )
180 call params%get_parameters( skey="Kcb_Feb", fvalues=kcb_feb )
181 call params%get_parameters( skey="Kcb_Mar", fvalues=kcb_mar )
182 call params%get_parameters( skey="Kcb_Apr", fvalues=kcb_apr )
183 call params%get_parameters( skey="Kcb_May", fvalues=kcb_may )
184 call params%get_parameters( skey="Kcb_Jun", fvalues=kcb_jun )
185 call params%get_parameters( skey="Kcb_Jul", fvalues=kcb_jul )
186 call params%get_parameters( skey="Kcb_Aug", fvalues=kcb_aug )
187 call params%get_parameters( skey="Kcb_Sep", fvalues=kcb_sep )
188 call params%get_parameters( skey="Kcb_Oct", fvalues=kcb_oct )
189 call params%get_parameters( skey="Kcb_Nov", fvalues=kcb_nov )
190 call params%get_parameters( skey="Kcb_Dec", fvalues=kcb_dec )
191
192
193 allocate( growth_stage_length_in_days( 5, inumberoflanduses ), stat=istat )
194 call assert( istat==0, "Failed to allocate memory for GROWTH_STAGE_LENGTH_IN_DAYS array", &
195 __file__, __line__ )
196
197 allocate( growth_stage_gdd( 5, inumberoflanduses ), stat=istat )
198 call assert( istat==0, "Failed to allocate memory for GROWTH_STAGE_GDD array", &
199 __file__, __line__ )
200
201 allocate( growth_stage_date( 6, inumberoflanduses ), stat=istat )
202 call assert( istat==0, "Failed to allocate memory for GROWTH_STAGE_DATE array", &
203 __file__, __line__ )
204
205 allocate( kcb_l( 16, inumberoflanduses ), stat=istat )
206 call assert( istat==0, "Failed to allocate memory for KCB_l array", &
207 __file__, __line__ )
208
209 allocate( kcb_method( inumberoflanduses ), stat=istat )
210 call assert( istat==0, "Failed to allocate memory for KCB_METHOD vector", &
211 __file__, __line__ )
212
213 kcb_method = -9999
214 kcb_l = -9999.
215 growth_stage_gdd = -9999.
217
218 if ( ubound(l_ini_l,1) == inumberoflanduses ) &
220
221 if ( ubound(l_dev_l,1) == inumberoflanduses ) &
223
224 if ( ubound(l_mid_l,1) == inumberoflanduses ) &
226
227 if ( ubound(l_late_l,1) == inumberoflanduses ) &
229
230 if ( ubound(l_fallow_l,1) == inumberoflanduses ) &
232
233 call logs%write(" ## Crop Kcb Curve Summary ##", ilinesafter=1)
234 call logs%write(" _only meaningful for landuses where the Kcb curve is defined " &
235 //"in terms of days _", ilinesafter=1)
236 call logs%write("Landuse Code | Planting Date | End of 'ini' | End of 'dev' " &
237 //"| End of 'mid' | End of 'late' | End of 'fallow' ")
238 call logs%write("-------------|---------------|--------------|--------------" &
239 //"|--------------|---------------|-----------------")
240
241 if ( slplantingdate%count == inumberoflanduses .and. slplantingdate%count > 0 ) then
242
243 do iindex=1, slplantingdate%count
244
245 plantingdate_str = slplantingdate%get( iindex )
246
247 ! if there is no planting date entry, assume user is specifying
248 ! a planting GDD for this landuse code
249 if ( len_trim(plantingdate_str) == 0 ) cycle
250
251 if ( plantingdate_str .contains. "/" ) then
252
253 ! append current year to the end of the user-entered planting date in mm/dd
254 smmddyyyy = trim(plantingdate_str)//"/"//ascharacter( sim_dt%start%iYear )
255 call growth_stage_date( planting_date, iindex)%parsedate( smmddyyyy, __file__, __line__ )
256
258
259 else
260 ! assume the value is a day-of-year value
261 dtplantingdate = sim_dt%start + asfloat( plantingdate_str)
262 call dtplantingdate%calcJulianDay()
263 growth_stage_date( planting_date, iindex) = dtplantingdate
264 endif
265
266 ! march forward through time calculating the various dates on the Kcb curve
267 ! GROWTH_STAGE_DATE( ENDDATE_INI, iIndex ) = GROWTH_STAGE_DATE( PLANTING_DATE, iIndex ) + L_ini_l( iIndex )
268 ! GROWTH_STAGE_DATE( ENDDATE_DEV, iIndex ) = GROWTH_STAGE_DATE( ENDDATE_INI, iIndex ) + L_dev_l( iIndex )
269 ! GROWTH_STAGE_DATE( ENDDATE_MID, iIndex ) = GROWTH_STAGE_DATE( ENDDATE_DEV, iIndex ) + L_mid_l( iIndex )
270 ! GROWTH_STAGE_DATE( ENDDATE_LATE, iIndex ) = GROWTH_STAGE_DATE( ENDDATE_MID, iIndex ) + L_late_l( iIndex )
271 ! GROWTH_STAGE_DATE( ENDDATE_FALLOW, iIndex ) = GROWTH_STAGE_DATE( ENDDATE_LATE, iIndex ) + L_fallow_l( iIndex )
272
273 ! if any of the L_* length values is missing, a value of zero will be used, resulting in a wierd looking Kcb curve
284
285 call logs%write( "| "//ascharacter( landuse_code( iindex ))//" | " &
286 //trim( growth_stage_date( planting_date, iindex )%prettydate() ) &
287 //" (doy:"//ascharacter( growth_stage_date( planting_date, iindex )%getDayOfYear() )//") | " &
288 //trim( growth_stage_date( enddate_ini, iindex )%prettydate() )//" | " &
289 //" (doy:"//ascharacter( growth_stage_date( enddate_ini, iindex )%getDayOfYear() )//") | " &
290 //trim( growth_stage_date( enddate_dev, iindex )%prettydate() )//" | " &
291 //" (doy:"//ascharacter( growth_stage_date( enddate_dev, iindex )%getDayOfYear() )//") | " &
292 //trim( growth_stage_date( enddate_mid, iindex )%prettydate() )//" | " &
293 //" (doy:"//ascharacter( growth_stage_date( enddate_mid, iindex )%getDayOfYear() )//") | " &
294 //trim( growth_stage_date( enddate_late, iindex )%prettydate() )//" | " &
295 //" (doy:"//ascharacter( growth_stage_date( enddate_late, iindex )%getDayOfYear() )//") | " &
296 //trim( growth_stage_date( enddate_fallow, iindex )%prettydate() ) &
297 //" (doy:"//ascharacter( growth_stage_date( enddate_fallow, iindex )%getDayOfYear() )//") | ")
298 enddo
299
300 endif
301
302 if (ubound(gdd_plant_l,1) == inumberoflanduses) growth_stage_gdd( gdd_plant, : ) = gdd_plant_l
303 if (ubound(gdd_ini_l,1) == inumberoflanduses) growth_stage_gdd( gdd_ini, : ) = gdd_ini_l
304 if (ubound(gdd_dev_l,1) == inumberoflanduses) growth_stage_gdd( gdd_dev, : ) = gdd_dev_l
305 if (ubound(gdd_mid_l,1) == inumberoflanduses) growth_stage_gdd( gdd_mid, : ) = gdd_mid_l
306 if (ubound(gdd_late_l,1) == inumberoflanduses) growth_stage_gdd( gdd_late, : ) = gdd_late_l
307
308 if (ubound(kcb_ini_l,1) == inumberoflanduses) kcb_l( kcb_ini, :) = kcb_ini_l
309 if (ubound(kcb_mid_l,1) == inumberoflanduses) kcb_l( kcb_mid, :) = kcb_mid_l
310 if (ubound(kcb_end_l,1) == inumberoflanduses) kcb_l( kcb_end, :) = kcb_end_l
311 if (ubound(kcb_min_l,1) == inumberoflanduses) kcb_l( kcb_min, :) = kcb_min_l
312
313 if (ubound(kcb_jan,1) == inumberoflanduses) kcb_l( jan, :) = kcb_jan
314 if (ubound(kcb_feb,1) == inumberoflanduses) kcb_l( feb, :) = kcb_feb
315 if (ubound(kcb_mar,1) == inumberoflanduses) kcb_l( mar, :) = kcb_mar
316 if (ubound(kcb_apr,1) == inumberoflanduses) kcb_l( apr, :) = kcb_apr
317 if (ubound(kcb_may,1) == inumberoflanduses) kcb_l( may, :) = kcb_may
318 if (ubound(kcb_jun,1) == inumberoflanduses) kcb_l( jun, :) = kcb_jun
319 if (ubound(kcb_jul,1) == inumberoflanduses) kcb_l( jul, :) = kcb_jul
320 if (ubound(kcb_aug,1) == inumberoflanduses) kcb_l( aug, :) = kcb_aug
321 if (ubound(kcb_sep,1) == inumberoflanduses) kcb_l( sep, :) = kcb_sep
322 if (ubound(kcb_oct,1) == inumberoflanduses) kcb_l( oct, :) = kcb_oct
323 if (ubound(kcb_nov,1) == inumberoflanduses) kcb_l( nov, :) = kcb_nov
324 if (ubound(kcb_dec,1) == inumberoflanduses) kcb_l( dec, :) = kcb_dec
325
326 ! go through the table values and try to figure out how Kcb curves should be constructed:
327 ! Monthly Kcb, GDD-based, or DOY-based
328 do iindex = lbound( kcb_method, 1), ubound( kcb_method, 1)
329
330 if ( all( kcb_l( jan:dec, iindex ) > 0.0_c_float ) ) then
332 kcb_l( kcb_min, iindex ) = minval( kcb_l(jan:dec, iindex) )
333 kcb_l( kcb_mid, iindex) = minval( kcb_l(jan:dec, iindex) )
334
335 elseif ( all( growth_stage_gdd( :, iindex ) >= 0.0_c_float ) &
336 .and. all( kcb_l( kcb_ini:kcb_min, iindex ) > 0.0_c_float ) ) then
337 kcb_method( iindex ) = kcb_method_gdd
338
339 elseif ( all( growth_stage_length_in_days( planting_date:, iindex ) >= 0.0_c_float ) &
340 .and. all( kcb_l( kcb_ini:kcb_min, iindex ) > 0.0_c_float ) ) then
341 kcb_method( iindex ) = kcb_method_fao56
342 endif
343
344 if ( kcb_method( iindex ) < 0 ) then
345 call warn("There are missing day-of-year (L_ini, L_dev, L_mid, L_late, L_fallow), " &
346 //"growing degree-day ~(GDD_plant, GDD_ini, GDD_dev, GDD_mid, GDD_late)," &
347 //" or monthly crop ~coefficients (Kcb_jan...Kcb_dec) for" &
348 //" landuse "//ascharacter( landuse_code( iindex ) ), lfatal=true )
349 endif
350
351 enddo
352
353! do iIndex = lbound( fSoilStorage, 1 ), ubound( fSoilStorage,1 )
354
355! fKcb_initial = update_crop_coefficient_date_as_threshold( iLanduseIndex( iIndex ) )
356
357 ! call calc_effective_root_depth( fRz_i=fRz_initial, iLanduseIndex=iLanduseIndex( iIndex ), &
358 ! fZr_max=fMax_Rooting_Depths( iLanduseIndex( iIndex ), &
359 ! iSoilGroup( iIndex ) ), &
360 ! Kcb=fKcb_initial )
361 !
362 ! fSoilStorage( iIndex ) = INITIAL_PERCENT_SOIL_MOISTURE( iIndex ) / 100.0_c_float &
363 ! * fRz_initial * fAvailable_Water_Content( iIndex )
364
365! enddo
366
367 !> @TODO Add more logic here to perform checks on the validity of this data.
368
369 !> @TODO Need to handle missing values. WHat do we do if an entire column of values
370 !! is missing?
371
373
374!------------------------------------------------------------------------------
375
376 !> Update the current basal crop coefficient (Kcb) for
377 !! a SINGLE irrigation table entry
378 !!
379 !! @param[inout] pIRRIGATION pointer to a single line of information in the irrigation file.
380 !! @param[in] iThreshold either the current day of year or the number of growing degree days.
381 !! @retval rKcb Basal crop coefficient given the irrigation table entries and the
382 !! current threshold values.
383
384 impure elemental function update_crop_coefficient_date_as_threshold( iLanduseIndex ) &
385 result(kcb)
386
387 integer (c_int), intent(in) :: ilanduseindex
388 real (c_float) :: kcb
389
390 ! [ LOCALS ]
391 real (c_double) :: ffrac
392
393 if ( kcb_method( ilanduseindex ) == kcb_method_monthly_values ) then
394
395 kcb = kcb_l( sim_dt%curr%iMonth, ilanduseindex )
396
397 else
398
399 ! define shorthand variable names for remainder of function
400 associate( date_ini => growth_stage_date( enddate_ini, ilanduseindex ), &
401 date_dev => growth_stage_date( enddate_dev, ilanduseindex ), &
402 date_mid => growth_stage_date( enddate_mid, ilanduseindex ), &
403 date_late => growth_stage_date( enddate_late, ilanduseindex ), &
404 date_fallow => growth_stage_date( enddate_fallow, ilanduseindex ), &
405 kcb_ini => kcb_l(kcb_ini, ilanduseindex), &
406 kcb_mid => kcb_l(kcb_mid, ilanduseindex), &
407 kcb_min => kcb_l(kcb_min, ilanduseindex), &
408 plantingdate => growth_stage_date( planting_date, ilanduseindex), &
409 kcb_end => kcb_l(kcb_end, ilanduseindex), &
410 current_date => sim_dt%curr )
411
412 ! now calculate Kcb for the given landuse
413
414 if( current_date > date_late ) then
415
416 kcb = kcb_min
417
418 elseif ( current_date > date_mid ) then
419
420 ffrac = ( current_date - date_mid ) / ( date_late - date_mid )
421
422 kcb = kcb_mid * (1.0_c_double - ffrac) + kcb_end * ffrac
423
424 elseif ( current_date > date_dev ) then
425
426 kcb = kcb_mid
427
428 elseif ( current_date > date_ini ) then
429
430 ffrac = ( current_date - date_ini ) / ( date_dev - date_ini )
431
432 kcb = kcb_ini * (1.0_c_double - ffrac) + kcb_mid * ffrac
433
434 elseif ( current_date >= plantingdate ) then
435
436 kcb = kcb_ini
437
438 else
439
440 kcb = kcb_min
441
442 endif
443
444 end associate
445
446 end if
447
449
450!------------------------------------------------------------------------------
451
452pure elemental function crop_coefficients_fao56_calculate_kcb_max(wind_speed_meters_per_sec, &
453 relative_humidity_min_pct, &
454 Kcb, &
455 plant_height_meters) result(kcb_max)
456
457 real (c_float), intent(in) :: wind_speed_meters_per_sec
458 real (c_float), intent(in) :: relative_humidity_min_pct
459 real (c_float), intent(in) :: kcb
460 real (c_float), intent(in) :: plant_height_meters
461
462 real (c_float) :: kcb_max
463 real (c_double) :: u2
464 real (c_double) :: rhmin
465 real (c_double) :: plant_height
466
467 ! Limits are as suggested on page 123 of FAO-56 with respect to
468 ! modifying mid-season KCB_mid values
469 rhmin = clip( relative_humidity_min_pct, minval=20., maxval=80. )
470 u2 = clip(wind_speed_meters_per_sec, minval=1., maxval=6.)
471 plant_height = clip(plant_height_meters, minval=1., maxval=10.)
472
473 ! equation 72, FAO-56, p 199
474 kcb_max = max( 1.2_c_double + ( (0.04_c_double * (u2 - 2._c_double) &
475 - 0.004_c_double * (rhmin - 45._c_double) ) ) &
476 * (plant_height_meters/3._c_double)**0.3_c_double, &
477 kcb + 0.05_c_double )
478
480
481!------------------------------------------------------------------------------
482
483 !> Update the current basal crop coefficient (Kcb), with GDD as the threhold
484 !!
485 !! @param[in] fGDD current growing degree day value associated with the cell.
486 !! @retval fKcb Basal crop coefficient given the irrigation table entries and the
487 !! current threshold values.
488
489 impure elemental function update_crop_coefficient_gdd_as_threshold( iLanduseIndex, fGDD ) &
490 result(fkcb)
491
492 integer (c_int), intent(in) :: ilanduseindex
493 real (c_float), intent(in) :: fgdd
494 real (c_float) :: fkcb
495
496 ! [ LOCALS ]
497 real (c_double) :: ffrac
498
499 ! define shorthand variable names for remainder of function
500 associate( gdd_ini_l => growth_stage_gdd( gdd_ini, ilanduseindex ), &
501 gdd_dev_l => growth_stage_gdd( gdd_dev, ilanduseindex ), &
502 gdd_mid_l => growth_stage_gdd( gdd_mid, ilanduseindex ), &
503 gdd_late_l => growth_stage_gdd( gdd_late, ilanduseindex ), &
504 kcb_ini => kcb_l(kcb_ini, ilanduseindex), &
505 kcb_mid => kcb_l(kcb_mid, ilanduseindex), &
506 kcb_min => kcb_l(kcb_min, ilanduseindex), &
507 plantingdoy => growth_stage_date( planting_date, &
508 ilanduseindex)%getDayOfYear(), &
509 gdd_plant_l => growth_stage_gdd( gdd_plant, ilanduseindex), &
510 kcb_end => kcb_l(kcb_end, ilanduseindex), &
511 current_doy => sim_dt%curr%getDayOfYear() )
512
513 ! now calculate Kcb for the given landuse
514 if( fgdd > gdd_late_l ) then
515
516 fkcb = kcb_min
517
518 elseif ( fgdd > gdd_mid_l ) then
519
520 ffrac = ( fgdd - gdd_mid_l ) / ( gdd_late_l - gdd_mid_l )
521
522 fkcb = kcb_mid * (1.0_c_double - ffrac) + kcb_end * ffrac
523
524 elseif ( fgdd > gdd_dev_l ) then
525
526 fkcb = kcb_mid
527
528 elseif ( fgdd > gdd_ini_l ) then
529
530 ffrac = ( fgdd - gdd_ini_l ) / ( gdd_dev_l - gdd_ini_l )
531
532 fkcb = kcb_ini * (1_c_double - ffrac) + kcb_mid * ffrac
533
534 elseif ( (plantingdoy > 0) .and. (current_doy >= plantingdoy) ) then
535
536 ! if there is a valid value for the Planting Date, use it
537 fkcb = kcb_ini
538
539 elseif ( (gdd_plant_l >= 0.) .and. (fgdd >= gdd_plant_l) ) then
540
541 fkcb = kcb_ini
542
543 else
544
545 fkcb = kcb_min
546
547 endif
548
549 end associate
550
552
553!------------------------------------------------------------------------------
554
556
557 ! [ LOCALS ]
558 integer (c_int) :: iindex
559 real (c_double) :: dtempdate
560 type (datetime_t) :: dttempdate
561 real (c_float) :: growing_cycle_length
562
563 do iindex=lbound(growth_stage_date,2), ubound(growth_stage_date,2)
564
565 ! print *, SIM_DT%curr%prettydate(), " | ", GROWTH_STAGE_DATE( ENDDATE_LATE, iIndex )%prettydate()
566
567 ! if we have not yet reached the enddate associated with the fallow period, skip
568 if ( sim_dt%curr < growth_stage_date( enddate_fallow, iindex ) ) cycle
569
570 if ( kcb_method( iindex ) /= kcb_method_fao56 ) cycle
571
572 ! current date is beyond the enddate associated with fallow period;
573 ! update Kcb curve and dates
574
575! print *, "a) ",SIM_DT%curr%prettydate(), " | ", GROWTH_STAGE_DATE( ENDDATE_LATE, iIndex )%prettydate()
576
577 ! it's possible that the planting date might be later in the current calendar year
578 call growth_stage_date( planting_date, iindex )%setYear( sim_dt%curr%iYear )
579
580! print *, "b) ",SIM_DT%curr%prettydate(), " | ", GROWTH_STAGE_DATE( PLANTING_DATE, iIndex )%prettydate()
581
582 ! however, if we are already past that point in the year, planting date must be
583 ! next calendar year
584 if ( sim_dt%iDOY > growth_stage_date( planting_date, iindex )%getDayOfYear() ) &
585 call growth_stage_date( planting_date, iindex )%addYear()
586
587! print *, "c) ",SIM_DT%curr%prettydate(), " | ", GROWTH_STAGE_DATE( PLANTING_DATE, iIndex )%prettydate()
588
589 ! now calculate dates associated with the rest of the Kcb curve
600
601 call logs%write("## Updating Kcb Date Values ##", ilinesafter=1, lecho=false )
602 call logs%write("Landuse Code | Planting Date | End of 'ini' | End of 'dev' " &
603 //"| End of 'mid' | End of 'late' | End of 'fallow' ", lecho=false )
604 call logs%write("-------------|---------------|--------------|--------------" &
605 //"|--------------|---------------|-----------------", lecho=false )
606
607 call logs%write( ascharacter( landuse_code( iindex ))//" | " &
608 //trim( growth_stage_date( planting_date, iindex )%prettydate() )//" | " &
609 //trim( growth_stage_date( enddate_ini, iindex )%prettydate() )//" | " &
610 //trim( growth_stage_date( enddate_dev, iindex )%prettydate() )//" | " &
611 //trim( growth_stage_date( enddate_mid, iindex )%prettydate() )//" | " &
612 //trim( growth_stage_date( enddate_late, iindex )%prettydate() )//" | " &
613 //trim( growth_stage_date( enddate_fallow, iindex )%prettydate() ), &
614 lecho=false, iloglevel=log_all )
615
616 enddo
617
618
620
621!--------------------------------------------------------------------------------------------------
622
623 impure elemental subroutine crop_coefficients_fao56_calculate( Kcb, landuse_index, GDD )
624
625 real (c_float), intent(inout) :: kcb
626 integer (c_int), intent(in) :: landuse_index
627 real (c_float), intent(in), optional :: gdd
628
629
630 if ( kcb_method( landuse_index ) == kcb_method_fao56 &
631 .or. kcb_method( landuse_index ) == kcb_method_monthly_values ) then
632
633 kcb = update_crop_coefficient_date_as_threshold( landuse_index )
634
635 else
636
637 kcb = update_crop_coefficient_gdd_as_threshold( landuse_index, gdd )
638
639 endif
640
642
643!--------------------------------------------------------------------------------------------------
644
646 landuse_index, &
647 Kcb, &
648 it_is_growing_season)
649
650 real (c_float), intent(in) :: kcb
651 integer (c_int), intent(in) :: landuse_index
652 logical (c_bool), intent(out) :: it_is_growing_season
653
654 if ( kcb > kcb_l( kcb_min, landuse_index) ) then
655 it_is_growing_season = true
656 else
657 it_is_growing_season = false
658 endif
659
661
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
real(c_double), parameter, public m_per_foot
integer(c_int), parameter, public itinyval
Update crop coefficients for crop types in simulation.
subroutine, public crop_coefficients_fao56_initialize()
impure elemental real(c_float) function, public update_crop_coefficient_gdd_as_threshold(ilanduseindex, fgdd)
Update the current basal crop coefficient (Kcb), with GDD as the threhold.
pure elemental real(c_float) function, public crop_coefficients_fao56_calculate_kcb_max(wind_speed_meters_per_sec, relative_humidity_min_pct, kcb, plant_height_meters)
real(c_float), dimension(:,:), allocatable, public growth_stage_length_in_days
real(c_float), dimension(:), allocatable growth_stage_shift_days
impure elemental subroutine, public crop_coefficients_fao56_update_growing_season(landuse_index, kcb, it_is_growing_season)
real(c_float), dimension(:,:), allocatable, public kcb_l
subroutine, public crop_coefficients_fao56_update_growth_stage_dates()
impure elemental real(c_float) function, public update_crop_coefficient_date_as_threshold(ilanduseindex)
Update the current basal crop coefficient (Kcb) for a SINGLE irrigation table entry.
real(c_float), dimension(:,:), allocatable growth_stage_gdd
integer(c_int), dimension(:), allocatable, public kcb_method
impure elemental subroutine, public crop_coefficients_fao56_calculate(kcb, landuse_index, gdd)
integer(c_int), dimension(:), allocatable landuse_code
type(datetime_t), dimension(:,:), allocatable, public growth_stage_date
Defines the DATA_CATALOG_T data type, which contains type-bound procedures to add,...
type(data_catalog_t), public dat
DAT is a global to hold data catalog entries.
This module contains the DATETIME_T class and associated time and date-related routines,...
Definition datetime.F90:9
subroutine, public warn(smessage, smodule, iline, shints, lfatal, iloglevel, lecho)
subroutine, public die(smessage, smodule, iline, shints, scalledby, icalledbyline)
type(logfile_t), public logs
Definition logfiles.F90:62
type(parameters_t), public params
type(date_range_t), public sim_dt