-
Notifications
You must be signed in to change notification settings - Fork 3
/
HMC_Module_Phys_HydraulicStructure.f90
479 lines (372 loc) · 25.6 KB
/
HMC_Module_Phys_HydraulicStructure.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
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
!------------------------------------------------------------------------------------------
! File: HMC_Module_Phys_HydraulicStructure.f90
! Author: fabio
!
! Created on February 24, 2016, 2:11 PM
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Module Header
module HMC_Module_Phys_HydraulicStructure
!------------------------------------------------------------------------------------------
! External module(s) for all subroutine in this module
use HMC_Module_Namelist, only: oHMC_Namelist
use HMC_Module_Vars_Loader, only: oHMC_Vars
use HMC_Module_Tools_Debug
! Implicit none for all subroutines in this module
implicit none
!------------------------------------------------------------------------------------------
contains
!------------------------------------------------------------------------------------------
! Subroutine for calculating dam volume spilling
subroutine HMC_Phys_Dam_Spilling(iID, iNDam, dDt, &
a1dVarVDam, a1dVarQoutDam, a1dVarCoeffDam, a1dVarLDam)
!------------------------------------------------------------------------------------------
! Variable(s)
integer(kind = 4) :: iID
integer(kind = 4) :: iI, iJ, iD, iK, iRank !aggiunte variabili
integer(kind = 4) :: iNDam
integer(kind = 4) :: iFlagD, iMaxC !aggiunte variabili
real(kind = 4) :: dDt, dQtmp, dQsprec !aggiunte variabili
real(kind = 4), dimension (iNDam) :: a1dVarVDam
real(kind = 4), dimension (iNDam) :: a1dVarQoutDam, a1dVarCoeffDam
real(kind = 4), dimension (iNDam) :: a1dVarHDam, a1dVarLDam !aggiunte variabili
character(len = 256) :: sQtmp, sD, sHDam, sHMaxDam !aggiunte variabili
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Initialize variable(s)
a1dVarQoutDam = 0.0;
dQtmp = 0.0;
iFlagD = 0;
sQtmp = ''; sHDam = ''; sHMaxDam = '';
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Re-initialize dam height (part added to compute Q of spilling WITHIN the convolution in [m^3/s])
a1dVarHDam = 0.0
! ! Update dam volume using observation(s)
! where (oHMC_Vars(iID)%a1dVDamObs .gt. -9999.0)
! a1dVarVDam = oHMC_Vars(iID)%a1dVDamObs
! endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Check dam availability
if (iNDam .gt. 0) then
!------------------------------------------------------------------------------------------
! Cycle on dam
do iD = 1, iNDam
!------------------------------------------------------------------------------------------
! Check limit(s)
if (a1dVarVDam(iD) .gt. oHMC_Namelist(iID)%dTV*oHMC_Vars(iID)%a1dVMaxDam(iD)) then
!------------------------------------------------------------------------------------------
! Get dam indexes
iI = 0; iJ = 0;
iI = oHMC_Vars(iID)%a2iXYDam(iD,2); iJ = oHMC_Vars(iID)%a2iXYDam(iD,1)
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Check dam length
if (a1dVarLDam(iD) .le. 0) then
a1dVarLDam(iD) = 40; ! [m]
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Compute H dam
iFlagD = 0
! Compute TV array length
!iRank = count(mask=oHMC_Vars(iID)%a2dVDam(iD,:).ge.0.0, dim=1)
iRank = SIZE (oHMC_Vars(iID)%a2dVDam,dim=2)
! Cycle on array dimension
do iK = iRank, 2,-1
! Check TV curve availability
if ( (oHMC_Vars(iID)%a2dVDam(iD, iK) .gt. a1dVarVDam(iD)) .and. &
(oHMC_Vars(iID)%a2dVDam(iD, iK-1) .lt. a1dVarVDam(iD)) ) then
! TV curve defined
iFlagD = 1
! Dam H (absolute) --> Linear interpolation
a1dVarHDam(iD) = oHMC_Vars(iID)%a2dLDam(iD,iK-1) + &
( oHMC_Vars(iID)%a2dLDam(iD,iK) - oHMC_Vars(iID)%a2dLDam(iD,iK-1) )* &
( a1dVarVDam(iD) - oHMC_Vars(iID)%a2dVDam(iD,iK-1) )/ &
( oHMC_Vars(iID)%a2dVDam(iD,iK) - oHMC_Vars(iID)%a2dVDam(iD,iK-1) )
! Dam h max (absolute)
!dHMaxDam = oHMC_Vars(iID)%a1dHMaxDam(iD) + oHMC_Vars(iID)%a2dLDam(iD,1)
elseif (oHMC_Vars(iID)%a2dVDam(iD, iK) .lt. 0.0) then
! TV curve undefined
iFlagD = 2
endif
enddo
! TV curve defined but V > VMax defined
if (iFlagD .eq. 0) then !Sono fuori dalla curva invaso volume
iMaxC=MAXLOC(oHMC_Vars(iID)%a2dLDam(iD,:),1)
! Dam H (absolute)
a1dVarHDam(iD) = oHMC_Vars(iID)%a2dLDam(iD,iMaxC) + &
( oHMC_Vars(iID)%a2dLDam(iD,iMaxC) - oHMC_Vars(iID)%a2dLDam(iD,iMaxC-1) )* &
( a1dVarVDam(iD) - oHMC_Vars(iID)%a2dVDam(iD,iMaxC) )/ &
( oHMC_Vars(iID)%a2dVDam(iD,iMaxC) - oHMC_Vars(iID)%a2dVDam(iD,iMaxC-1) )
! Dam h max (absolute)
!dHMaxDam = oHMC_Vars(iID)%a1dHMaxDam(iD) + oHMC_Vars(iID)%a2dLDam(iD,1)
endif
! TV curve undefined
if (iFlagD .eq. 2) then !Non ho la curva invaso volume uso relazione lineare
! Dam h (relative)
a1dVarHDam(iD) = oHMC_Vars(iID)%a1dHMaxDam(iD)*a1dVarVDam(iD)/oHMC_Vars(iID)%a1dVMaxDam(iD)
! Dam h max (relative)
!dHMaxDam = oHMC_Vars(iID)%a1dHMaxDam(iD)
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Check Hdam > Hdammax-DiffDamSpill (if true starting with dam outgoing flow) !
if ( a1dVarHDam(iD) .gt. (oHMC_Vars(iID)%a1dHMaxDam(iD) - oHMC_Namelist(iID)%dDamSpillH) ) then
! Outgoing dam discharge [m^3/s]
a1dVarQoutDam(iD) = 0.0;
a1dVarQoutDam(iD) = 0.385*a1dVarLDam(iD)*( (2*9.81)**0.5)*(a1dVarHDam(iD) - &
( oHMC_Vars(iID)%a1dHMaxDam(iD) - oHMC_Namelist(iID)%dDamSpillH) )**1.5
! Check Qdam <= Qdamoutmax
if (a1dVarQoutDam(iD) .gt. oHMC_Vars(iID)%a1dQcSLDam(iD)) then
a1dVarQoutDam(iD) = oHMC_Vars(iID)%a1dQcSLDam(iD)
endif
! Check Hdam > Hdammax
!!! if ( a1dVarHDam(iD) .gt. oHMC_Vars(iID)%a1dHMaxDam(iD) ) then
! Exceeded volume used for outgoing dam discharge [m^3/s]
!!! a1dVarQoutDam(iD) = (a1dVarVDam(iD) - oHMC_Vars(iID)%a1dVMaxDam(iD))/dDt
! Check max between a1dVarQoutDam(iD) and max outgoing dam flow [comment: to insert a
! a1dVarQoutDam(iD) = max(a1dVarQoutDam(iD), oHMC_Vars(iID)%a1dQcSLDam(iD));
!endif
else
a1dVarQoutDam(iD) = 0.0;
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Impose Vdam not less than dTV * Vdam_max --> static spillage
if ((a1dVarQoutDam(iD)*dDt).gt.(a1dVarVDam(iD) - (oHMC_Namelist(iID)%dTV*oHMC_Vars(iID)%a1dVMaxDam(iD)))) then
! Discharge in [m^3/s] and update dam volume
a1dVarQoutDam(iD) = (a1dVarVDam(iD) - oHMC_Namelist(iID)%dTV*oHMC_Vars(iID)%a1dVMaxDam(iD))/dDt
a1dVarVDam(iD) = oHMC_Namelist(iID)%dTV*oHMC_Vars(iID)%a1dVMaxDam(iD)
else
a1dVarVDam(iD) = a1dVarVDam(iD) - a1dVarQoutDam(iD)*dDt
endif
! If volume exceed the maximum volume that the dam can contain spill add to the spilled water
! the amount of water needed to take the level to and acceptable one
if (a1dVarVDam(iD).gt.oHMC_Vars(iID)%a1dVMaxDam(iD)) then
a1dVarQoutDam(iD) = a1dVarQoutDam(iD) + (a1dVarVDam(iD) - oHMC_Vars(iID)%a1dVMaxDam(iD))/dDt
a1dVarVDam(iD) = oHMC_Vars(iID)%a1dVMaxDam(iD)
endif
! Discharge in [mm]
a1dVarQoutDam(iD) = a1dVarQoutDam(iD)*1000*dDt/(oHMC_Vars(iID)%a2dAreaCell(iI,iJ))
!write(*,*) a1dVarQoutDam(iD)
if (a1dVarVDam(iD) .lt. 0.0) then
a1dVarVDam(iD) = 0.0
a1dVarQoutDam(iD) = 0.0
endif
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
enddo
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Dam_Spilling
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Subroutine to compute dam discharge
subroutine HMC_Phys_Dam_Discharge(iID, iNDam, dDt, &
a1dVarVDam, a1dVarHDam, a1dVarLDam, a1dVarCoeffDam)
!------------------------------------------------------------------------------------------
! Variable(s)
integer(kind = 4) :: iID
integer(kind = 4) :: iD, iK, iRank
integer(kind = 4) :: iNDam
real(kind = 4) :: dDt, dQtmp
real(kind = 4) :: dHMaxDam
character(len = 256) :: sQtmp, sD, sCoeffDam, sHDam, sHMaxDam
integer(kind = 4) :: iFlagD
real(kind = 4), dimension (iNDam) :: a1dVarVDam, a1dVarHDam, a1dVarLDam, a1dVarCoeffDam
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Initialize variable(s)
iFlagD = 0; dQtmp = 0.0;
sQtmp = ''; sCoeffDam = ''; sHDam = ''; sHMaxDam = '';
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Check dam availability
if (iNDam .gt. 0) then
!------------------------------------------------------------------------------------------
! Cycle on dam
do iD = 1, iNDam
!------------------------------------------------------------------------------------------
! Re-initialize dam height
a1dVarHDam(iD) = 0.0
! Update dam volume using observation(s)
where (oHMC_Vars(iID)%a1dVDamObs .gt. -9999.0)
a1dVarVDam = oHMC_Vars(iID)%a1dVDamObs
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Check considering 0.95 total volume
if ( a1dVarVDam(iD) .gt. oHMC_Namelist(iID)%dTV*oHMC_Vars(iID)%a1dVMaxDam(iD) ) then
!------------------------------------------------------------------------------------------
! Check dam length
if (a1dVarLDam(iD) .le. 0) then
a1dVarLDam(iD) = 40; ! [m]
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Compute H dam
iFlagD = 0
! Compute TV array length
iRank = count(mask=oHMC_Vars(iID)%a2dVDam(iD,:).ge.0.0, dim=1)
! Cycle on array dimension
do iK = iRank, 1, -1
! Check TV curve availability
if ( (oHMC_Vars(iID)%a2dVDam(iD, iK) .gt. a1dVarVDam(iD)) .and. &
(oHMC_Vars(iID)%a2dVDam(iD, iK-1) .lt. a1dVarVDam(iD)) ) then
! TV curve defined
iFlagD = 1
! Dam H (absolute) --> Linear interpolation
a1dVarHDam(iD) = oHMC_Vars(iID)%a2dLDam(iD,iK-1) + &
( oHMC_Vars(iID)%a2dLDam(iD,iK) - oHMC_Vars(iID)%a2dLDam(iD,iK-1) )* &
( a1dVarVDam(iD) - oHMC_Vars(iID)%a2dVDam(iD,iK-1) )/ &
( oHMC_Vars(iID)%a2dVDam(iD,iK) - oHMC_Vars(iID)%a2dVDam(iD,iK-1) )
! Dam h max (absolute)
dHMaxDam = oHMC_Vars(iID)%a1dHMaxDam(iD) + oHMC_Vars(iID)%a2dLDam(iD,1)
elseif (oHMC_Vars(iID)%a2dVDam(iD, iK) .lt. 0.0) then
! TV curve undefined
iFlagD = 2
endif
enddo
! TV curve defined but V > VMax defined
if (iFlagD .eq. 0) then !Sono fuori dalla curva invaso volume
! Dam H (absolute)
a1dVarHDam(iD) = oHMC_Vars(iID)%a2dLDam(iD,iRank) + &
( oHMC_Vars(iID)%a2dLDam(iD,iRank) - oHMC_Vars(iID)%a2dLDam(iD,iRank-1) )* &
( a1dVarVDam(iD) - oHMC_Vars(iID)%a2dVDam(iD,iRank) )/ &
( oHMC_Vars(iID)%a2dVDam(iD,iRank) - oHMC_Vars(iID)%a2dVDam(iD,iRank-1) )
! Dam h max (absolute)
dHMaxDam = oHMC_Vars(iID)%a1dHMaxDam(iD) + oHMC_Vars(iID)%a2dLDam(iD,1)
endif
! TV curve undefined
if (iFlagD .eq. 2) then !Non ho la curva invaso volume uso relazione lineare
! Dam h (relative)
a1dVarHDam(iD) = oHMC_Vars(iID)%a1dHMaxDam(iD)*a1dVarVDam(iD)/oHMC_Vars(iID)%a1dVMaxDam(iD)
! Dam h max (relative)
dHMaxDam = oHMC_Vars(iID)%a1dHMaxDam(iD)
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Check Hdam > Hdammax-DiffDamSpill (if true starting with dam outgoing flow)
if ( a1dVarHDam(iD) .gt. (dHMaxDam - oHMC_Namelist(iID)%dDamSpillH) ) then
! Outgoing dam discharge [m^3/s]
dQtmp = 0.0;
dQtmp = 0.385*a1dVarLDam(iD)*( (2*9.81)**0.5)*(a1dVarHDam(iD) - &
( dHMaxDam - oHMC_Namelist(iID)%dDamSpillH) )**1.5
! Check Qdam <= Qdamoutmax
if (dQtmp .gt. oHMC_Vars(iID)%a1dQcSLDam(iD)) then
dQtmp = oHMC_Vars(iID)%a1dQcSLDam(iD)
endif
! Check Hdam > Hdammax
if ( a1dVarHDam(iD) .gt. dHMaxDam ) then
! Exceeded volume used for outgoing dam discharge [m^3/s]
dQtmp = (a1dVarVDam(iD) - oHMC_Vars(iID)%a1dVMaxDam(iD))/3600
! Check max between dQtmp and max outgoing dam flow
dQtmp = max(dQtmp, oHMC_Vars(iID)%a1dQcSLDam(iD));
endif
else
dQtmp = 0.0;
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Compute dam outgoing flow coefficient [1/s]
if (dQtmp .gt. 0.0) then
a1dVarCoeffDam(iD) = dQtmp/(a1dVarVDam(iD) - &
oHMC_Namelist(iID)%dTV*oHMC_Vars(iID)%a1dVMaxDam(iD))
else
a1dVarCoeffDam(iD) = 0.000
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Debug
!write(*,*) '==============================='
!write(*,*) iRank, iD
!write(*,*) oHMC_Vars(iID)%a2dLDam(iD,iRank)
!write(*,*) oHMC_Vars(iID)%a2dVDam(iD,iRank)
!write(*,*) oHMC_Vars(iID)%a2dLDam(iD,iRank-1)
!write(*,*) oHMC_Vars(iID)%a2dVDam(iD,iRank-1)
!write(*,*) a1dVarLDam(iD)
!write(*,*) a1dVarHDam(iD)
!write(*,*) oHMC_Vars(iID)%a1dHMaxDam(iD)
!write(*,*) a1dVarVDam(iD)
!write(*,*) oHMC_Vars(iID)%a1dVMaxDam(iD)
!write(*,*) dQtmp
!write(*,*) a1dVarCoeffDam(iD)
!write(*,*) '==============================='
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Info dam coefficient updating
write(sQtmp, *) dQtmp; write(sD, *) iD; write(sCoeffDam, *) a1dVarCoeffDam(iD);
write(sHDam, *) a1dVarHDam(iD); write(sHMaxDam, *) oHMC_Vars(iID)%a1dHMaxDam(iD)
call mprintf(.true., iINFO_Verbose, &
' Phys :: Convolution :: Discharge Dam :: NDam: '//trim(sD)//' [-]'// &
' QDam: '//trim(sQtmp)//' [m^3/s]'// &
' CoeffDam: '//trim(sCoeffDam)//' [-]'// &
' HDam: '//trim(sHDam)//' [m]'// &
' HMaxDam: '//trim(sHMaxDam)//' [m]')
!------------------------------------------------------------------------------------------
enddo
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Dam_Discharge
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Subroutine to compute lake tank
subroutine HMC_Phys_Lake_Tank(iID, iRows, iCols, dDtDataForcing, iNLake, &
a1dVarVLake, a1dVarQoutLake, &
a2dVarHydro, a2dVarFlowDeep)
!------------------------------------------------------------------------------------------
! Variable(s)
integer(kind = 4) :: iID
integer(kind = 4) :: iL, iRows, iCols, iNLake
integer(kind = 4) :: iI, iII, iIII, iJ, iJJ, iJJJ
real(kind = 4) :: dDtDataForcing
real(kind = 4), dimension (iRows, iCols) :: a2dVarFlowDeep, a2dVarHydro
real(kind = 4), dimension (iNLake) :: a1dVarVLake, a1dVarQoutLake
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Initialize variable(s)
a1dVarQoutLake = 0.0
! Check lake availability
if (iNLake .gt. 0) then
! Cycle on lake
iL = 0
do iL = 1, iNLake
! Check lake volume
if ( a1dVarVLake(iL) .gt. oHMC_Vars(iID)%a1dVMinLake(iL) ) then
iI = 0; iJ = 0;
iI = oHMC_Vars(iID)%a2iXYLake(iL,2); iJ = oHMC_Vars(iID)%a2iXYLake(iL,1);
! Output lake discharge [mm/h]
a1dVarQoutLake(iL) = (a1dVarVLake(iL) - oHMC_Vars(iID)%a1dVMinLake(iL))* &
(oHMC_Vars(iID)%a1dCostLake(iL)/3600)/(oHMC_Vars(iID)%a2dAreaCell(iI,iJ))*3600*1000
! Lake outgoing flow
a1dVarVLake(iL) = a1dVarVLake(iL) - &
(a1dVarVLake(iL) - oHMC_Vars(iID)%a1dVMinLake(iL))* &
(oHMC_Vars(iID)%a1dCostLake(iL)/3600)*dDtDataForcing
! Defining flow directions
iII = int((int(oHMC_Vars(iID)%a2iPNT(iI,iJ)) - 1)/3) - 1;
iJJ = oHMC_Vars(iID)%a2iPNT(iI,iJ) - 5 - 3*iII
iIII = iI + iII; iJJJ = iJ + iJJ
! Lake equation to fill lake [mm]
if (oHMC_Vars(iID)%a1dCodeLake(iL).gt.0) then
where (oHMC_Vars(iID)%a2iChoice.eq.oHMC_Vars(iID)%a1dCodeLake(iL) .and. oHMC_Vars(iID)%a2dDem.gt.0.0)
! Lake volume to lake mean level
a2dVarHydro = a1dVarVLake(iL)/(oHMC_Vars(iID)%a1dCodeLake(iL)*oHMC_Vars(iID)%a2dAreaCell(iI,iJ))*1000
endwhere
endif
! Add lake flow to deep flow (summed to intensity in Horton subroutine) [mm/dt]
a2dVarFlowDeep(iIII, iJJJ) = a2dVarFlowDeep(iIII, iJJJ) + a1dVarQoutLake(iL)*dDtDataForcing/3600
endif
enddo
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Lake_Tank
!------------------------------------------------------------------------------------------
end module HMC_Module_Phys_HydraulicStructure
!------------------------------------------------------------------------------------------