-
Notifications
You must be signed in to change notification settings - Fork 3
/
flow_sim.f90
executable file
·382 lines (297 loc) · 11.9 KB
/
flow_sim.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
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
!***************************************************************
! Copyright (c) 2017 Battelle Memorial Institute
! Licensed under modified BSD License. A copy of this license can be
! found in the LICENSE file in the top level directory of this
! distribution.
!***************************************************************
!
! NAME: flow_sim
!
! VERSION and DATE: MASS1 v0.75 3/25/98
!
! PURPOSE: hydraulics solver for St. Venant Equations
!
! RETURNS:
!
! REQUIRED:
!
! LOCAL VARIABLES:
!
! COMMENTS: based on an old homework problem from Holly.
!
!
! MOD HISTORY: added hydrad to CALL section and added
! calculation of froude_num, friction_slope, etc. ; mcr 11/21/1997
! added lateral inflow; mcr 3/25/98
!
!***************************************************************
!
! ----------------------------------------------------------------
! SUBROUTINE depth_check
! ----------------------------------------------------------------
SUBROUTINE depth_check(thalweg, y, q)
USE general_vars, ONLY: depth_minimum
IMPLICIT NONE
DOUBLE PRECISION, INTENT(IN) :: thalweg
DOUBLE PRECISION, INTENT(INOUT) :: y, q
DOUBLE PRECISION :: depth
depth = y - thalweg
IF (depth .LT. depth_minimum) THEN
y = thalweg + depth_minimum
! q = 0.0
END if
END SUBROUTINE depth_check
SUBROUTINE flow_sim
! $DEBUG
USE general_vars
USE link_vars
USE bctable
USE point_vars
USE fluvial_coeffs
USE flow_coeffs
USE logicals , ONLY : do_latflow
IMPLICIT NONE
DOUBLE PRECISION :: a,b,c,d,g,ap,bp,cp,dp,gp,denom
DOUBLE PRECISION :: latq_old,latq_new
DOUBLE PRECISION :: depth,area_temp,width,conveyance,dkdy,hydrad
DOUBLE PRECISION :: sum,sum2,bcval,dy,dq,y_new_time, q_new_time
DOUBLE PRECISION :: temp
DOUBLE PRECISION :: delta_x
INTEGER :: i,j,point_num,link,point,table_type
LOGICAL :: fluvial
CHARACTER (LEN=1024) :: msg
! run through links top down according to computational order
links_forward: DO i=1,maxlinks
link = comporder(i)
SELECT CASE(linktype(link))
CASE(1,20,21)
fluvial = .TRUE.
CASE(2,3,4,5,6,7,12,13)
fluvial = .FALSE.
END SELECT
! set upstream bc q(t) or junction condition for first point
point = 1
IF(num_con_links(link) == 0)THEN ! must be an upstream most link
SELECT CASE(linktype(link))
CASE(1,20)
call bc_table_interpolate(linkbc, linkbc_table(link), time/time_mult)
bcval = bc_table_current(linkbc, linkbc_table(link), 1)
CASE(21)
call bc_table_interpolate(hydrobc, linkbc_table(link), time/time_mult)
temp = bc_table_current(hydrobc, linkbc_table(link), 1)
bcval = bc_table_current(hydrobc, linkbc_table(link), 2)
bcval = bcval + temp ! total flow rate at the dam
END SELECT
q1 = q(link,point)
e(link,point) = 0.0
f(link,point) = bcval - q1
ELSE
sum = 0.0
sum2 = 0.0
DO j=1,num_con_links(link)
sum = sum + e(con_links(link,j),maxpoints(con_links(link,j)))
sum2 = sum2 + q(con_links(link,j),maxpoints(con_links(link,j))) + &
&f(con_links(link,j),maxpoints(con_links(link,j))) + &
&e(con_links(link,j),maxpoints(con_links(link,j)))*&
&(y(link,point) - y(con_links(link,j),maxpoints(con_links(link,j))))
END DO
e(link,point) = sum
f(link,point) = -q(link,point) + sum2
END IF
! do internal links
points: DO point=1,maxpoints(link)-1
IF( fluvial )THEN
! fluvial links
! set geometric data for points i, i+1
point_num = point
depth = y(link,point) - thalweg(link,point) !remember y is ELEVATION
depth = MAX(depth, depth_minimum)
CALL section(link,point_num,depth,area_temp,width,conveyance,dkdy,hydrad)
!CALL section(link,point_num,depth,area_temp,width,conveyance,dkdy)
d1 = depth
fr1 = froude_num(link,point)
y1 = y(link,point)
q1 = q(link,point)
a1 = area_temp
b1 = width
k1 = conveyance
ky1 = dkdy
IF (a1 .GT. 0.0D00) THEN
vel(link,point_num) = q1/a1
area_old(link,point_num) = a1
ELSE
vel(link,point_num) = 0.0
area_old(link,point_num) = 0.0
END IF
q_old(link,point_num) = q1
y_old(link,point_num) = y1
point_num = point + 1
depth = y(link,point+1) - thalweg(link,point+1)
depth = MAX(depth, depth_minimum)
CALL section(link,point_num,depth,area_temp,width,conveyance,dkdy,hydrad)
!CALL section(link,point_num,depth,area_temp,width,conveyance,dkdy)
d2 = depth
fr2 = froude_num(link,point+1)
y2 = y(link,point+1)
q2 = q(link,point+1)
a2 = area_temp
b2 = width
k2 = conveyance
ky2 = dkdy
IF (a2 .GT. 0.0) THEN
vel(link,point_num) = q2/a2
area_old(link,point_num) = a2
ELSE
vel(link,point_num) = 0.0
area_old(link,point_num) = 0.0
END IF
q_old(link,point_num) = q2
y_old(link,point_num) = y2
delta_x = ABS(x(link,point+1) - x(link,point))
! uniform lateral inflow per unit length
IF(do_latflow)THEN
IF(latflowbc_table(link) /= 0)THEN
latq_old = lateral_inflow(link,point)
lateral_inflow_old(link,point) = latq_old
call bc_table_interpolate(latflowbc, latflowbc_table(link), time/time_mult)
lateral_inflow(link,point) = bc_table_current(latflowbc, latflowbc_table(link), 1)
latq_new = lateral_inflow(link,point)
ELSE
latq_old = 0.0
latq_new = 0.0
ENDIF
ELSE
latq_old = 0.0
latq_new = 0.0
ENDIF
CALL fluvial_coeff(a,b,c,d,g,ap,bp,cp,dp,gp,delta_x,delta_t,grav,&
&latq_old,latq_new,lpiexp(link))
! nonfluvial internal links ----------------------------
ELSE
! nonfluvial links also need q_old for
! transport
IF(linktype(link) == 6)THEN ! hydropower plant
call bc_table_interpolate(hydrobc, linkbc_table(link), time/time_mult)
temp = bc_table_current(hydrobc, linkbc_table(link), 1)
bcval = bc_table_current(hydrobc, linkbc_table(link), 2)
bcval = bcval + temp ! total flow rate at the dam
ELSE ! other non-fluvial links
call bc_table_interpolate(linkbc, linkbc_table(link), time/time_mult)
bcval = bc_table_current(linkbc, linkbc_table(link), 1)
ENDIF
CALL nonfluvial_coeff(link,point,bcval,a,b,c,d,g,ap,bp,cp,dp,gp)
q_old(link, :) = q(link, :)
END IF
denom = (c*dp - cp*d)
l(link,point) = (a*dp - ap*d)/denom
m(link,point) = (b*dp - bp*d)/denom
n(link,point) = (d*gp - dp*g)/denom
denom = b - m(link,point)*(c + d*e(link,point))
e(link,point+1) = (l(link,point)*(c + d*e(link,point)) - a)/denom
f(link,point+1) = (n(link,point)*(c + d*e(link,point)) +d*f(link,point) + g)/denom
END DO points
END DO links_forward
!------------------------------------------------------------------------------
! run through links bottom to top
links_backward: DO i=maxlinks,1,-1
link = comporder(i)
point = maxpoints(link)
! set downstream bc y(t) or Q(t) OR junction conditions
IF(ds_conlink(link) == 0)THEN
SELECT CASE(dsbc_type)
CASE(1)
! given downstream stage y(t)
y1 = y(link,point)
call bc_table_interpolate(linkbc, dsbc_table(link), time/time_mult)
y_new_time = bc_table_current(linkbc, dsbc_table(link), 1)
dy = y_new_time - y1
dq = e(link,point)*dy + f(link,point)
CASE(2)
! given Q(t)
q1 = q(link,point)
call bc_table_interpolate(linkbc, dsbc_table(link), time/time_mult)
q_new_time = bc_table_current(linkbc, dsbc_table(link), 1)
dq = q_new_time - q1
dy = (dq - f(link,point))/e(link,point)
END SELECT
!update stage and discharge at last point on the link
y(link,point) = y(link,point) + dy
q(link,point) = q(link,point) + dq
ELSE
! junction conditions
dy = y(ds_conlink(link),1) - y(link,point)
dq = e(link,point)*dy + f(link,point)
y(link,point) = y(link,point) + dy
q(link,point) = q(link,point) + dq
END IF
DO point=maxpoints(link)-1,1,-1
IF(linktype(link) == 2)THEN
call bc_table_interpolate(linkbc, linkbc_table(link), time/time_mult)
bcval = bc_table_current(linkbc, linkbc_table(link), 1)
dq = bcval - q(link,point)
dy = (dq - f(link,point))/e(link,point)
ELSEIF(linktype(link) == 6)THEN
call bc_table_interpolate(hydrobc, linkbc_table(link), time/time_mult)
temp = bc_table_current(hydrobc, linkbc_table(link), 1)
bcval = bc_table_current(hydrobc, linkbc_table(link), 2)
bcval = bcval + temp ! total flow rate at the dam
dq = bcval - q(link,point)
dy = (dq - f(link,point))/e(link,point)
ELSE
dy = l(link,point)*dy + m(link,point)*dq + n(link,point)
dq = e(link,point)*dy + f(link,point)
ENDIF
y(link,point) = y(link,point) + dy
q(link,point) = q(link,point) + dq
END DO
END DO links_backward
! Initialize these so unused parts of the array do not affect later
! calls to max/min
froude_num = 10.0
courant_num = 0.0
diffuse_num = 0.0
!------------------------------------------------------------------------------
! computes additional data after hydraulics have been
! updated for this time
!-------------------------------------------------------------------------------
DO link = 1,maxlinks
SELECT CASE(linktype(link))
CASE(1,20,21)
DO point = 1,maxpoints(link)
CALL depth_check(thalweg(link, point), y(link,point), q(link,point))
depth = y(link,point) - thalweg(link,point)
CALL section(link,point,depth,area_temp,width,conveyance,dkdy,hydrad)
!CALL section(link,point,depth,area_temp,width,conveyance,dkdy)
IF (point .GE. maxpoints(link)) THEN
delta_x = ABS(x(link,point-1) - x(link,point))
ELSE
delta_x = ABS(x(link,point+1) - x(link,point))
END IF
top_width(link,point) = width
hyd_radius(link,point) = hydrad
IF (area_temp .GT. 0.0) THEN
area(link,point) = area_temp
froude_num(link,point) = &
&SQRT((q(link,point)**2*width)/(grav*area_temp**3))
friction_slope(link,point) =&
& ((q(link,point)*manning(link,point))/&
& (res_coeff*area_temp*(hydrad**2.0)**0.3333333))**2.0
courant_num(link, point) = ABS(q(link,point))/area_temp*delta_t/delta_x
ELSE
area(link,point) = 0.0
froude_num(link,point) = 0.0
friction_slope(link,point) = 0.0
courant_num(link, point) = 0.0
END IF
bed_shear(link,point) = &
&unit_weight_h2o*hydrad*friction_slope(link,point)
diffuse_num(link, point) = 2.0*k_diff(link,point)*delta_t/delta_x/delta_x
IF (froude_num(link, point) .GE. 1.0) THEN
WRITE (msg, '("warning: supercritial (Fr=", F5.1, ") indicated at link ", I3, ", point ", I3)')&
&froude_num(link, point), link, point
CALL status_message(msg)
END IF
END DO
END SELECT
END DO
END SUBROUTINE flow_sim