Merge remote-tracking branch 'origin/release-v4.5'
[WRF.git] / hydro / Routing / Reservoirs / reservoir_tests.F
blob934e49116902f7d1e26e74c13590471723c5ed67
1 ! This program unit tests the various reservoir implementations along with edge cases.
2 ! It is important to run these unit tests whenever making any changes to code in this
3 ! module to ensure nothing is broken. If nothing is broken, the user will see
4 ! "All Reservoir Tests Passed". To compile and run these tests, go to the Reservoir
5 ! directory in the terminal and type "make" and then "make test". Then type "./reservoir_tests".
7 program reservoir_unit_tests
8     use module_levelpool_tests
9     use module_persistence_levelpool_hybrid_tests
10     use module_rfc_forecasts_tests
12     implicit none
14     logical :: rv1 = .false.
15     logical :: rv2 = .false.
16     logical :: rv3 = .false.
17     logical :: rv4 = .false.
18     logical :: rv5 = .false.
19     logical :: rv6 = .false.
20     logical :: rv7 = .false.
21     logical :: rv8 = .false.
22     logical :: rv9 = .false.
23     logical :: rv10 = .false.
24     logical :: rv11 = .false.
25     logical :: rv12 = .false.
26     logical :: rv13 = .false.
27     logical :: rv14 = .false.
28     logical :: rv15 = .false.
29     logical :: rv16 = .false.
30     logical :: rv17 = .false.
32     real, dimension(120) :: inflow_array
34     inflow_array = (/189.22899, 189.27005, 189.31049, 189.35042, 189.38965, 189.42819, 189.46588, 189.50273, &
35     189.53859, 189.57346, 189.60719, 189.63979, 189.6711, 189.7011, 189.72968, &
36     189.75679, 189.7823, 189.80617, 189.82822, 189.84842, 189.86653, 189.88255, 189.89622, &
37     189.90752, 189.91612, 189.922, 189.92482, 189.92447, 189.92067, 189.91319, 189.90175, &
38     189.88611, 189.86592, 189.84088, 189.81064, 189.77487, 189.73317, 189.6852, 189.63051, &
39     189.56873, 189.49939, 189.42207, 189.33635, 189.24176, 189.13782, 189.02408, &
40     188.90009, 188.76535, 188.61945, 188.46188, 188.29224, 188.11006, 187.91493, 187.70644, &
41     187.48419, 187.24779, 186.9969, 186.73119, 186.45035, 186.15407, 185.84213, 185.51424, &
42     185.17023, 184.80989, 184.43312, 184.03975, 183.62973, 183.20296, 182.75943, 182.29909, &
43     181.82205, 181.32828, 180.81792, 80.29099, 179.74774, 179.1882, 178.61267, 178.02129, &
44     177.41437, 176.79207, 176.15475, 175.50269, 174.83627, 174.15576, 173.46162, &
45     172.75417, 172.03389, 171.3011, 170.55634, 169.79997, 169.03255, 168.25441, 167.46616, &
46     166.66815, 165.86099, 165.04509, 164.22101, 163.38913, 162.55011, 161.70428, 160.85229, &
47     159.99452, 159.13156, 158.26382, 157.39188, 156.51611, 155.63715, 154.75531, 153.8712, 152.98517, &
48     152.09779, 151.2094, 150.32057, 149.43166, 148.54315, 147.6554, 146.76892, 145.88405, 145.00128, 144.12091/)
50     rv1 = test_levelpool()
52     rv2 = test_levelpool_overflow_max_height()
54     rv3 = test_persistence_levelpool_hybrid_usgs()
56     rv4 = test_persistence_levelpool_hybrid_usace()
58     rv5 = test_rfc_forecasts_object()
60     rv6 = test_rfc_forecasts_time_series_object()
62     rv7 = test_rfc_forecasts_levelpool_fallback()
64     rv8 = test_rfc_forecasts_time_series_output_with_lookback_and_offset()
66     rv9 = test_rfc_forecasts_time_series_output_with_negative_data()
68     rv10 = test_rfc_forecasts_time_series_output_all_synthetic_data()
70     rv11 = test_rfc_forecasts_data_exceeding_max_range()
72     rv12 = test_rfc_forecasts_with_offset_for_extended_AnA()
74     rv13 = test_persistence_levelpool_hybrid_usace_overtop()
76     rv14 = test_persistence_levelpool_hybrid_usace_no_timeslice_available()
78     rv15 = test_ak_rfc_forecasts_pass_through_fallback()
80     rv16 = test_ak_rfc_forecasts_time_series_with_lookback_and_offset()
82     rv17 = test_rfc_forecasts_over_max_water_elevation()
84     if (rv1 .and. rv2 .and. rv3 .and. rv4 .and. rv5 .and. rv6 .and. rv7 .and. rv8 .and. rv9 .and. &
85     rv10 .and. rv11 .and. rv12 .and. rv13 .and. rv14 .and. rv15 .and. rv16 .and. rv17) then
86         print *, "========================================================================"
87         print *, 'All Reservoir Tests Passed'
88         print *, "========================================================================"
90     else
91         print *, "========================================================================"
92         print *, 'Not All Reservoir Tests Passed'
93         print *, "========================================================================"
94     end if
96     contains
98     !------------------------------------------------------------------------------!
99     !                              test_levelpool()                                !
100     ! this function verifies that the constructor for the levelpool type correctly !
101     ! initializes all data members                                                 !
102     !------------------------------------------------------------------------------!
104     function test_levelpool() result(rv)
105         implicit none
107         logical rv                        ! test result
108         logical :: call_status
110         type (levelpool) :: levelpool_reservoir_data
111         real :: water_elevation = 2.
113         integer(kind=int64):: lake_number = 20
115         call_status = .false.
117         print *, "calling init for levelpool"
118         call levelpool_reservoir_data%init(water_elevation, 4., 6., 8., 10., 11., 12., 14., 16., 18., lake_number)
120         print *, "testing data in levelpool"
121         call_status = levelpool_data_info(levelpool_reservoir_data)
123         if (call_status) then
124             rv = .true.
125         end if
127     end function test_levelpool
130     ! This tests the Persistence Levelpool Hybrid Module run reservoir function for USGS reservoirs.
131     ! It also reads the USGS persistence parameters from the provided reservoir parameter file and
132     ! a gage discharge from the provided USGS Timeslice file.
133     function test_persistence_levelpool_hybrid_usgs() result(rv)
134         implicit none
135         logical rv, rv1, rv2                        ! test result
136         type (persistence_levelpool_hybrid) :: persistence_levelpool_hybrid_reservoir_data
137         real :: outflow, inflow
138         real :: water_elevation
139         real :: prev_time_inflow
140         real :: lake_area, weir_elevation, weir_coefficient
141         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
142         real :: orifice_area, max_depth, initial_fractional_depth
143         integer(kind=int64):: lake_number
144         integer :: reservoir_type
145         integer :: timestep_count
146         character*256 :: cwd_full
147         integer :: dynamic_reservoir_type
148         real :: assimilated_value
149         character(len=256) :: assimilated_source_file
151         prev_time_inflow = 0.0
152         timestep_count = 0
153         water_elevation = 0.0
154         rv = .false.
155         rv1 = .false.
156         rv2 = .false.
158         lake_area = 2.096320037841796875e+02
159         weir_elevation = 1.332074047851562455e+03
160         weir_coefficient = 4.000000000000000222e-01
161         weir_length = 1.000000000000000000e+01
162         dam_length = 10.0
163         orifice_elevation = 1.314473347981770758e+03
164         orifice_coefficient = 1.000000000000000056e-01
165         orifice_area = 1.0
166         max_depth = 1.335180053710937500e+03
167         initial_fractional_depth = 8.999999761581420898e-01
168         lake_number = 402142
169         reservoir_type = 2
171         print *, "calling init for persistence_levelpool_hybrid USGS type reservoir"
173         cwd_full = "../../../../tests/local/reservoir_testing_files/"
175         call persistence_levelpool_hybrid_reservoir_data%init(water_elevation, lake_area, weir_elevation, &
176         weir_coefficient, weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, &
177         initial_fractional_depth, lake_number, reservoir_type, "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", &
178         "2010-10-01_07:00:00", cwd_full, cwd_full, 12, 1000000000)
180         print *, "testing data in persistence_levelpool_hybrid"
181         rv1 = persistence_levelpool_hybrid_data_info(persistence_levelpool_hybrid_reservoir_data)
183         print *, "calling reservoir run for persistence_levelpool_hybrid USGS type reservoir"
185         water_elevation = 1331.18005
187         do timestep_count = 1, 120
189             inflow = inflow_array(timestep_count)
190             call persistence_levelpool_hybrid_reservoir_data%run(inflow, &
191             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
192             assimilated_value, assimilated_source_file)
194             prev_time_inflow = inflow
196             print *, outflow
197         end do
199         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
201         if (outflow .ge. 13.73367 - epsilon(13.73367) .and. outflow .le. 13.73367 + epsilon(13.73367) &
202         .and. dynamic_reservoir_type == 2 .and. assimilated_value .ge. 13.73367 - epsilon(13.73367) &
203         .and. assimilated_value .le. 13.73367 + epsilon(13.73367) .and. &
204         assimilated_source_file == "2010-10-01_06:00:00.15min.usgsTimeSlice.ncdf") then
205             rv2 = .true.
206             print *, "========================================================================"
207             print *, 'Persistence Levelpool Hybrid Run USGS Reservoir Test Passed'
208             print *, "========================================================================"
209         else
210             print *, "========================================================================"
211             print *, 'Persistence Levelpool Hybrid Run USGS reservoir Test Failed'
212             print *, 'Outflow should be 13.7336712'
213             print *, 'dynamic_reservoir_type needs to be 2 for USGS persistence output'
214             print *, "========================================================================"
215             print *, outflow
216         end if
218         if (rv1 .and. rv2) then
219             rv = .true.
220         end if
222     end function test_persistence_levelpool_hybrid_usgs
225     ! This tests the Persistence Levelpool Hybrid Module run reservoir function for U.S. Army Corps of Engineers (USACE)
226     ! type reservoirs. It also reads the USACE persistence parameters from the provided reservoir parameter file
227     ! and a gage discharge from the provided USACE Timeslice file.
228     function test_persistence_levelpool_hybrid_usace() result(rv)
229         implicit none
230         logical rv, rv1, rv2                        ! test result
231         type (persistence_levelpool_hybrid) :: persistence_levelpool_hybrid_reservoir_data
232         real :: outflow, inflow
233         real :: water_elevation
234         real :: prev_time_inflow
235         real :: lake_area, weir_elevation, weir_coefficient
236         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
237         real :: orifice_area, max_depth, initial_fractional_depth
238         integer(kind=int64):: lake_number
239         integer :: reservoir_type
240         integer :: timestep_count
241         character*256 :: cwd_full
242         integer :: dynamic_reservoir_type
243         real :: assimilated_value
244         character(len=256) :: assimilated_source_file
246         prev_time_inflow = 0.0
247         timestep_count = 0
248         water_elevation = 0.0
249         rv = .false.
250         rv2 = .false.
252         lake_area = 2.096320037841796875e+02
253         weir_elevation = 1.332074047851562455e+03
254         weir_coefficient = 4.000000000000000222e-01
255         weir_length = 1.000000000000000000e+01
256         dam_length = 10.0
257         orifice_elevation = 1.314473347981770758e+03
258         orifice_coefficient = 1.000000000000000056e-01
259         orifice_area = 1.0
260         max_depth = 1.335180053710937500e+03
261         initial_fractional_depth = 8.999999761581420898e-01
262         lake_number = 4672717
263         reservoir_type = 3
265         print *, "calling init for persistence_levelpool_hybrid USACE type reservoir"
267         cwd_full = "../../../../tests/local/reservoir_testing_files/"
269         call persistence_levelpool_hybrid_reservoir_data%init(water_elevation, lake_area, weir_elevation, &
270         weir_coefficient, weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, &
271         max_depth, initial_fractional_depth, lake_number, reservoir_type, &
272         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", &
273         "2016-10-10_02:00:00", cwd_full, cwd_full, 12, 1000000000)
275         print *, "calling reservoir run for persistence_levelpool_hybrid for USACE type reservoir"
277         water_elevation = 1331.18005
279         do timestep_count = 1, 120
281             inflow = inflow_array(timestep_count)
282             call persistence_levelpool_hybrid_reservoir_data%run(inflow, &
283             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
284             assimilated_value, assimilated_source_file)
286             prev_time_inflow = inflow
288             print *, outflow
289         end do
291         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
293         if (outflow .ge. 16.0820713 - epsilon(16.0820713) .and. outflow .le. 16.0820713 + epsilon(16.0820713) &
294         .and. dynamic_reservoir_type == 3 .and. assimilated_value .ge. 16.082071 - epsilon(16.082071) &
295         .and. assimilated_value .le. 16.082071 + epsilon(16.082071) .and. &
296         assimilated_source_file == "2016-10-10_00:00:00.15min.usaceTimeSlice.ncdf") then
297             rv2 = .true.
298             print *, "========================================================================"
299             print *, 'Persistence Levelpool Hybrid Run USACE Reservoir Test Passed'
300             print *, "========================================================================"
301         else
302             print *, "========================================================================"
303             print *, 'Persistence Levelpool Hybrid Run USACE Reservoir Test Failed'
304             print *, 'Outflow should be 16.08207'
305             print *, 'dynamic_reservoir_type needs to be 3 for USACE persistence output'
306             print *, "========================================================================"
307             print *, outflow
308         end if
310         if (rv2) then
311             rv = .true.
312         end if
314     end function test_persistence_levelpool_hybrid_usace
317     ! This tests the Persistence Levelpool Hybrid Module run reservoir function USACE
318     ! type reservoirs to output level pool whenever it reaches the overtop condition.
319     function test_persistence_levelpool_hybrid_usace_overtop() result(rv)
321         implicit none
322         logical rv, rv1, rv2                        ! test result
323         type (persistence_levelpool_hybrid) :: persistence_levelpool_hybrid_reservoir_data
324         real :: outflow, inflow
325         real :: water_elevation
326         real :: prev_time_inflow
327         real :: lake_area, weir_elevation, weir_coeffecient
328         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
329         real :: orifice_area, max_depth, initial_fractional_depth
330         integer(kind=int64):: lake_number
331         integer :: reservoir_type
332         integer :: timestep_count
333         character*256 :: cwd_full
334         integer :: dynamic_reservoir_type
335         real :: assimilated_value
336         character(len=256) :: assimilated_source_file
338         prev_time_inflow = 0.0
339         timestep_count = 0
340         water_elevation = 0.0
341         rv = .false.
342         rv2 = .false.
344         lake_area = 2.096320037841796875e+02
345         weir_elevation = 1.332074047851562455e+03
346         weir_coeffecient = 4.000000000000000222e-01
347         weir_length = 1.000000000000000000e+01
348         dam_length = 10.0
349         orifice_elevation = 1.314473347981770758e+03
350         orifice_coefficient = 1.000000000000000056e-01
351         orifice_area = 1.0
352         max_depth = 1.335180053710937500e+03
353         initial_fractional_depth = 8.999999761581420898e-01
354         lake_number = 4672717
355         reservoir_type = 3
357         print *, "calling init for persistence_levelpool_hybrid USACE type reservoir"
359         cwd_full = "../../../../tests/local/reservoir_testing_files/"
361         call persistence_levelpool_hybrid_reservoir_data%init(water_elevation, lake_area, weir_elevation, &
362         weir_coeffecient, weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, &
363         max_depth, initial_fractional_depth, lake_number, reservoir_type, &
364         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", &
365         "2016-10-10_02:00:00", cwd_full, cwd_full, 12, 1000000000)
367         print *, "calling reservoir run for persistence_levelpool_hybrid for USACE type reservoir"
369         ! This initial water elevation sets to reservoir overtop
370         water_elevation = 1338.18005
372         do timestep_count = 1, 120
374             inflow = inflow_array(timestep_count)
375             call persistence_levelpool_hybrid_reservoir_data%run(inflow, &
376             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
377             assimilated_value, assimilated_source_file)
379             prev_time_inflow = inflow
381             print *, outflow
382         end do
384         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
386         if (outflow .ge. 229.226059 - epsilon(229.226059) .and. outflow .le. 229.226059 + epsilon(229.226059) &
387         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
388         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
389         assimilated_source_file == "") then
390             rv2 = .true.
391             print *, "========================================================================"
392             print *, 'Persistence Levelpool Hybrid Run USACE Reservoir Overflow Max Height Test Passed'
393             print *, "========================================================================"
394         else
395             print *, "========================================================================"
396             print *, 'Persistence Levelpool Hybrid Run USACE Reservoir Overflow Max Height Test Failed'
397             print *, 'Outflow should be 16.08207'
398             print *, 'dynamic_reservoir_type needs to be 1 because reservoir outputs level pool instead of USACE persistence'
399             print *, "========================================================================"
400             print *, outflow
401         end if
403         if (rv2) then
404             rv = .true.
405         end if
407     end function test_persistence_levelpool_hybrid_usace_overtop
410     ! This tests the Persistence Levelpool Hybrid Module run reservoir function USACE
411     ! type reservoirs to output level pool whenever no timeslices are available.
412     function test_persistence_levelpool_hybrid_usace_no_timeslice_available() result(rv)
414         implicit none
415         logical rv, rv1, rv2                        ! test result
416         type (persistence_levelpool_hybrid) :: persistence_levelpool_hybrid_reservoir_data
417         real :: outflow, inflow
418         real :: water_elevation
419         real :: prev_time_inflow
420         real :: lake_area, weir_elevation, weir_coeffecient
421         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
422         real :: orifice_area, max_depth, initial_fractional_depth
423         integer(kind=int64):: lake_number
424         integer :: reservoir_type
425         integer :: timestep_count
426         character*256 :: cwd_full
427         integer :: dynamic_reservoir_type
428         real :: assimilated_value
429         character(len=256) :: assimilated_source_file
431         prev_time_inflow = 0.0
432         timestep_count = 0
433         water_elevation = 0.0
434         rv = .false.
435         rv2 = .false.
437         lake_area = 2.096320037841796875e+02
438         weir_elevation = 1.332074047851562455e+03
439         weir_coeffecient = 4.000000000000000222e-01
440         weir_length = 1.000000000000000000e+01
441         dam_length = 10.0
442         orifice_elevation = 1.314473347981770758e+03
443         orifice_coefficient = 1.000000000000000056e-01
444         orifice_area = 1.0
445         max_depth = 1.335180053710937500e+03
446         initial_fractional_depth = 8.999999761581420898e-01
447         lake_number = 347987  ! the corresponding USACE gage id is CA00265
448         reservoir_type = 3
450         print *, "calling init for persistence_levelpool_hybrid USACE type reservoir"
452         cwd_full = "../../../../tests/local/reservoir_testing_files/"
454         call persistence_levelpool_hybrid_reservoir_data%init(water_elevation, lake_area, weir_elevation, &
455         weir_coeffecient, weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, &
456         max_depth, initial_fractional_depth, lake_number, reservoir_type, &
457         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", &
458         "2016-10-12_02:00:00", cwd_full, cwd_full, 12, 1000000000)
460         print *, "calling reservoir run for persistence_levelpool_hybrid for USACE type reservoir"
462         water_elevation = 1331.18005
464         do timestep_count = 1, 120
466             inflow = inflow_array(timestep_count)
467             call persistence_levelpool_hybrid_reservoir_data%run(inflow, &
468             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
469             assimilated_value, assimilated_source_file)
471             prev_time_inflow = inflow
473             print *, outflow
474         end do
476         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
478         if (outflow .ge. 1.81549609 - epsilon(1.81549609) .and. outflow .le. 1.81549609 + epsilon(1.81549609) &
479         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
480         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
481         assimilated_source_file == "") then
482             rv2 = .true.
483             print *, "========================================================================"
484             print *, 'Persistence Levelpool Hybrid Run USACE Reservoir No Timeslice Available Test Passed'
485             print *, "========================================================================"
486         else
487             print *, "========================================================================"
488             print *, 'Persistence Levelpool Hybrid Run USACE Reservoir No Timeslice Available Test Failed'
489             print *, 'Outflow should be 1.81549609'
490             print *, 'dynamic_reservoir_type needs to be 1 because reservoir outputs level pool instead of USACE persistence'
491             print *, "========================================================================"
492             print *, outflow
493         end if
495         if (rv2) then
496             rv = .true.
497         end if
499     end function test_persistence_levelpool_hybrid_usace_no_timeslice_available
502     ! This tests the reservoir function of the level pool module under the specific condition
503     ! where the water elevation reaches the max height.
504     function test_levelpool_overflow_max_height() result(rv)
506         implicit none
507         logical rv                       ! test result
508         type (levelpool) :: levelpool_reservoir_data
509         real :: lake_area, weir_elevation, weir_coefficient
510         real :: weir_length, dam_length, orifice_elevation
511         real :: orifice_coefficient, orifice_area, max_depth
512         integer(kind=int64):: lake_number
513         real :: inflow, prev_time_inflow, outflow, water_elevation
514         real, dimension(108) :: inflow_array_overflow
515         integer :: timestep_count
516         integer :: dynamic_reservoir_type
517         real :: assimilated_value
518         character(len=256) :: assimilated_source_file
520         rv = .false.
521         prev_time_inflow = 0.0
522         timestep_count = 0
523         water_elevation = 0.0
525         lake_area = 1.509490013122558594e+01
526         weir_elevation = 9.626000022888183238e+00
527         weir_coefficient = 0.4
528         weir_length = 1.000000000000000000e+01
529         dam_length = 10.0
530         orifice_elevation = 7.733333269755045869e+00
531         orifice_coefficient = 1.000000000000000056e-01
532         orifice_area = 1.0
533         max_depth = 9.960000038146972656e+00
534         lake_number = 16944276
536         inflow_array_overflow = (/91.27196, 91.7394, 92.15904, 92.1518, 91.84663, &
537         91.38554, 90.86131, 90.32736, 89.81273, 89.3325, 88.89427, 88.5025, 88.16228, &
538         87.41539, 86.80043, 86.03979, 85.3849, 85.33451, 86.84274, 91.6084, 101.81398, &
539         118.85916, 143.99232, 177.7355, 219.2348, 267.22351, 319.90402, 374.54324, 428.86066, &
540         480.92096, 529.23584, 572.77673, 610.93237, 643.4389, 670.28516, 691.67767, 707.96088, &
541         719.57312, 726.96997, 730.63269, 731.03186, 728.61438, 723.79578, 716.9549, 708.43268, &
542         698.53247, 687.52112, 675.63123, 663.06421, 649.99976, 636.57898, 622.92926, 609.1745, &
543         595.40369, 581.68799, 568.08588, 554.64484, 541.4032, 528.39185, 515.63513, 503.14838, &
544         490.95123, 479.05109, 467.45493, 456.16663, 445.18753, 434.51706, 424.15311,414.0921, &
545         404.32956, 394.86014, 385.67789, 376.77621, 368.14966, 359.78958, 351.68875, 343.83972, &
546         336.23505, 328.86719, 321.7287, 314.81219, 308.11047, 301.61646, 295.32312, 289.22369, &
547         283.31207, 277.5813, 272.02521, 266.63776, 261.41315, 256.34564, 251.42978, 246.66023, &
548         242.03192, 237.53989, 233.17944, 228.94595, 224.83511, 220.84265, 216.96449, 213.19672, &
549         209.53554, 205.97734, 202.51857, 199.1559, 195.88605, 192.70595, 189.61255 /)
551         call levelpool_reservoir_data%init(water_elevation, lake_area, weir_elevation, &
552         weir_coefficient, weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, lake_number)
554         water_elevation = 9.73733330
556         print *, "outflow"
558         do timestep_count = 1, 108
559             inflow = inflow_array_overflow(timestep_count)
560             call levelpool_reservoir_data%run(inflow, inflow, 0.0, water_elevation, outflow, 300.0, dynamic_reservoir_type, &
561             assimilated_value, assimilated_source_file)
563             prev_time_inflow = inflow
565             print *, outflow
566         end do
568         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
570         if (outflow .ge. 17.0451450 - epsilon(17.0451450) .and. outflow .le. 17.0451450 + epsilon(17.0451450) &
571         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
572         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
573         assimilated_source_file == "") then
574             rv = .true.
575             print *, "========================================================================"
576             print *, 'Levelpool Overflow Max Height Test Passed'
577             print *, "========================================================================"
578         else
579             print *, "========================================================================"
580             print *, 'Levelpool Overflow Max Height Test Failed'
581             print *, 'Final outflow should be 17.0451450'
582              print *, 'dynamic_reservoir_type needs to be 1 for level pool type'
583             print *, "========================================================================"
584         end if
586     end function test_levelpool_overflow_max_height
589     ! This tests the basic object creation of an RFC Forecast type reservoir.
590     function test_rfc_forecasts_object() result(rv)
591         implicit none
593         logical rv                        ! test result
594         logical :: call_status
596         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
597         type (time_series_data_type) :: time_series_data
598         character*256 :: cwd_full
599         integer(kind=int64):: lake_number
601         real :: water_elevation = 2.
603         lake_number = 3745478
605         cwd_full = "../../../../tests/local/reservoir_testing_files/"
607         call_status = .false.
609         print *, "calling init for rfc_forecasts"
610         call rfc_forecasts_reservoir_data%init(water_elevation, 4., 6., 8., 10., 11., 12., 14., 16., 18., 0.9, lake_number, 4, &
611         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", "2019-08-18_07:00:00", cwd_full, 24)
613         print *, "testing data in rfc_forecasts"
615         call_status = rfc_forecasts_data_info(rfc_forecasts_reservoir_data)
617         if (call_status) then
618             rv = .true.
619         end if
621     end function test_rfc_forecasts_object
624     ! This tests the object creation and function of the RFC Forecasts Time Series Reader.
625     function test_rfc_forecasts_time_series_object() result(rv)
626         implicit none
628         logical rv                        ! test result
629         logical :: call_status
630         type (time_series_data_type) :: time_series_data, time_series_data_second
631         character*256 :: cwd_full
632         integer :: lookback_seconds, total_counts, observed_counts, time_step_seconds
633         integer :: timeslice_offset_hours
634         real, allocatable, dimension(:) :: forecast_discharges
635         logical :: forecast_found
636         logical :: all_discharges_synthetic
637         character(len=256):: time_series_file_name
638         integer :: lake_number
640         lake_number = 3745478
642         call_status = .false.
644         cwd_full = "../../../../tests/local/reservoir_testing_files/"
646         forecast_found = .false.
648         all_discharges_synthetic = .false.
650         timeslice_offset_hours = 3
652         call time_series_data%init("2019-08-18_07:00:00", cwd_full, 24, timeslice_offset_hours, "CCHC1", lake_number, &
653         lookback_seconds, total_counts, observed_counts, time_step_seconds, forecast_discharges, forecast_found, &
654         all_discharges_synthetic, time_series_file_name)
656         print *, "Checking reading of time series file"
657         if ( forecast_discharges(8) .ge. 0.8 - epsilon(0.8) .and. forecast_discharges(8) .le. 0.8 + epsilon(0.8) ) then
658             print *, "========================================================================"
659             print *, "All RFC Forecast Reservoir Time Series Object Tests Passed"
660             print *, "========================================================================"
661             rv = .true.
662         else
663             print *, "========================================================================"
664             print *, "Not All RFC Forecast Reservoir Time Series Object Tests Passed"
665             print *, '8th discharge in the discharge array should be 0.8'
666             print *, "========================================================================"
667         end if
669     end function test_rfc_forecasts_time_series_object
672     ! This tests an RFC Forecasts Reservoir functionality to run levelpool whenever it does
673     ! not find its corresponding time series file within a given lookback window.
674     function test_rfc_forecasts_levelpool_fallback() result(rv)
675         implicit none
676         logical rv, rv1, rv2                        ! test result
677         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
678         real :: outflow, inflow
679         real :: water_elevation
680         real :: prev_time_inflow
681         real :: lake_area, weir_elevation, weir_coefficient
682         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
683         real :: orifice_area, max_depth, initial_fractional_depth
684         integer(kind=int64):: lake_number
685         integer :: reservoir_type
686         integer :: timestep_count
687         character*256 :: cwd_full
688         integer :: dynamic_reservoir_type
689         real :: assimilated_value
690         character(len=256) :: assimilated_source_file
692         prev_time_inflow = 0.0
693         timestep_count = 0
694         water_elevation = 0.0
695         rv = .false.
696         rv2 = .false.
698         lake_area = 2.096320037841796875e+02
699         weir_elevation = 1.332074047851562455e+03
700         weir_coefficient = 4.000000000000000222e-01
701         weir_length = 1.000000000000000000e+01
702         dam_length = 10.0
703         orifice_elevation = 1.314473347981770758e+03
704         orifice_coefficient = 1.000000000000000056e-01
705         orifice_area = 1.0
706         max_depth = 1.335180053710937500e+03
707         initial_fractional_depth = 8.999999761581420898e-01
708         lake_number = 3745478
709         reservoir_type = 4
711         cwd_full = "../../../../tests/local/reservoir_testing_files/"
713         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
714         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
715         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", "2019-05-01_22:00:00", cwd_full, 24)
717         water_elevation = 1331.18005
719         do timestep_count = 1, 120
721             inflow = inflow_array(timestep_count)
722             call rfc_forecasts_reservoir_data%run(inflow, &
723             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
724             assimilated_value, assimilated_source_file)
726             prev_time_inflow = inflow
728             print *, outflow
729         end do
731         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
733         if (outflow .ge. 1.81549609 - epsilon(1.81549609) .and. outflow .le. 1.81549609 + epsilon(1.81549609) &
734         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
735         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
736         assimilated_source_file == "") then
737             rv = .true.
738             print *, "========================================================================"
739             print *, 'RFC Forecasts Levelpool Fallback Test Passed'
740             print *, "========================================================================"
741         else
742             print *, "========================================================================"
743             print *, 'RFC Forecasts Levelpool Fallback Test Failed'
744             print *, 'Outflow should be 1.81549609'
745             print *, 'dynamic_reservoir_type needs to be 1 because reservoir outputs level pool instead of RFC forecast'
746             print *, "========================================================================"
748         end if
750     end function test_rfc_forecasts_levelpool_fallback
753     ! This tests an RFC Forecast Reservoirs functionality to offset 3 hours in the future
754     ! and look back up to 24 hours to find a time series file and output the appropriate
755     ! discharge in the array that matches up with model time as would be done in a standard
756     ! analysis and assimilation run.
757     function test_rfc_forecasts_time_series_output_with_lookback_and_offset() result(rv)
758         implicit none
759         logical rv, rv1, rv2                        ! test result
760         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
761         real :: outflow, inflow
762         real :: water_elevation
763         real :: prev_time_inflow
764         real :: lake_area, weir_elevation, weir_coefficient
765         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
766         real :: orifice_area, max_depth, initial_fractional_depth
767         integer(kind=int64):: lake_number
768         integer :: reservoir_type
769         integer :: timestep_count
770         character*256 :: cwd_full
771         integer :: dynamic_reservoir_type
772         real :: assimilated_value
773         character(len=256) :: assimilated_source_file
775         prev_time_inflow = 0.0
776         timestep_count = 0
777         water_elevation = 0.0
778         rv = .false.
780         lake_area = 2.096320037841796875e+02
781         weir_elevation = 1.332074047851562455e+03
782         weir_coefficient = 4.000000000000000222e-01
783         weir_length = 1.000000000000000000e+01
784         dam_length = 10.0
785         orifice_elevation = 1.314473347981770758e+03
786         orifice_coefficient = 1.000000000000000056e-01
787         orifice_area = 1.0
788         max_depth = 1.335180053710937500e+03
789         initial_fractional_depth = 8.999999761581420898e-01
790         lake_number = 17609317
791         reservoir_type = 4
793         cwd_full = "../../../../tests/local/reservoir_testing_files/"
795         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
796         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
797         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", "2019-08-18_09:00:00", cwd_full, 24)
799         water_elevation = 1331.18005
801         do timestep_count = 1, 80
803             inflow = inflow_array(timestep_count)
804             call rfc_forecasts_reservoir_data%run(inflow, &
805             inflow, 0.0, water_elevation, outflow, 3600.0, dynamic_reservoir_type, &
806             assimilated_value, assimilated_source_file)
808             prev_time_inflow = inflow
810             print *, outflow
811         end do
813         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
815         if (outflow .ge. 3.6 - epsilon(3.6) .and. outflow .le. 3.6 + epsilon(3.6) .and. &
816         water_elevation .ge. 1331.43604 - epsilon(1331.43604) .and. water_elevation .le. 1331.43604 + epsilon(1331.43604) .and. &
817         dynamic_reservoir_type == 4 .and. assimilated_value .ge. 3.6 - epsilon(3.6) &
818         .and. assimilated_value .le. 3.6 + epsilon(3.6) .and. &
819         assimilated_source_file == "2019-08-18_00.60min.CCHC1.RFCTimeSeries.ncdf") then
820             rv = .true.
821             print *, "========================================================================"
822             print *, 'RFC Forecasts Time Series Output Test Passed'
823             print *, "========================================================================"
824         else
825             print *, "========================================================================"
826             print *, 'RFC Forecasts Time Series Output Test Failed'
827             print *, 'Outflow should be 3.6'
828             print *, 'dynamic_reservoir_type needs to be 4 for RFC forecast output'
829             print *, "========================================================================"
831         end if
833     end function test_rfc_forecasts_time_series_output_with_lookback_and_offset
836     ! This tests an RFC Forecast Reservoir's functionality to run levelpool if
837     ! it receives any values from a Time Series Forecast file that are negative.
838     function test_rfc_forecasts_time_series_output_with_negative_data() result(rv)
839         implicit none
840         logical rv, rv1, rv2                        ! test result
841         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
842         real :: outflow, inflow
843         real :: water_elevation
844         real :: prev_time_inflow
845         real :: lake_area, weir_elevation, weir_coefficient
846         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
847         real :: orifice_area, max_depth, initial_fractional_depth
848         integer(kind=int64):: lake_number
849         integer :: reservoir_type
850         integer :: timestep_count
851         character*256 :: cwd_full
852         integer :: dynamic_reservoir_type
853         real :: assimilated_value
854         character(len=256) :: assimilated_source_file
855         prev_time_inflow = 0.0
856         timestep_count = 0
857         water_elevation = 0.0
858         rv = .false.
859         lake_area = 2.096320037841796875e+02
860         weir_elevation = 1.332074047851562455e+03
861         weir_coefficient = 4.000000000000000222e-01
862         weir_length = 1.000000000000000000e+01
863         dam_length = 10.0
864         orifice_elevation = 1.314473347981770758e+03
865         orifice_coefficient = 1.000000000000000056e-01
866         orifice_area = 1.0
867         max_depth = 1.335180053710937500e+03
868         initial_fractional_depth = 8.999999761581420898e-01
869         lake_number = 17609317
870         reservoir_type = 4
871         cwd_full = "../../../../tests/local/reservoir_testing_files/"
872         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
873         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
874         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", "2019-09-18_09:00:00", cwd_full, 24)
875         water_elevation = 1331.18005
876         do timestep_count = 1, 120
877             inflow = inflow_array(timestep_count)
878             call rfc_forecasts_reservoir_data%run(inflow, &
879             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
880             assimilated_value, assimilated_source_file)
881             prev_time_inflow = inflow
882             print *, outflow
883         end do
884         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
885         if (outflow .ge. 1.81549609 - epsilon(1.81549609) .and. outflow .le. 1.81549609 + epsilon(1.81549609) &
886         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
887         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
888         assimilated_source_file == "") then
889             rv = .true.
890             print *, "========================================================================"
891             print *, 'RFC Forecasts Time Series Output Negative Data Test Passed'
892             print *, "========================================================================"
893         else
894             print *, "========================================================================"
895             print *, 'RFC Forecasts Time Series Output Negative Data Test Failed'
896             print *, 'Outflow should be 1.81549609'
897             print *, 'dynamic_reservoir_type needs to be 1 because reservoir outputs level pool instead of RFC forecast'
898             print *, "========================================================================"
899         end if
900     end function test_rfc_forecasts_time_series_output_with_negative_data
901     ! This tests an RFC Forecast Reservoir's functionality to run levelpool if
902     ! it receives all synthetic values from a Time Series Forecast file.
903     function test_rfc_forecasts_time_series_output_all_synthetic_data() result(rv)
904         implicit none
905         logical rv, rv1, rv2                        ! test result
906         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
907         real :: outflow, inflow
908         real :: water_elevation
909         real :: prev_time_inflow
910         real :: lake_area, weir_elevation, weir_coefficient
911         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
912         real :: orifice_area, max_depth, initial_fractional_depth
913         integer(kind=int64):: lake_number
914         integer :: reservoir_type
915         integer :: timestep_count
916         character*256 :: cwd_full
917         integer :: dynamic_reservoir_type
918         real :: assimilated_value
919         character(len=256) :: assimilated_source_file
921         prev_time_inflow = 0.0
922         timestep_count = 0
923         water_elevation = 0.0
924         rv = .false.
926         lake_area = 2.096320037841796875e+02
927         weir_elevation = 1.332074047851562455e+03
928         weir_coefficient = 4.000000000000000222e-01
929         weir_length = 1.000000000000000000e+01
930         dam_length = 10.0
931         orifice_elevation = 1.314473347981770758e+03
932         orifice_coefficient = 1.000000000000000056e-01
933         orifice_area = 1.0
934         max_depth = 1.335180053710937500e+03
935         initial_fractional_depth = 8.999999761581420898e-01
936         lake_number = 17609317
937         reservoir_type = 4
939         cwd_full = "../../../../tests/local/reservoir_testing_files/"
941         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
942         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
943         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", "2019-10-18_09:00:00", cwd_full, 24)
945         water_elevation = 1331.18005
947         do timestep_count = 1, 120
949             inflow = inflow_array(timestep_count)
950             call rfc_forecasts_reservoir_data%run(inflow, &
951             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
952             assimilated_value, assimilated_source_file)
954             prev_time_inflow = inflow
956             print *, outflow
957         end do
959         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
961         if (outflow .ge. 1.81549609 - epsilon(1.81549609) .and. outflow .le. 1.81549609 + epsilon(1.81549609) &
962         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
963         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
964         assimilated_source_file == "") then
965             rv = .true.
966             print *, "========================================================================"
967             print *, 'RFC Forecasts Time Series Output All Synthetic Data Test Passed'
968             print *, "========================================================================"
969         else
970             print *, "========================================================================"
971             print *, 'RFC Forecasts Time Series Output All Synthetic Data Test Failed'
972             print *, 'Outflow should be 1.81549609'
973             print *, 'dynamic_reservoir_type needs to be 1 because reservoir outputs level pool instead of RFC forecast'
974             print *, "========================================================================"
976         end if
978     end function test_rfc_forecasts_time_series_output_all_synthetic_data
981     ! This tests an RFC Forecast Reservoir's functionality to run levelpool if
982     ! it receives any values from a Time Series Forecast file that are greater than the max
983     ! historical MS river flow.
984     function test_rfc_forecasts_data_exceeding_max_range() result(rv)
985         implicit none
986         logical rv, rv1, rv2                        ! test result
987         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
988         real :: outflow, inflow
989         real :: water_elevation
990         real :: prev_time_inflow
991         real :: lake_area, weir_elevation, weir_coefficient
992         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
993         real :: orifice_area, max_depth, initial_fractional_depth
994         integer(kind=int64):: lake_number
995         integer :: reservoir_type
996         integer :: timestep_count
997         character*256 :: cwd_full
998         integer :: dynamic_reservoir_type
999         real :: assimilated_value
1000         character(len=256) :: assimilated_source_file
1001         prev_time_inflow = 0.0
1002         timestep_count = 0
1003         water_elevation = 0.0
1004         rv = .false.
1005         lake_area = 2.096320037841796875e+02
1006         weir_elevation = 1.332074047851562455e+03
1007         weir_coefficient = 4.000000000000000222e-01
1008         weir_length = 1.000000000000000000e+01
1009         dam_length = 10.0
1010         orifice_elevation = 1.314473347981770758e+03
1011         orifice_coefficient = 1.000000000000000056e-01
1012         orifice_area = 1.0
1013         max_depth = 1.335180053710937500e+03
1014         initial_fractional_depth = 8.999999761581420898e-01
1015         lake_number = 17609317
1016         reservoir_type = 4
1017         cwd_full = "../../../../tests/local/reservoir_testing_files/"
1018         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
1019         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
1020         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", "2019-11-18_09:00:00", cwd_full, 24)
1021         water_elevation = 1331.18005
1022         do timestep_count = 1, 120
1023             inflow = inflow_array(timestep_count)
1024             call rfc_forecasts_reservoir_data%run(inflow, &
1025             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
1026             assimilated_value, assimilated_source_file)
1027             prev_time_inflow = inflow
1028             print *, outflow
1029         end do
1030         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
1031         if (outflow .ge. 1.81549609 - epsilon(1.81549609) .and. outflow .le. 1.81549609 + epsilon(1.81549609) &
1032         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
1033         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
1034         assimilated_source_file == "") then
1035             rv = .true.
1036             print *, "========================================================================"
1037             print *, 'RFC Forecasts Time Series Output Data Exceeding Max Test Passed'
1038             print *, "========================================================================"
1039         else
1040             print *, "========================================================================"
1041             print *, 'RFC Forecasts Time Series Output Data Exceeding Max Test Failed'
1042             print *, 'Outflow should be 1.81549609'
1043             print *, 'dynamic_reservoir_type needs to be 1 because reservoir outputs level pool instead of RFC forecast'
1044             print *, "========================================================================"
1045         end if
1046     end function test_rfc_forecasts_data_exceeding_max_range
1047     ! This tests an RFC Forecast Reservoirs functionality to offset 28 hours in the future
1048     ! and look back up to 24 hours to find a time series file and output the appropriate
1049     ! discharge in the array that matches up with model time as would be done in an extended
1050     ! analysis and assimilation run.
1051     function test_rfc_forecasts_with_offset_for_extended_AnA() result(rv)
1052         implicit none
1053         logical rv, rv1, rv2                        ! test result
1054         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
1055         real :: outflow, inflow
1056         real :: water_elevation
1057         real :: prev_time_inflow
1058         real :: lake_area, weir_elevation, weir_coefficient
1059         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
1060         real :: orifice_area, max_depth, initial_fractional_depth
1061         integer(kind=int64):: lake_number
1062         integer :: reservoir_type
1063         integer :: timestep_count
1064         character*256 :: cwd_full
1065         integer :: dynamic_reservoir_type
1066         real :: assimilated_value
1067         character(len=256) :: assimilated_source_file
1068         prev_time_inflow = 0.0
1069         timestep_count = 0
1070         water_elevation = 0.0
1071         rv = .false.
1072         lake_area = 2.096320037841796875e+02
1073         weir_elevation = 1.332074047851562455e+03
1074         weir_coefficient = 4.000000000000000222e-01
1075         weir_length = 1.000000000000000000e+01
1076         dam_length = 10.0
1077         orifice_elevation = 1.314473347981770758e+03
1078         orifice_coefficient = 1.000000000000000056e-01
1079         orifice_area = 1.0
1080         max_depth = 1.335180053710937500e+03
1081         initial_fractional_depth = 8.999999761581420898e-01
1082         lake_number = 17609317
1083         reservoir_type = 4
1084         cwd_full = "../../../../tests/local/reservoir_testing_files/"
1085         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
1086         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
1087         "../../../../tests/local/reservoir_testing_files/reservoir_index_AnA_represents_Extended_MR.nc", &
1088         "2019-12-18_03:00:00", cwd_full, 24)
1089         water_elevation = 1331.18005
1090         do timestep_count = 1, 35
1091             inflow = inflow_array(timestep_count)
1092             call rfc_forecasts_reservoir_data%run(inflow, &
1093             inflow, 0.0, water_elevation, outflow, 3600.0, dynamic_reservoir_type, &
1094             assimilated_value, assimilated_source_file)
1095             prev_time_inflow = inflow
1096             print *, outflow
1097         end do
1098         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
1099         if (outflow .ge. 7.8 - epsilon(7.8) .and. outflow .le. 7.8 + epsilon(7.8) &
1100         .and. dynamic_reservoir_type == 4 .and. assimilated_value .ge. 7.8 - epsilon(7.8) &
1101         .and. assimilated_value .le. 7.8 + epsilon(7.8) .and. &
1102         assimilated_source_file == "2019-12-19_05.60min.CCHC1.RFCTimeSeries.ncdf") then
1103             rv = .true.
1104             print *, "========================================================================"
1105             print *, 'RFC Forecasts Time Series Output For Extended AnA Test Passed'
1106             print *, "========================================================================"
1107         else
1108             print *, "========================================================================"
1109             print *, 'RFC Forecasts Time Series Output For Extended AnA Test Failed'
1110             print *, 'Outflow should be 7.8'
1111             print *, 'dynamic_reservoir_type needs to be 4 for RFC forecast output'
1112             print *, "========================================================================"
1113         end if
1114     end function test_rfc_forecasts_with_offset_for_extended_AnA
1116     ! This tests an Alaska RFC Forecasts Reservoir/Glacier functionality to pass inflow directly to
1117     ! outflow whenever it does not find its corresponding time series file within a given
1118     ! lookback window.
1119     function test_ak_rfc_forecasts_pass_through_fallback() result(rv)
1120         implicit none
1121         logical rv, rv1, rv2                        ! test result
1122         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
1123         real :: outflow, inflow
1124         real :: water_elevation
1125         real :: prev_time_inflow
1126         real :: lake_area, weir_elevation, weir_coefficient
1127         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
1128         real :: orifice_area, max_depth, initial_fractional_depth
1129         integer(kind=int64):: lake_number
1130         integer :: reservoir_type
1131         integer :: timestep_count
1132         character*256 :: cwd_full
1133         integer :: dynamic_reservoir_type
1134         real :: assimilated_value
1135         character(len=256) :: assimilated_source_file
1136         prev_time_inflow = 0.0
1137         timestep_count = 0
1138         water_elevation = 0.0
1139         rv = .false.
1140         rv2 = .false.
1141         lake_area = 0.0
1142         weir_elevation = 0.0
1143         weir_coefficient = 0.0
1144         weir_length = 0.0
1145         dam_length = 0.0
1146         orifice_elevation = 0.0
1147         orifice_coefficient = 0.0
1148         orifice_area = 0.0
1149         max_depth = 0.0
1150         initial_fractional_depth = 0.0
1151         lake_number = 1
1152         reservoir_type = 5
1153         cwd_full = "../../../../tests/local/reservoir_testing_files/"
1154         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
1155         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.0, lake_number, reservoir_type, &
1156         "../../../../tests/local/reservoir_testing_files/reservoir_index_Standard_AnA_APRFC_GDLs.nc", "2019-05-01_22:00:00", cwd_full, 24)
1157         water_elevation = 0.0
1158         do timestep_count = 1, 120
1159             inflow = inflow_array(timestep_count)
1160             call rfc_forecasts_reservoir_data%run(inflow, &
1161             inflow, 0.0, water_elevation, outflow, 900.0, dynamic_reservoir_type, &
1162             assimilated_value, assimilated_source_file)
1163             prev_time_inflow = inflow
1164             print *, outflow
1165         end do
1166         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
1167         if (outflow .ge. inflow - epsilon(inflow) .and. outflow .le. inflow + epsilon(inflow) &
1168         .and. water_elevation .ge. 0.0 - epsilon(0.0) .and. water_elevation .le. 0.0 + epsilon(0.0) &
1169         .and. dynamic_reservoir_type == 1 .and. assimilated_value .ge. -9999.0 - epsilon(-9999.0) &
1170         .and. assimilated_value .le. -9999.0 + epsilon(-9999.0) .and. &
1171         assimilated_source_file == "") then
1172             rv = .true.
1173             print *, "========================================================================"
1174             print *, 'AK RFC Forecasts Pass Through Fallback Test Passed'
1175             print *, "========================================================================"
1176         else
1177             print *, "========================================================================"
1178             print *, 'AK RFC Forecasts Pass Throug Fallback Test Failed'
1179             print *, 'Outflow should be 144.120911'
1180             print *, 'dynamic_reservoir_type needs to be 1 because reservoir outputs level pool instead of RFC forecast'
1181             print *, "========================================================================"
1182         end if
1183     end function test_ak_rfc_forecasts_pass_through_fallback
1185     ! This tests an Alaska RFC Forecast Reservoir/Glacier functionality to offset 3 hours in the future
1186     ! and look back up to 24 hours to find a time series file and output the appropriate
1187     ! discharge in the array that matches up with model time as would be done in a standard
1188     ! analysis and assimilation run.
1189     function test_ak_rfc_forecasts_time_series_with_lookback_and_offset() result(rv)
1190         implicit none
1191         logical rv, rv1, rv2                        ! test result
1192         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
1193         real :: outflow, inflow
1194         real :: water_elevation
1195         real :: prev_time_inflow
1196         real :: lake_area, weir_elevation, weir_coefficient
1197         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
1198         real :: orifice_area, max_depth, initial_fractional_depth
1199         integer(kind=int64):: lake_number
1200         integer :: reservoir_type
1201         integer :: timestep_count
1202         character*256 :: cwd_full
1203         integer :: dynamic_reservoir_type
1204         real :: assimilated_value
1205         character(len=256) :: assimilated_source_file
1206         prev_time_inflow = 0.0
1207         timestep_count = 0
1208         water_elevation = 0.0
1209         rv = .false.
1211         lake_area = 0.0
1212         weir_elevation = 0.0
1213         weir_coefficient = 0.0
1214         weir_length = 0.0
1215         dam_length = 0.0
1216         orifice_elevation = 0.0
1217         orifice_coefficient = 0.0
1218         orifice_area = 0.0
1219         max_depth = 0.0
1220         initial_fractional_depth = 0.0
1221         lake_number = 1
1222         reservoir_type = 5
1223         cwd_full = "../../../../tests/local/reservoir_testing_files/"
1224         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
1225         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.0, lake_number, reservoir_type, &
1226         "../../../../tests/local/reservoir_testing_files/reservoir_index_Standard_AnA_APRFC_GDLs.nc", "2019-08-18_09:00:00", cwd_full, 24)
1227         water_elevation = 0.0
1228         do timestep_count = 1, 80
1229             inflow = inflow_array(timestep_count)
1230             call rfc_forecasts_reservoir_data%run(inflow, &
1231             inflow, 0.0, water_elevation, outflow, 3600.0, dynamic_reservoir_type, &
1232             assimilated_value, assimilated_source_file)
1233             prev_time_inflow = inflow
1234             print *, outflow
1235         end do
1236         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
1237         if (outflow .ge. 180.392075 - epsilon(180.392075) .and. outflow .le. 180.392075 + epsilon(180.392075) &
1238         .and. water_elevation .ge. 0.0 - epsilon(0.0) .and. water_elevation .le. 0.0 + epsilon(0.0) &
1239         .and. dynamic_reservoir_type == 5 .and. assimilated_value .ge. 3.6 - epsilon(3.6) &
1240         .and. assimilated_value .le. 3.6 + epsilon(3.6) .and. &
1241         assimilated_source_file == "2019-08-18_00.60min.SNOA2.RFCTimeSeries.ncdf") then
1242             rv = .true.
1243             print *, "========================================================================"
1244             print *, 'AK RFC Forecasts Time Series Output Test Passed'
1245             print *, "========================================================================"
1246         else
1247             print *, "========================================================================"
1248             print *, 'AK RFC Forecasts Time Series Output Test Failed'
1249             print *, 'Outflow should be 180.392075'
1250             print *, 'dynamic_reservoir_type needs to be 5 for AK RFC forecast output'
1251             print *, "========================================================================"
1252         end if
1253     end function test_ak_rfc_forecasts_time_series_with_lookback_and_offset
1255     ! This tests an RFC Forecast Reservoirs functionality to correct water elevation going
1256     ! over max water elevation
1257     function test_rfc_forecasts_over_max_water_elevation() result(rv)
1258         implicit none
1259         logical rv, rv1, rv2                        ! test result
1260         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
1261         real :: outflow, inflow
1262         real :: water_elevation
1263         real :: prev_time_inflow
1264         real :: lake_area, weir_elevation, weir_coefficient
1265         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
1266         real :: orifice_area, max_depth, initial_fractional_depth
1267         integer :: lake_number, reservoir_type
1268         integer :: timestep_count
1269         character*256 :: cwd_full
1270         integer :: dynamic_reservoir_type
1271         real :: assimilated_value
1272         character(len=256) :: assimilated_source_file
1274         prev_time_inflow = 0.0
1275         timestep_count = 0
1276         water_elevation = 0.0
1277         rv = .false.
1279         lake_area = 2.096320037841796875e+02
1280         weir_elevation = 1.332074047851562455e+03
1281         weir_coefficient = 4.000000000000000222e-01
1282         weir_length = 1.000000000000000000e+01
1283         dam_length = 10.0
1284         orifice_elevation = 1.314473347981770758e+03
1285         orifice_coefficient = 1.000000000000000056e-01
1286         orifice_area = 1.0
1287         max_depth = 1.335180053710937500e+03
1288         initial_fractional_depth = 8.999999761581420898e-01
1289         lake_number = 17609317
1290         reservoir_type = 4
1292         cwd_full = "../../../../tests/local/reservoir_testing_files/"
1294         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
1295         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
1296         "../../../../tests/local/reservoir_testing_files/reservoir_index_short_range.nc", "2019-08-18_09:00:00", cwd_full, 24)
1298         !water_elevation = 1331.18005
1300         water_elevation = 1.335180053710937500e+03
1302         do timestep_count = 1, 80
1304             inflow = inflow_array(timestep_count)
1305             call rfc_forecasts_reservoir_data%run(inflow, &
1306             inflow, 0.0, water_elevation, outflow, 3600.0, dynamic_reservoir_type, &
1307             assimilated_value, assimilated_source_file)
1309             prev_time_inflow = inflow
1311             print *, outflow
1312         end do
1314         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
1316         if (outflow .ge. 3.6 - epsilon(3.6) .and. outflow .le. 3.6 + epsilon(3.6) .and. &
1317         water_elevation .ge. 1335.18005 - epsilon(1335.18005) .and. water_elevation .le. 1335.18005 + epsilon(1335.18005) .and. &
1318         dynamic_reservoir_type == 4 .and. assimilated_value .ge. 3.6 - epsilon(3.6) &
1319         .and. assimilated_value .le. 3.6 + epsilon(3.6) .and. &
1320         assimilated_source_file == "2019-08-18_00.60min.CCHC1.RFCTimeSeries.ncdf") then
1321             rv = .true.
1322             print *, "========================================================================"
1323             print *, 'RFC Forecasts Over Max Water Elevation Test Passed'
1324             print *, "========================================================================"
1325         else
1326             print *, "========================================================================"
1327             print *, 'RFC Forecasts Over Max Water Elevation Test Failed'
1328             print *, 'Outflow should be 3.6'
1329             print *, 'Water Elevation should be 1335.18005'
1330             print *, 'dynamic_reservoir_type needs to be 4 for RFC forecast output'
1331             print *, "========================================================================"
1333         end if
1335     end function test_rfc_forecasts_over_max_water_elevation
1338     ! This tests an RFC Forecast Reservoirs functionality to offset 12 hours in the future
1339     ! and look back up to 24 hours to find a time series file and output the appropriate
1340     ! discharge in the array that matches up with model time as would be done in a long range
1341     ! analysis and assimilation run.
1342     function test_rfc_forecasts_with_offset_for_long_range_AnA() result(rv)
1343         implicit none
1344         logical rv, rv1, rv2                        ! test result
1345         type (rfc_forecasts) :: rfc_forecasts_reservoir_data
1346         real :: outflow, inflow
1347         real :: water_elevation
1348         real :: prev_time_inflow
1349         real :: lake_area, weir_elevation, weir_coefficient
1350         real :: weir_length, dam_length, orifice_elevation, orifice_coefficient
1351         real :: orifice_area, max_depth, initial_fractional_depth
1352         integer :: lake_number, reservoir_type
1353         integer :: timestep_count
1354         character*256 :: cwd_full
1355         integer :: dynamic_reservoir_type
1356         real :: assimilated_value
1357         character(len=256) :: assimilated_source_file
1359         prev_time_inflow = 0.0
1360         timestep_count = 0
1361         water_elevation = 0.0
1362         rv = .false.
1364         lake_area = 2.096320037841796875e+02
1365         weir_elevation = 1.332074047851562455e+03
1366         weir_coefficient = 4.000000000000000222e-01
1367         weir_length = 1.000000000000000000e+01
1368         dam_length = 10.0
1369         orifice_elevation = 1.314473347981770758e+03
1370         orifice_coefficient = 1.000000000000000056e-01
1371         orifice_area = 1.0
1372         max_depth = 1.335180053710937500e+03
1373         initial_fractional_depth = 8.999999761581420898e-01
1374         lake_number = 17609317
1375         reservoir_type = 4
1377         cwd_full = "../../../../tests/local/reservoir_testing_files/"
1379         call rfc_forecasts_reservoir_data%init(water_elevation, lake_area, weir_elevation, weir_coefficient, &
1380         weir_length, dam_length, orifice_elevation, orifice_coefficient, orifice_area, max_depth, 0.9, lake_number, reservoir_type, &
1381         "../../../../tests/local/reservoir_testing_files/reservoir_index_Long_Range_AnA.nc", &
1382         "2019-12-18_19:00:00", cwd_full, 24)
1384         water_elevation = 1331.18005
1386         do timestep_count = 1, 19
1388             inflow = inflow_array(timestep_count)
1389             call rfc_forecasts_reservoir_data%run(inflow, &
1390             inflow, 0.0, water_elevation, outflow, 3600.0, dynamic_reservoir_type, &
1391             assimilated_value, assimilated_source_file)
1393             prev_time_inflow = inflow
1395             print *, outflow
1396         end do
1398         print *, 'dynamic_reservoir_type: ', dynamic_reservoir_type
1400         if (outflow .ge. 7.8 - epsilon(7.8) .and. outflow .le. 7.8 + epsilon(7.8) &
1401         .and. dynamic_reservoir_type == 4 .and. assimilated_value .ge. 7.8 - epsilon(7.8) &
1402         .and. assimilated_value .le. 7.8 + epsilon(7.8) .and. &
1403         assimilated_source_file == "2019-12-19_05.60min.CCHC1.RFCTimeSeries.ncdf") then
1404             rv = .true.
1405             print *, "========================================================================"
1406             print *, 'RFC Forecasts Time Series Output For Long Range AnA Test Passed'
1407             print *, "========================================================================"
1408         else
1409             print *, "========================================================================"
1410             print *, 'RFC Forecasts Time Series Output For Long Range AnA Test Failed'
1411             print *, 'Outflow should be 7.8'
1412             print *, 'dynamic_reservoir_type needs to be 4 for RFC forecast output'
1413             print *, "========================================================================"
1415         end if
1417     end function test_rfc_forecasts_with_offset_for_long_range_AnA
1419 end program