-
Notifications
You must be signed in to change notification settings - Fork 0
/
ch_rtmusk.f90
344 lines (302 loc) · 14.1 KB
/
ch_rtmusk.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
subroutine ch_rtmusk
!! ~ ~ ~ PURPOSE ~ ~ ~
!! this subroutine routes a daily flow through a reach using the
!! Muskingum method
!! code provided by Dr. Valentina Krysanova, Pottsdam Institute for
!! Climate Impact Research, Germany
!! Modified by Balaji Narasimhan
!! Spatial Sciences Laboratory, Texas A&M University
use basin_module
use channel_data_module
use channel_module
use hydrograph_module !, only : ob, icmd, jrch, isdch, fp_stor, ch_stor, wet
use time_module
use channel_velocity_module
use sd_channel_module
use climate_module
use reservoir_module
use reservoir_data_module
use water_body_module
use hru_module, only : hru
use conditional_module
implicit none
integer :: ii !none |current day of simulation
integer :: ihru
integer :: iihru
integer :: icha
integer :: irtstep
integer :: isubstep
integer :: ires
integer :: ihyd
integer :: irel
real :: ch_stor_init !m3 |storage in channel at beginning of day
real :: fp_stor_init !m3 |storage in flood plain above wetlands emergency spillway at beginning of day
real :: wet_stor_init !m3 |storage in flood plain wetlands at beginning of day
real :: inout !m3 |inflow - outflow for day
real :: del_stor !m3 |change in storage of channel + flood plain + wetlands
real :: c1 !units |description
real :: c2 !units |description
real :: c3 !units |description
real :: c4 !m^3 H2O |
real :: p !m |wetted perimeter
real :: c !none |inverse of channel side slope
real :: rh !m |hydraulic radius
real :: topw !m |top width of main channel
real :: qinday !units |description
real :: qoutday !units |description
real :: volrt !units |description
real :: maxrt !units |description
real :: adddep !units |description
real :: addp !units |description
real :: addarea !units |description
real :: rttlc1 !units |description
real :: rttlc2 !units |description
real :: rtevp1 !units |description
real :: rtevp2 !units |description
real :: qman !m^3/s or m/s |flow rate or flow velocity
real :: vc !m/s |flow velocity in reach
real :: aaa !units |description
real :: inflo !m^3 |inflow water volume
real :: inflo_rate !m^3/s |inflow rate
real :: xs_area !m^2 |cross section area of channel
real :: dep_flo = 0. !m |depth of flow
real :: wet_perim !m |wetted perimeter
real :: ttime !hr |travel time through the reach
real :: t_inc !hr |time in routing step - 1/time%step
real :: outflo !m^3 |outflow water volume
real :: tl !m^3 |transmission losses during time step
real :: trans_loss = 0. !m^3 |transmission losses during day
real :: rate_flo !m^3/s |flow rate
real :: ev !m^3 |evaporation during time step
real :: evap = 0. !m^3 |evaporation losses during day
real :: above_prin_fr
real :: ch_out_fr
real :: fp_out_fr
real :: ch_fp_fr
real :: fp_ch_fr
real :: rto
real :: rto1
real :: rto_w
real :: rto_emer
real :: outflo_rate
real :: ch_st !m^3 |water storage in and above channel
real :: fp_st !m^3 |water storage in flood plain
real :: wet_st !m^3 |wetland storage above principal
real :: dts !seconds |time step interval for substep
real :: dthr
real :: scoef
real :: vol_ch
real :: vol_fp_av
real :: vol_tot_av
real :: sum_inflo, sum_outflo
real :: dep
jrch = isdch
jhyd = sd_dat(jrch)%hyd
qinday = 0
qoutday = 0
ob(icmd)%hd(1) = hz
ob(icmd)%hyd_flo = 0.
hyd_rad = 0.
trav_time = 0.
flo_dep = 0.
trans_loss = 0.
evap = 0.
sum_inflo = sum (ob(icmd)%tsin)
!! total wetland volume at start of day
wet_stor(jrch) = hz
do ihru = 1, sd_ch(jrch)%fp%hru_tot
wet_stor(jrch) = wet_stor(jrch) + wet(ihru)
end do
wet_stor_init = wet_stor(jrch)%flo
ch_stor_init = ch_stor(jrch)%flo
fp_stor_init = fp_stor(jrch)%flo
irtstep = 1
isubstep = 0
dts = time%dtm / sd_ch(jrch)%msk%substeps * 60.
dthr = dts / 3600.
if (jrch == 3) ht1%flo = 480000.
!! subdaily time step
do ii = 1, sd_ch(jrch)%msk%nsteps
!! water entering reach during time step - substeps for stability
isubstep = isubstep + 1
if (isubstep > sd_ch(jrch)%msk%substeps) then
irtstep = irtstep + 1
isubstep = 1
end if
inflo = ob(icmd)%tsin(irtstep) / sd_ch(jrch)%msk%substeps
if (jrch == 3) inflo = 10000.
inflo_rate = inflo / dts
!! interpolate rating curve using inflow rates
icha = jrch
call rcurv_interp_flo (icha, inflo_rate)
ch_rcurv(jrch)%in2 = rcurv
!! save variables at each routing time step for sediment routing
if (isubstep == 1 .and. rcurv%wet_perim > 1.e-6) then
hyd_rad(irtstep) = rcurv%xsec_area / rcurv%wet_perim
trav_time(irtstep) = rcurv%ttime
flo_dep(irtstep) = rcurv%dep
end if
!! add inflow to total storage
if (ht1%flo > 1.e-6) then
rto = inflo / ht1%flo
tot_stor(jrch) = tot_stor(jrch) + rto * ht1
end if
!! partition channel and flood plain based on bankfull volume
if (tot_stor(jrch)%flo > ch_rcurv(jrch)%elev(2)%vol_ch) then
!! fill channel to bank full if below
rto = (ch_rcurv(jrch)%elev(2)%vol_ch - ch_stor(jrch)%flo) / tot_stor(jrch)%flo
vol_ch = ch_stor(jrch)%flo
if (rto > 0.) then
ch_stor(jrch) = ch_stor(jrch) + rto * tot_stor(jrch)
end if
!! if flood plain link - adjust flood plain and wetland storages for inflow
if (bsn_cc%i_fpwet == 2) then
!! add flood plain inflow to wetland
do ihru = 1, sd_ch(jrch)%fp%hru_tot
iihru = sd_ch(jrch)%fp%hru(ihru)
!! distribute water by hru fraction of the flood plain
rto1 = (inflo - (ch_rcurv(jrch)%elev(2)%vol_ch - vol_ch)) / tot_stor(jrch)%flo
rto_w = rto1 * sd_ch(jrch)%fp%hru_fr(ihru)
wet(iihru) = wet(iihru) + rto_w * tot_stor(jrch)
!! if above emergency - move back to flood plain storage
rto_emer = (wet(iihru)%flo - wet_ob(iihru)%evol) / wet(iihru)%flo
if (rto_emer > 0.) then
fp_stor(jrch) = fp_stor(jrch) + rto_emer * wet(iihru)
wet(iihru) = wet(iihru) - rto_emer * wet(iihru)
end if
end do
else
!! if no flood plain link - add rest to flood plain
fp_stor(jrch) = fp_stor(jrch) + (1. - rto) * tot_stor(jrch)
end if
else
!! total volume below bankfull
ch_stor(jrch) = tot_stor(jrch)
fp_stor(jrch) = hz
end if ! flow rate above bank full
tot_stor(jrch) = ch_stor(jrch) + fp_stor(jrch)
!! if no water in channel - skip routing and set rating curves to zero
if (tot_stor(jrch)%flo < 1.e-6) then
ch_rcurv(jrch)%in1 = rcz
ch_rcurv(jrch)%out1 = rcz
sd_ch(jrch)%in1_vol = 0.
sd_ch(jrch)%out1_vol = 0.
else
!! Muskingum flood routing method
outflo = sd_ch(jrch)%msk%c1 * inflo + sd_ch(jrch)%msk%c2 * sd_ch(jrch)%in1_vol + &
sd_ch(jrch)%msk%c3 * sd_ch(jrch)%out1_vol
!! save inflow/outflow volumes for next time step (and day) for Muskingum
sd_ch(jrch)%in1_vol = inflo
sd_ch(jrch)%out1_vol = outflo
!! Variable Storage Coefficent method - sc=2*dt/(2*ttime+dt) - ttime=(in2+out1)/2
scoef = 2. * dthr / (ch_rcurv(jrch)%in2%ttime + ch_rcurv(jrch)%out1%ttime + dthr)
scoef = Min (scoef, 1.)
!outflo = scoef * tot_stor(jrch)%flo
outflo = Min (outflo, tot_stor(jrch)%flo)
outflo = Max (outflo, 0.)
!! compute outflow rating curve for next time step
outflo_rate = outflo / dts !convert to cms
call rcurv_interp_flo (jrch, outflo_rate)
ch_rcurv(jrch)%out2 = rcurv
!! add outflow to daily hydrograph and subdaily flow
rto = outflo / tot_stor(jrch)%flo
rto = amin1 (1., rto)
ob(icmd)%hd(1) = ob(icmd)%hd(1) + rto * tot_stor(jrch)
ob(icmd)%hyd_flo(1,irtstep) = ob(icmd)%hyd_flo(1,irtstep) + outflo
!! subtract outflow from total storage
rto = outflo / tot_stor(jrch)%flo
tot_stor(jrch) = (1. - rto) * tot_stor(jrch)
!! readjust channel and flood plain volumes after outflow
if (tot_stor(jrch)%flo <= ch_rcurv(jrch)%elev(2)%vol_ch) then
ch_stor(jrch) = tot_stor(jrch)
fp_stor(jrch) = hz
else
rto = ch_rcurv(jrch)%elev(2)%vol_ch / tot_stor(jrch)%flo
ch_stor(jrch) = rto * tot_stor(jrch)
rto = 1. - rto
fp_stor(jrch) = rto * tot_stor(jrch)
end if
!! if flood plain link - adjust wetland storages for outflow
if (bsn_cc%i_fpwet == 2) then
!! compute release from flood plain wetlands
do ihru = 1, sd_ch(jrch)%fp%hru_tot
iihru = sd_ch(jrch)%fp%hru(ihru)
ires= hru(iihru)%dbs%surf_stor
ihyd = wet_dat(ires)%hyd
irel = wet_dat(ires)%release
!! calc release from decision table
d_tbl => dtbl_res(irel)
wbody => wet(iihru)
wbody_wb => wet_wat_d(iihru)
if (wet_ob(ihru)%area_ha > 1.e-6) then
dep = wbody%flo / wet_ob(ihru)%area_ha / 10000. !m = m3 / ha / 10000m2/ha
else
dep = 0.
end if
call conditions (iihru, irel)
call res_hydro (iihru, irel, ihyd, wet_ob(ihru)%pvol, wet_ob(ihru)%evol, dep, wet_ob(ihru)%weir_hgt)
!! subtract outflow from wetland and add to flood plain storage
wet(iihru) = wet(iihru) - ht2
fp_stor(jrch) = fp_stor(jrch) + ht2
tot_stor(jrch) = tot_stor(jrch) + ht2
end do
end if ! bsn_cc%i_fpwet == 2
!! set rating curve for next time step
ch_rcurv(jrch)%in1 = ch_rcurv(jrch)%in2
ch_rcurv(jrch)%out1 = ch_rcurv(jrch)%out2
!! calculate transmission losses
if (tot_stor(jrch)%flo > 1.e-6) then
ttime = Min(24., rcurv%ttime)
tl = sd_ch(jrch)%chk * sd_ch(jrch)%chl * rcurv%wet_perim * ttime !mm/hr * km * mm * hr = m3
tl = Min(tl, tot_stor(jrch)%flo)
trans_loss = trans_loss + tl
!! subtract evap and transmission loses from channel and flood plain storage
if (tot_stor(jrch)%flo > 1.e-6) then
rto = tl / tot_stor(jrch)%flo
rto = amin1 (1., rto)
ch_stor(jrch) = (1. - rto) * ch_stor(jrch)
fp_stor(jrch) = (1. - rto) * fp_stor(jrch)
tot_stor(jrch) = (1. - rto) * tot_stor(jrch)
else
ch_stor(jrch) = hz
fp_stor(jrch) = hz
tot_stor(jrch) = hz
end if
end if
!! calculate evaporation
if (tot_stor(jrch)%flo > 1.e-6) then
!! calculate width of channel at water level - flood plain evap calculated in wetlands
if (dep_flo <= sd_ch(jrch)%chd) then
topw = ch_rcurv(jrch)%out2%surf_area
else
topw = 1000. * sd_ch(jrch)%chl * sd_ch(jrch)%chw
end if
iwst = ob(icmd)%wst
!! mm/day * m2 / (1000. * sd_ch(jrch)%msk%nsteps)
ev = bsn_prm%evrch * wst(iwst)%weat%pet * topw / (1000. * sd_ch(jrch)%msk%nsteps)
if (ev < 0.) ev = 0.
ev = Min(ev, tot_stor(jrch)%flo)
evap = evap + ev
!! subtract evap and transmission loses from channel and flood plain storage
if (tot_stor(jrch)%flo > 1.e-6) then
rto = ev / tot_stor(jrch)%flo
rto = amin1 (1., rto)
ch_stor(jrch) = (1. - rto) * ch_stor(jrch)
fp_stor(jrch) = (1. - rto) * fp_stor(jrch)
tot_stor(jrch) = (1. - rto) * tot_stor(jrch)
else
ch_stor(jrch) = hz
fp_stor(jrch) = hz
tot_stor(jrch) = hz
end if
end if
tot_stor(jrch) = ch_stor(jrch) + fp_stor(jrch)
end if !! tot_stor(jrch)%flo < 1.e-6 loop
end do !! end of sub-daily loop
!! check water balance at end of day
sum_outflo = sum (ob(icmd)%hyd_flo(1,:))
inout = sum_inflo - sum_outflo - trans_loss - evap
del_stor = (ch_stor(jrch)%flo - ch_stor_init) + (fp_stor(jrch)%flo - fp_stor_init) + &
(wet_stor(jrch)%flo - wet_stor_init)
return
end subroutine ch_rtmusk