Soil Water Balance (SWB2)
Loading...
Searching...
No Matches
precipitation__method_of_fragments.F90
Go to the documentation of this file.
1!> @file
2!! Contains the module @ref precipitation__method_of_fragments.
3
4!>
5!! Module @ref precipitation__method_of_fragments
6!! provides support for creating synthetic daily precipitation
7!! given grids of monthly sum precipitation and a "fragments" file.
8!! The fragments file is generated from observations at discrete locations,
9!! the values of which range from 0 to 1, and the sum of which is 1.
10!! The fragment value is simply the daily observed precipitation value divided
11!! by the monthly sum of all observed precipitation values for that station.
12!!
13!! In addition, this routine accepts another set of rainfall adjustment grids,
14!! needed in order to ensure the resulting precipitation totals fall in line with
15!! other published values, in the development case, the Rainfall Atlas of Hawaii.
16
18
19 use iso_c_binding
21 use data_catalog
23 use dictionary
24 use exceptions
27 use logfiles, only : logs, log_all, log_debug
28 use parameters
29 use fstring
30 use fstring_list
32 use grid
33 implicit none
34
35 private
36
38 public :: read_daily_fragments
40
41 !> Module variable that holds the rainfall gage (zone) number
42 integer (c_int), allocatable, public :: rain_gage_id(:)
43
44 !> Module variable that holds the current day's rainfall fragment value
45 real (c_float), allocatable, public :: fragment_value(:)
46
47 !> Module variable indicating which "simulation number" is active
48 !! Only has meaning if the rainfall fragments are being applied via a predetermined
49 !! sequence file
50 integer (c_int), public :: simulation_number = 1
51
52 !> Module variable that holds the rainfall adjustment factor
53 real (c_float), allocatable, public :: rainfall_adjust_factor(:)
54
55 !> Module variable that holds a sequence of random numbers associated with the selection
56 !! of the fragment set to use
57 real (c_double), allocatable :: random_values(:,:)
58
59 !> Module level variable used to create subsets of the FRAGMENT_SEQUENCES file
60 logical (c_bool), allocatable :: sequence_selection(:)
61
62 !> Module variable detemining whether fragment sequences are chosen at random or
63 !! selected from an external file
64 logical (c_bool) :: random_fragment_sequences = .true._c_bool
65
66 !> Data structure that holds a single line of data from the input rainfall fragments file.
67 type, public :: fragments_t
68 integer (c_int) :: imonth
69 integer (c_int) :: iraingagezone
70 integer (c_int) :: ifragmentset
71 real (c_float) :: ffragmentvalue(31)
72 end type fragments_t
73
74 !> Pointer to a rainfall fragments data structure.
75 type, public :: ptr_fragments_t
76 type (fragments_t), pointer :: pfragment => null()
77 end type ptr_fragments_t
78
79 !> Data structure to hold the current active rainfall fragments for
80 !! a particular rain gage zone.
81 type, public :: fragments_set_t
82 integer (c_int) :: iraingagezone
83 integer (c_int) :: inumberoffragments(12)
84 integer (c_int) :: istartrecord(12)
85 end type fragments_set_t
86
87 !> Array of all fragments read in from the rainfall fragments file.
88 type (fragments_t), allocatable, target, public :: fragments(:)
89
90 !> Subset of rainfall fragments file pointing to the currently active fragments.
91 type (ptr_fragments_t), allocatable :: current_fragments(:,:)
92
93 !> Array of fragments sets; fragments sets include indices to the start record
94 !! associated with the fragment for each month; FRAGMENTS_SETS will have a
95 !! number of elements equal to the number of rainfall gages in the model domain
96 type (fragments_set_t), allocatable, public :: fragments_sets(:)
97
98 !> Data structure to hold static (pre-calculated) fragment selection numbers
99 type, public :: fragments_sequence_t
100 integer (c_int) :: sim_number
101 integer (c_int) :: sim_month
102 integer (c_int) :: sim_rainfall_zone
103 integer (c_int) :: sim_year
104 real (c_float) :: sim_random_number
105 integer (c_int) :: sim_selected_set
106 end type fragments_sequence_t
107
108 !> Pointer to all or some of the FRAGMENTS_SEQUENCE array
110
111 !> Array of fragment sequence sets
112 type (fragments_sequence_t), allocatable, public :: fragments_sequence(:)
113
115
116 integer (c_int) :: lu_fragments_echo
117
118contains
119
120 !> Initialize method of fragments.
121 !!
122 !! This routine accesses the "RAINFALL_ZONE" gridded data object and
123 !! calls the routine to read in the rainfall fragments file. Values of RAINFALL_ZONE are stored
124 !! in a module variable @ref RAIN_GAGE_ID for future reference.
125 !!
126 !! @params[in] lActive 2-D boolean array defining active and inactive cells
127
129
130 logical (c_bool), intent(in) :: lactive(:,:)
131
132 ! [ LOCALS ]
133 integer (c_int) :: istat
134 type (data_catalog_entry_t), pointer :: prainfall_zone
135 type (fstring_list_t) :: slstring
136 integer (c_int) :: imaxrainzones
137 integer (c_int), allocatable :: isimulationnumbers(:)
138 character (len=256) :: error_str
139
140 ! look up the simulation number associated with the desired fragment sequence set
141 call cf_dict%get_values( skey="FRAGMENTS_SEQUENCE_SIMULATION_NUMBER", ivalues=isimulationnumbers )
142 if ( isimulationnumbers(1) > 0 ) simulation_number = isimulationnumbers(1)
143
144 ! locate the data structure associated with the gridded rainfall zone entries
145 prainfall_zone => dat%find("RAINFALL_ZONE")
146 if ( .not. associated(prainfall_zone) ) &
147 call die("A RAINFALL_ZONE grid must be supplied in order to make use of this option.", &
148 __file__, __line__)
149
150 allocate( rain_gage_id( count(lactive) ), stat=istat )
151 call assert( istat == 0, "Problem allocating memory", __file__, __line__ )
152
153 call prainfall_zone%getvalues()
154
155 ! map the 2D array of RAINFALL_ZONE values to the vector of active cells
156 rain_gage_id = pack( prainfall_zone%pGrdBase%iData, lactive )
157
158 allocate( rainfall_adjust_factor( count(lactive) ), stat=istat )
159 call assert( istat == 0, "Problem allocating memory", __file__, __line__ )
160
161 allocate( fragment_value( count(lactive) ), stat=istat )
162 call assert( istat == 0, "Problem allocating memory", __file__, __line__ )
163
164
165 ! look up the name of the fragments file in the control file dictionary
166 call cf_dict%get_values( skey="FRAGMENTS_DAILY_FILE", slstring=slstring )
167
168 ! use the first entry in the string list slString as the filename to open for
169 ! use with the daily fragments routine
170
171
172 call read_daily_fragments( slstring%get(1) )
173 call slstring%clear()
174
175 ! look up the name of the fragments SEQUENCE file in the control file dictionary
176 call cf_dict%get_values( skey="FRAGMENTS_SEQUENCE_FILE", slstring=slstring )
177
178 if ( .not. ( slstring%get(1) .strequal. "<NA>" ) ) then
179 call read_fragments_sequence( slstring%get(1) )
181 allocate ( sequence_selection( count(fragments_sequence%sim_month > 0) ), stat=istat )
182 call assert( istat == 0, "Problem allocating memory", __file__, __line__ )
183 endif
184
185 !> Now the fragments file is in memory. Create an ancillary data structure
186 !> to keep track of which records correspond to various rain zones
187
188 imaxrainzones = maxval(fragments%iRainGageZone)
189
190 allocate ( fragments_sets( imaxrainzones ), stat=istat )
191 call assert( istat == 0, "Problem allocating memory", __file__, __line__ )
192
193 if ( .not. allocated(current_fragments) ) then
194 allocate (current_fragments( imaxrainzones, 1 ), stat=istat, errmsg=error_str )
195 call assert( istat == 0, "Problem allocating memory, stat="//ascharacter(istat) &
196 //"; msg: "//trim(error_str), __file__, __line__ )
197 endif
198
199 if ( .not. allocated( random_values) ) then
200 allocate (random_values( imaxrainzones, 1 ), stat=istat )
201 call assert( istat == 0, "Problem allocating memory, stat="//ascharacter(istat) &
202 //"; msg: "//trim(error_str), __file__, __line__ )
203 endif
204
206
208
209 open( newunit=lu_fragments_echo, file="Fragments_as_implemented_by_SWB.csv")
210 write( lu_fragments_echo, fmt="(a, 30('fragment, '),'fragment')") &
211 "Simulation_Number, Month, Rain_Zone, Year, Random_Number, Fragment_Set,"
212
214
215 !--------------------------------------------------------------------------------------------------
216
217 subroutine read_daily_fragments( sFilename )
218
219 character (len=*), intent(in) :: sfilename
220
221 ! [ LOCALS ]
222 character (len=512) :: srecord, ssubstring
223 integer (c_int) :: istat
224 integer (c_int) :: icount
225 integer (c_int) :: iindex
226 integer (c_int) :: last_zone
227 integer (c_int) :: last_fragment
228 integer (c_int) :: last_month
229 integer (c_int) :: inumlines
230 real (c_float) :: ftempvalue
231 type (ascii_file_t) :: fragments_file
232
233
234 call fragments_file%open( sfilename = sfilename, &
235 scommentchars = "#%!", &
236 sdelimiters = "WHITESPACE", &
237 lhasheader = .false._c_bool )
238
239 inumlines = fragments_file%numLines()
240
241 allocate( fragments( inumlines ), stat=istat )
242 call assert( istat == 0, "Problem allocating memory for fragments table", __file__, __line__ )
243
244 icount = 0
245 last_zone = 1
246 last_fragment = 0
247 last_month = 1
248
249 do
250
251 ! read in next line of file
252 srecord = fragments_file%readLine()
253
254 if ( fragments_file%isEOF() ) exit
255
256 icount = icount + 1
257
258 ! read in month number
259 call chomp(srecord, ssubstring, fragments_file%sDelimiters )
260
261 if ( len_trim(ssubstring) == 0 ) &
262 call die( "Missing month number in the daily fragments file", &
263 __file__, __line__, "Problem occured on line number " &
264 //ascharacter(fragments_file%currentLineNum() ) &
265 //" of file "//dquote(sfilename) )
266
267 fragments(icount)%iMonth = asint(ssubstring)
268
269 ! reset the counter tracking the previous fragment ID
270 if ( fragments(icount)%iMonth == (last_month + 1) ) then
271 last_fragment = 0
272 last_zone = 0
273 endif
274
275 if ( fragments(icount)%iMonth < last_month ) &
276 call die( "Out-of-order month value in the daily fragments file", &
277 __file__, __line__, "Problem occured on line number " &
278 //ascharacter(fragments_file%currentLineNum() ) &
279 //" of file "//dquote(sfilename) )
280
281 last_month = fragments(icount)%iMonth
282
283 ! read in rain gage zone
284 call chomp(srecord, ssubstring, fragments_file%sDelimiters )
285
286 if ( len_trim(ssubstring) == 0 ) &
287 call die( "Missing rain gage zone number in the daily fragments file", &
288 __file__, __line__, "Problem occured on line number " &
289 //ascharacter(fragments_file%currentLineNum() ) &
290 //" of file "//dquote(sfilename) )
291
292 fragments(icount)%iRainGageZone = asint(ssubstring)
293
294 if ( fragments(icount)%iRainGageZone < last_zone ) &
295 call die( "Rain gage zone number out of order in the daily fragments file", &
296 __file__, __line__, "Problem occured on line number " &
297 //ascharacter(fragments_file%currentLineNum() ) &
298 //" of file "//dquote(sfilename) )
299
300 ! reset the counter tracking the previous fragment ID
301 if ( fragments(icount)%iRainGageZone == (last_zone + 1) ) &
302 last_fragment = 0
303
304 last_zone = fragments(icount)%iRainGageZone
305
306 ! read in fragment set number for this zone
307 call chomp(srecord, ssubstring, fragments_file%sDelimiters )
308
309 if ( len_trim(ssubstring) == 0 ) &
310 call die( "Missing fragment set number in the daily fragments file", &
311 __file__, __line__, "Problem occured on line number " &
312 //ascharacter(fragments_file%currentLineNum() ) &
313 //" of file "//dquote(sfilename) )
314
315 fragments(icount)%iFragmentSet = asint(ssubstring)
316
317 if ( fragments(icount)%iFragmentSet /= (last_fragment + 1) ) &
318 call die( "Missing or out-of-order fragment value in the daily fragments file", &
319 __file__, __line__, "Problem occured on line number " &
320 //ascharacter(fragments_file%currentLineNum() ) &
321 //" of file "//dquote(sfilename) )
322
323 last_fragment = fragments(icount)%iFragmentSet
324
325 ! read in fragments for each day of the current month
326 do iindex = 1, 31
327
328 ! read in fragment for given day of month
329 call chomp(srecord, ssubstring, fragments_file%sDelimiters )
330
331 if ( len_trim(ssubstring) == 0 ) &
332 call die( "Missing fragment value in the daily fragments file", &
333 __file__, __line__, "Problem occured on line number " &
334 //ascharacter(fragments_file%currentLineNum() ) &
335 //" of file "//dquote(sfilename) )
336
337 ftempvalue = asfloat( ssubstring )
338
339 ! This substitution is needed to prevent "-9999" or "9999" values embedded in the fragments file
340 ! from creeping into calculations. In the event of a "9999", the appropriate substitution is
341 ! zero, since the previous fragments for the month to this point should already sum to 1.0
342 if ( ( ftempvalue < 0.0_c_float ) .or. ( ftempvalue > 1.0_c_float ) ) then
343 fragments(icount)%fFragmentValue(iindex) = 0.0_c_float
344 else
345 fragments(icount)%fFragmentValue(iindex) = ftempvalue
346 endif
347
348 enddo
349
350 if ( fragments(icount)%iMonth == 2 ) call normalize_february_fragment_sequence( icount )
351
352 enddo
353
354 call logs%write("Maximum rain gage zone number: "//ascharacter(maxval(fragments%iRainGageZone)), &
355 itab=31, ilinesafter=1, iloglevel=log_all)
356
357 end subroutine read_daily_fragments
358
359!--------------------------------------------------------------------------------------------------
360
361 !> after fragments file has been read in, iterate over a set of rainfall fragments
362 ! and keep track of the index values that correspond with changes
363 ! in month and rain gage numbers
365
366 integer (c_int) :: iCount
367 integer (c_int) :: iIndex
368 integer (c_int) :: iRainGageZone
369 integer (c_int) :: iPreviousRainGageZone
370 integer (c_int) :: iFragmentChunk
371 integer (c_int) :: iMonth
372 integer (c_int) :: iPreviousMonth
373 character (len=10) :: sBuf0
374 character (len=10) :: sBuf1
375 character (len=12) :: sBuf2
376 character (len=10) :: sBuf3
377 character (len=52) :: sBuf4
378
379 ! this counter is used to accumulate the number of fragments associated with the
380 ! current raingage zone/month combination
381 icount = 1
382
383 iraingagezone = fragments( lbound( fragments, 1) )%iRainGageZone
384 ipreviousraingagezone = iraingagezone
385 ipreviousmonth = fragments( lbound( fragments, 1) )%iMonth
386
387 ! at this point, iRainGageZone should be 1, and iPreviousMonth should be 1,
388 ! assuming that the fragments were sorted properly upon input
389
390 ! populate the first record of FRAGMENT_SETS
391 fragments_sets( iraingagezone )%iRainGageZone = iraingagezone
392 fragments_sets( iraingagezone )%iStartRecord(ipreviousmonth) = lbound( fragments, 1)
393
394 ! now iterate through *all* fragments, keeping track of the starting record for each new rainfall gage
395 ! zone number
396 do iindex = lbound( fragments, 1) + 1, ubound( fragments, 1 )
397
398 iraingagezone = fragments(iindex)%iRainGageZone
399 imonth = fragments(iindex)%iMonth
400
401 if ( iraingagezone /= ipreviousraingagezone ) then
402 ! the previous record was the last one associated with the previous
403 ! rainfall gage zone; do not count the current record as part of the
404 ! collection of records associated with previous zone
405
406 fragments_sets( ipreviousraingagezone )%iNumberOfFragments(ipreviousmonth) = icount
407
408
409 fragments_sets( iraingagezone )%iRainGageZone = iraingagezone
410 fragments_sets( iraingagezone )%iStartRecord(imonth) = iindex
411 ! ! need to handle the last fragment set as a special case
412 ! FRAGMENTS_SETS( iRainGageZone )%iNumberOfFragments(iMonth) = iCount
413 icount = 1
414
415 else
416 icount = icount + 1
417 endif
418
419 ipreviousmonth = imonth
420 ipreviousraingagezone = iraingagezone
421
422 enddo
423
424 ! need to handle the last month of the last fragment set as a special case
425 fragments_sets( iraingagezone )%iNumberOfFragments(imonth) = icount
426
427 call logs%write("### Summary of fragment sets in memory ###", &
428 iloglevel=log_all, ilinesbefore=1, ilinesafter=1, lecho=false )
429 call logs%write("gage number | month | start index | num records ")
430 call logs%write("----------- | ---------- | ------------ | ------------")
431 do iindex=1, ubound( fragments_sets, 1)
432 do imonth=1,12
433 write (sbuf0, fmt="(i10)") iindex
434 write (sbuf1, fmt="(i10)") imonth
435 write (sbuf2, fmt="(i12)") fragments_sets(iindex)%iStartRecord(imonth)
436 write (sbuf3, fmt="(i10)") fragments_sets(iindex)%iNumberOfFragments(imonth)
437 write (sbuf4, fmt="(a10,' | ', a10,' | ', a12,' | ',a10)") adjustl(sbuf0), &
438 adjustl(sbuf1), adjustl(sbuf2), adjustl(sbuf3)
439 call logs%write( sbuf4 )
440 enddo
441 end do
442
443 end subroutine process_fragment_sets
444
445!--------------------------------------------------------------------------------------------------
446
447 !> eliminate rainfall on the 29th day of February; bump up all other values to ensure sum = 1
449
450 integer (c_int), intent(in) :: iCount
451
452 ! [ LOCALS ]
453 real (c_float) :: sum_fragments
454
455 ! we only want to correct the fragment if it was actually generated during a
456 ! leap year
457 if (fragments(icount)%fFragmentValue(29) > 0.0_c_float) then
458 sum_fragments = sum( fragments(icount)%fFragmentValue(1:28) )
459
460 ! bump up all February fragments so that their 28-day sum is 1
461 fragments(icount)%fFragmentValue(1:28) = fragments(icount)%fFragmentValue(1:28) &
462 / sum_fragments
463 ! zero out remaining values
464 fragments(icount)%fFragmentValue(29:31) = 0.0_c_float
465 endif
466
468
469!--------------------------------------------------------------------------------------------------
470
471 ! in order to compare Hawaii Water Budget results to SWB results,
472 ! it is necessary to force SWB to use the same sequence of fragments
473 ! that was used in the HWB simulations
474 subroutine read_fragments_sequence( sFilename )
475
476 character (len=*), intent(in) :: sFilename
477
478 ! [ LOCALS ]
479 character (len=512) :: sRecord, sSubstring
480 integer (c_int) :: iStat
481 integer (c_int) :: iCount
482 integer (c_int) :: iIndex
483 integer (c_int) :: iNumLines
484 type (ASCII_FILE_T) :: SEQUENCE_FILE
485 character (len=10) :: sBuf0
486 character (len=10) :: sBuf1
487 character (len=12) :: sBuf2
488 character (len=10) :: sBuf3
489 character (len=10) :: sBuf4
490 character (len=256) :: sBuf5
491 character (len=256) :: error_str
492 type (FSTRING_LIST_T) :: slHeader
493 integer (c_int) :: max_rain_gage_number
494 integer (c_int) :: max_simulation_number
495
496
497 call sequence_file%open( sfilename = sfilename, &
498 scommentchars = "#%!", &
499 sdelimiters = "WHITESPACE", &
500 lhasheader = .true._c_bool )
501
502 slheader = sequence_file%readHeader()
503
504 inumlines = sequence_file%numLines()
505
506 allocate( fragments_sequence( inumlines ), stat=istat )
507 call assert( istat == 0, "Problem allocating memory for fragments sequence table", &
508 __file__, __line__ )
509
510 icount = 0
511
512 do
513
514 ! read in next line of file
515 srecord = sequence_file%readLine()
516
517 if ( sequence_file%isEOF() ) exit
518
519 icount = icount + 1
520
521 ! read in simulation number
522 call chomp(srecord, ssubstring, sequence_file%sDelimiters )
523
524 if ( len_trim(ssubstring) == 0 ) &
525 call die( "Missing simulation number in the fragments sequence file", &
526 __file__, __line__, "Problem occured on line number " &
527 //ascharacter(sequence_file%currentLineNum() ) &
528 //" of file "//dquote(sfilename) )
529
530 fragments_sequence(icount)%sim_number = asint(ssubstring)
531
532 ! read in month
533 call chomp(srecord, ssubstring, sequence_file%sDelimiters )
534
535 if ( len_trim(ssubstring) == 0 ) &
536 call die( "Missing month number in the fragments sequence file", &
537 __file__, __line__, "Problem occured on line number " &
538 //ascharacter(sequence_file%currentLineNum() ) &
539 //" of file "//dquote(sfilename) )
540
541 fragments_sequence(icount)%sim_month = asint(ssubstring)
542
543 ! read in rainfall zone
544 call chomp(srecord, ssubstring, sequence_file%sDelimiters )
545
546 if ( len_trim(ssubstring) == 0 ) &
547 call die( "Missing rainfall zone number in the fragments sequence file", &
548 __file__, __line__, "Problem occured on line number " &
549 //ascharacter(sequence_file%currentLineNum() ) &
550 //" of file "//dquote(sfilename) )
551
552 fragments_sequence(icount)%sim_rainfall_zone = asint(ssubstring)
553
554
555 ! read in sim_year
556 call chomp(srecord, ssubstring, sequence_file%sDelimiters )
557
558 if ( len_trim(ssubstring) == 0 ) &
559 call die( "Missing year number in the fragments sequence file", &
560 __file__, __line__, "Problem occured on line number " &
561 //ascharacter(sequence_file%currentLineNum() ) &
562 //" of file "//dquote(sfilename) )
563
564 fragments_sequence(icount)%sim_year = asint(ssubstring)
565
566 ! read in sim_random_number
567 call chomp(srecord, ssubstring, sequence_file%sDelimiters )
568
569 if ( len_trim(ssubstring) == 0 ) &
570 call die( "Missing simulation random number in the fragments sequence file", &
571 __file__, __line__, "Problem occured on line number " &
572 //ascharacter(sequence_file%currentLineNum() ) &
573 //" of file "//dquote(sfilename) )
574
575 fragments_sequence(icount)%sim_random_number = asfloat(ssubstring)
576
577 ! read in simulation selected set (fragment set selected by HWB)
578 call chomp(srecord, ssubstring, sequence_file%sDelimiters )
579
580 if ( len_trim(ssubstring) == 0 ) &
581 call die( "Missing selected fragment set number in the fragments sequence file", &
582 __file__, __line__, "Problem occured on line number " &
583 //ascharacter(sequence_file%currentLineNum() ) &
584 //" of file "//dquote(sfilename) )
585
586 fragments_sequence(icount)%sim_selected_set = asint(ssubstring)
587
588 enddo
589
590 max_rain_gage_number = maxval(fragments_sequence(:)%sim_rainfall_zone,1)
591 max_simulation_number = maxval(fragments_sequence(:)%sim_number,1)
592
593 ! idea here is that if we are reading in a sequence file, there is a good chance
594 ! that the users is running multiple simulations; this will allow for a
595 ! separate pointer to be established for each rain gage/simulation number combination
596 ! as well as a means to keep the random value sequences separate by simulation
597
598 if ( allocated(current_fragments)) deallocate(current_fragments, stat=istat, &
599 errmsg=error_str)
600 call assert( istat == 0, "Problem deallocating memory, stat="//ascharacter(istat) &
601 //"; msg: "//trim(error_str), __file__, __line__ )
602
603 allocate(current_fragments(max_rain_gage_number, max_simulation_number), stat=istat, &
604 errmsg=error_str)
605 call assert( istat == 0, "Problem allocating memory, stat="//ascharacter(istat) &
606 //"; msg: "//trim(error_str), __file__, __line__ )
607
608
609 if ( allocated(random_values)) deallocate(random_values, stat=istat, errmsg=error_str)
610 call assert( istat == 0, "Problem deallocating memory, stat="//ascharacter(istat) &
611 //"; msg: "//trim(error_str), __file__, __line__ )
612
613 allocate(random_values(max_rain_gage_number, max_simulation_number), stat=istat, &
614 errmsg=error_str )
615 call assert( istat == 0, "Problem allocating memory, stat="//ascharacter(istat) &
616 //"; msg: "//trim(error_str), __file__, __line__ )
617
618 call logs%write("### Summary of fragment sequence sets in memory ###", &
619 iloglevel=log_debug, ilinesbefore=1, ilinesafter=1, lecho=false )
620 call logs%write("sim number | rainfall zone | month | year | selected set ")
621 call logs%write("----------- | ---------- | ------------ | ------------|------------")
622 do iindex=1, ubound( fragments_sequence, 1)
623 write (sbuf0, fmt="(i10)") fragments_sequence( iindex )%sim_number
624 write (sbuf1, fmt="(i10)") fragments_sequence( iindex )%sim_rainfall_zone
625 write (sbuf2, fmt="(i12)") fragments_sequence( iindex )%sim_month
626 write (sbuf3, fmt="(i10)") fragments_sequence( iindex )%sim_year
627 write (sbuf4, fmt="(i10)") fragments_sequence( iindex )%sim_selected_set
628
629 write (sbuf5, fmt="(a,' | ', a,' | ', a,' | ',a,' | ',a)") &
630 adjustl(sbuf0), adjustl(sbuf1), adjustl(sbuf2), adjustl(sbuf3), adjustl(sbuf4)
631 call logs%write( trim( sbuf5 ) )
632 end do
633
634 end subroutine read_fragments_sequence
635
636!--------------------------------------------------------------------------------------------------
637
638 !> Update rainfall fragments on daily basis.
639 !!
640 !! If called when lShuffle is TRUE:
641 !! 1) update random values
642 !! 2) random values are used to select the next active fragment set
643 !! for the current RainGageZone
644 !!
645 !! *Each* time the routine is called, the appropriate fragment is
646 !! selected from the current active fragment set and is assigned
647 !! to all cells that share a common RainGageZone
648
649 subroutine update_fragments( lShuffle )
650
651 logical (c_bool), intent(in) :: lShuffle
652
653 ! [ LOCALS ]
654 integer (c_int) :: rain_zone
655 integer (c_int) :: iMaxRainZones
656 integer (c_int) :: iMonth
657 integer (c_int) :: iDay
658 integer (c_int) :: iYearOfSimulation
659
660 integer (c_int) :: iNumberOfFragments
661 integer (c_int) :: iStartRecord
662 integer (c_int) :: iEndRecord
663 integer (c_int) :: iTargetRecord
664 integer (c_int) :: iStat
665 integer (c_int) :: iUBOUND_FRAGMENTS
666 integer (c_int) :: iUBOUND_CURRENT_FRAGMENTS
667 character (len=512) :: sBuf
668
669
670 imaxrainzones = maxval(fragments%iRainGageZone)
671 imonth = sim_dt%curr%iMonth
672 iday = sim_dt%curr%iDay
673 iyearofsimulation=sim_dt%iYearOfSimulation
674
675 ! equal the number of fragments are in memory
676 iubound_fragments = ubound( fragments, 1)
677
678 ! should equal the number of rainfall gages in model domain
679 iubound_current_fragments = ubound( current_fragments, 1)
680
681 ! if by chance a mismatch in shape-to-grid results in an active cell with *NO* valid
682 ! rain gage ID, we need to set the entire array to zero to quash any spurious values getting in
683 fragment_value = 0.0_c_float
684
685 do rain_zone = 1, imaxrainzones
686
687 if ( lshuffle ) then
688 ! find next fragment *record*
689
690 ! update the module variable RANDOM_VALUES
692
693 istartrecord = fragments_sets( rain_zone )%iStartRecord(imonth)
694 inumberoffragments = fragments_sets(rain_zone)%iNumberOfFragments(imonth)
695 iendrecord = istartrecord + inumberoffragments - 1
696 itargetrecord = istartrecord &
697 + int(random_values(rain_zone, simulation_number) &
698 * real( inumberoffragments ))
699
700 if ( ( rain_zone > iubound_current_fragments ) .or. ( itargetrecord > iubound_fragments ) &
701 .or. ( rain_zone < 1 ) .or. ( itargetrecord < 1) ) then
702 call logs%write("Error detected in method of fragments routine; dump of current" &
703 //" variables follows:", ilinesbefore=1)
704 call logs%write("rain_zone : "//ascharacter(rain_zone), itab=3 )
705 call logs%write("simulation_number : "//ascharacter(rain_zone), itab=3 )
706 call logs%write("iStartRecord : "//ascharacter(istartrecord), itab=3 )
707 call logs%write("iNumberOfFragments: "//ascharacter(inumberoffragments), itab=3 )
708 call logs%write("iEndRecord : "//ascharacter(iendrecord), itab=3 )
709 call logs%write("iTargetRecord : "//ascharacter(itargetrecord), itab=3 )
710 call logs%write("ubound(CURRENT_FRAGMENTS, 1): "//ascharacter(iubound_current_fragments), &
711 itab=3 )
712 call logs%write("ubound(FRAGMENTS, 1): "//ascharacter(iubound_fragments), itab=3 )
713 call logs%write("RANDOM_VALUES(rain_zone,SIMULATION_NUMBER): " &
714 //ascharacter(random_values(rain_zone,simulation_number)), itab=3 )
715 call die( "Miscalculation in target record: calculated record index is out of bounds", &
716 __file__, __line__ )
717 endif
718
719 ! reassign fragment pointer for this rain zone to the newly selected record
720 current_fragments(rain_zone, simulation_number)%pFragment => fragments( itargetrecord )
721
722 write(lu_fragments_echo,fmt="(4(i5,','),f10.6,',',i5,',',30(f8.3,','),f8.3)") &
724 fragments( itargetrecord)%iMonth, &
725 fragments( itargetrecord)%iRainGageZone, &
726 iyearofsimulation, &
727 random_values(rain_zone, simulation_number), &
728 fragments( itargetrecord)%iFragmentSet, &
729 fragments( itargetrecord)%fFragmentValue
730
731 ! call LOGS%write( trim(sBuf), iLogLevel=LOG_DEBUG, lEcho=FALSE )
732
733 endif
734
735 if ( ( current_fragments( rain_zone, simulation_number )%pFragment%fFragmentValue( iday ) < 0.0 ) &
736 .or. ( current_fragments( rain_zone, simulation_number )%pFragment%fFragmentValue( iday ) > 1.0 ) ) then
737
738 call logs%write("Error detected in method of fragments routine; dump of current variables" &
739 //" follows:", ilinesbefore=1, iloglevel=log_all )
740 call logs%write("rain_zone:"//ascharacter(rain_zone), itab=3 )
741 call logs%write("iDay: "//ascharacter(iday), itab=3 )
742 call logs%write("iRainGageZone: "//ascharacter(fragments( itargetrecord)%iRainGageZone), itab=3 )
743 call logs%write("iFragmentSet: "//ascharacter(fragments( itargetrecord)%iFragmentSet), itab=3 )
744 call logs%write("fFragmentValue: "//ascharacter(fragments( itargetrecord)%fFragmentValue(iday) ), itab=3 )
745
746 endif
747
748 ! call LOGS%write("frag: "//asCharacter(rain_zone)//" day: "//asCharacter(iDay) &
749 ! //" value: "//asCharacter( CURRENT_FRAGMENTS( rain_zone )%pFragment%fFragmentValue( iDay ) ), &
750 ! lEcho=FALSE )
751
752 ! now place current days' fragment value into the matching cells
753 where ( rain_gage_id == rain_zone )
754
755 fragment_value = current_fragments( rain_zone, simulation_number )%pFragment%fFragmentValue( iday )
756
757 endwhere
758
759 enddo
760
761 end subroutine update_fragments
762
763!--------------------------------------------------------------------------------------------------
764
766
767 ! [ LOCALS ]
768 integer (c_int) :: iIndex, iIndex2
769 logical (c_bool) :: lSequenceSelection
770
771 if ( random_fragment_sequences ) then
772
773 do iindex2=1,size(random_values,1)
774
775 !call random_number( RANDOM_VALUES )
777
778 enddo
779
780 else
781
782 random_values(:,simulation_number) = -9999999.9
783
784 do iindex=1, size(fragments_sequence%sim_month, 1)
785
786 lsequenceselection = ( fragments_sequence(iindex)%sim_month == sim_dt%curr%iMonth ) &
787 .and. ( fragments_sequence(iindex)%sim_year == sim_dt%iYearOfSimulation ) &
788 .and. ( fragments_sequence(iindex)%sim_number == simulation_number )
789
790 if ( .not. lsequenceselection ) cycle
791
792 do iindex2=1,size(random_values,1)
793
794 if ( fragments_sequence( iindex )%sim_rainfall_zone == iindex2 ) then
795
796 random_values( iindex2, simulation_number ) = fragments_sequence( iindex )%sim_random_number
797 exit
798
799 endif
800
801 enddo
802
803 enddo
804
805 endif
806
807 if (any( random_values(:, simulation_number) < 0.0 ) ) then
808
809 call logs%write("Error detected in method of fragments routine - random values " &
810 //"not found in sequence file for rainfall zone(s):", ilinesbefore=1)
811 do iindex=1,size(random_values, 1)
812 if ( random_values(iindex, simulation_number) < 0.0 ) &
813 call logs%write("simulation number, rainfall zone: " &
814 //trim(ascharacter(simulation_number))//", "//trim(ascharacter(iindex)), itab=3 )
815 enddo
816
817 endif
818
819 end subroutine update_random_values
820
821!--------------------------------------------------------------------------------------------------
822
824
825 logical (c_bool), intent(in) :: lactive(:,:)
826
827 ! [ LOCALS ]
828 integer (c_int) :: iindex
829 integer (c_int) :: imaxrainzones
830 integer (c_int) :: istat
831 logical (c_bool), save :: lfirstcall = true
832
834
835
836 !! if it is the first day of the month, update the rainfall adjustment factor grid
837 !! and update the fragments
838 if ( sim_dt%curr%iDay == 1 .or. lfirstcall ) then
839
840 ! locate the data structure associated with the gridded rainfall adjustment factor
841 prainfall_adjust_factor => dat%find("RAINFALL_ADJUST_FACTOR")
842
843 if ( .not. associated(prainfall_adjust_factor) ) &
844 call die("A RAINFALL_ADJUST_FACTOR grid must be supplied in order to make use" &
845 //" of this option.", __file__, __line__)
846
847 call prainfall_adjust_factor%getvalues( sim_dt%curr )
848
849 ! map the 2D array of RAINFALL_ADJUST_FACTOR values to the vector of active cells
850 rainfall_adjust_factor = pack( prainfall_adjust_factor%pGrdBase%rData, lactive )
851
852 call update_fragments( lshuffle = true)
853 lfirstcall = false
854
855 else
856
857 call update_fragments( lshuffle = false )
858
859 endif
860
862
863!--------------------------------------------------------------------------------------------------
864
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
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.
type(dict_t), public cf_dict
subroutine, public die(smessage, smodule, iline, shints, scalledby, icalledbyline)
Provides support for input and output of gridded ASCII data, as well as for creation and destruction ...
Definition grid.F90:8
type(logfile_t), public logs
Definition logfiles.F90:62
Module precipitation__method_of_fragments provides support for creating synthetic daily precipitation...
type(fragments_sequence_t), dimension(:), allocatable, public fragments_sequence
Array of fragment sequence sets.
real(c_float), dimension(:), allocatable, public rainfall_adjust_factor
Module variable that holds the rainfall adjustment factor.
subroutine, public precipitation_method_of_fragments_calculate(lactive)
subroutine update_fragments(lshuffle)
Update rainfall fragments on daily basis.
type(fragments_sequence_t), pointer pfragments_sequence
Pointer to all or some of the FRAGMENTS_SEQUENCE array.
subroutine, public precipitation_method_of_fragments_initialize(lactive)
Initialize method of fragments.
type(data_catalog_entry_t), pointer prainfall_adjust_factor
logical(c_bool) random_fragment_sequences
Module variable detemining whether fragment sequences are chosen at random or selected from an extern...
type(ptr_fragments_t), dimension(:,:), allocatable current_fragments
Subset of rainfall fragments file pointing to the currently active fragments.
subroutine process_fragment_sets()
after fragments file has been read in, iterate over a set of rainfall fragments
integer(c_int), public simulation_number
Module variable indicating which "simulation number" is active Only has meaning if the rainfall fragm...
integer(c_int), dimension(:), allocatable, public rain_gage_id
Module variable that holds the rainfall gage (zone) number.
real(c_double), dimension(:,:), allocatable random_values
Module variable that holds a sequence of random numbers associated with the selection of the fragment...
type(fragments_set_t), dimension(:), allocatable, public fragments_sets
Array of fragments sets; fragments sets include indices to the start record associated with the fragm...
real(c_float), dimension(:), allocatable, public fragment_value
Module variable that holds the current day's rainfall fragment value.
logical(c_bool), dimension(:), allocatable sequence_selection
Module level variable used to create subsets of the FRAGMENT_SEQUENCES file.
type(fragments_t), dimension(:), allocatable, target, public fragments
Array of all fragments read in from the rainfall fragments file.
subroutine normalize_february_fragment_sequence(icount)
eliminate rainfall on the 29th day of February; bump up all other values to ensure sum = 1
type(date_range_t), public sim_dt
Data structure to hold static (pre-calculated) fragment selection numbers.
Data structure to hold the current active rainfall fragments for a particular rain gage zone.
Data structure that holds a single line of data from the input rainfall fragments file.