-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathHMC_Module_Phys_Convolution_Apps_Discharge.f90
395 lines (313 loc) · 21.8 KB
/
HMC_Module_Phys_Convolution_Apps_Discharge.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
!------------------------------------------------------------------------------------
! File: HMC_Module_Phys_Convolution_Apps_Discharge.f90
!
! Author(s): Fabio Delogu, Francesco Silvestro, Simone Gabellani
! Date: 20190410
!
! Convolution Apps Discharge subroutine(s) for HMC model
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Module Header
module HMC_Module_Phys_Convolution_Apps_Discharge
!------------------------------------------------------------------------------------
! External module(s)
use HMC_Module_Namelist, only: oHMC_Namelist
use HMC_Module_Vars_Loader, only: oHMC_Vars
use HMC_Module_Tools_Debug
use HMC_Module_Phys_HydraulicStructure, only: HMC_Phys_Dam_Spilling, &
HMC_Phys_Dam_Discharge, &
HMC_Phys_Lake_Tank
use HMC_Module_Tools_Generic, only: getIntValue, getIntRange
! Implicit none for all subroutines in this module
implicit none
!------------------------------------------------------------------------------------------
contains
!------------------------------------------------------------------------------------
! Subroutine for calculating discharge channel network type
subroutine HMC_Phys_Convolution_Apps_Discharge_ChannelNetwork(iID, iRows, iCols, iNSection, &
dDtDataForcing, dDtAct, iNTime, iTq, dDtMax)
!------------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iID, iRows, iCols, iTAct
real(kind = 4) :: dDtDataForcing, dDtMax, dDtAct, dDtDischarge
integer(kind = 4), intent(inout) :: iNTime, iTq
integer(kind = 4) :: iI, iJ, iT, iS, iD
integer(kind = 4) :: iNSection
!integer(kind = 4) :: iFlagVarUc
real(kind = 4) :: dRate
!real(kind = 4), dimension (iRows, iCols) :: a2dVarHydroPrev, a2dVarHydro, a2dVarIntensity
real(kind = 4), dimension (iRows, iCols) :: a2dVarQDisOut, a2dVarQVolOut, a2dVarQTot, a2dVarQout
real(kind = 4), dimension (iNSection) :: a1dVarQoutSection
character(len = 10) :: sVarQoutSectionFirst, sVarQoutSectionMax
character(len = 10), parameter :: sFMTVarQoutSection = "(F10.2)"
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Variable(s) initialization
a2dVarQDisOut = 0.0; a2dVarQVolOut = 0.0; a2dVarQTot = 0.0; a2dVarQout = 0.0;
a1dVarQoutSection = 0.0;
! Re-initializing Q distributed discharge matrix (in global memory)
oHMC_Vars(iID)%a2dQout = 0.0; ! Distributed discharge
oHMC_Vars(iID)%a1dQoutSection = 0.0; ! Section discharge OUT
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Flags
!iFlagVarUc = oHMC_Namelist(iID)%iFlagVarUc
! Temporal model step
iT = int(oHMC_Vars(iID)%iTime)
! Variable(s) time dependent from global declaration
a2dVarQDisOut = oHMC_Vars(iID)%a2dQDisOut ! Qout
a2dVarQVolOut = oHMC_Vars(iID)%a2dQVolOut ! Qtmp
a2dVarQTot = oHMC_Vars(iID)%a2dQTot
a2dVarQout = oHMC_Vars(iID)%a2dQout
a1dVarQoutSection = oHMC_Vars(iID)%a1dQoutSection
! Info start
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: Discharge ... ' )
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Debug
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' ========= DISCHARGE START =========== ')
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQDisOut, oHMC_Vars(iID)%a2iMask, 'QDIS OUT START (Qout)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQVolOut, oHMC_Vars(iID)%a2iMask, 'QVOL OUT START (Qtmp)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQTot, oHMC_Vars(iID)%a2iMask, 'QTOT START (Q)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQout, oHMC_Vars(iID)%a2iMask, 'QOUT START (Qmap)') )
call mprintf(.true., iINFO_Extra, checkvar(oHMC_Vars(iID)%a2dUc, oHMC_Vars(iID)%a2iMask, 'UC ') )
call mprintf(.true., iINFO_Extra, ' ')
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Integrating step (SurfaceFlow)
dDtDischarge = dDtAct
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Discharge in channel cells [mm/s]
where( (oHMC_Vars(iID)%a2iChoice.eq.1) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarQTot = a2dVarQTot + a2dVarQVolOut
endwhere
! Discharge in hills cells [mm/s]
where( (oHMC_Vars(iID)%a2iChoice.eq.0) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarQTot = a2dVarQTot + a2dVarQVolOut
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Compute discharge from mm/s to m^3/s (constant is in 1/h)
if( (real(iTq)*dDtDischarge) .ge. (dDtDataForcing - dDtMax*1.001) ) then ! 1.001 for numerical approx
!------------------------------------------------------------------------------------------
! Compute distributed discharge (total and for each step)
where( oHMC_Vars(iID)%a2iMask.gt.0.0 )
a2dVarQTot = a2dVarQTot/(real(iTq)*1000)*oHMC_Vars(iID)%a2dAreaCell
a2dVarQout = a2dVarQTot
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Calculating discharge in selected outlet(s)
a1dVarQoutSection = 0.0
do iS = 1, iNSection
! Extracting outlet(s) indexes
iI = 0; iJ = 0;
iI = oHMC_Vars(iID)%a2iXYSection(iS, 2); iJ = oHMC_Vars(iID)%a2iXYSection(iS, 1)
! Get section discharge
a1dVarQoutSection(iS) = a2dVarQTot(iI, iJ)
! Subsurface flux
dRate = sin(oHMC_Vars(iID)%a2dBeta(iI, iJ))
if (dRate.gt.1) dRate = 0.99
if (dRate.lt.0) dRate = 0.1
if (oHMC_Vars(iID)%a2dWTableMax(iI,iJ).eq.0.0) dRate = 1.0 !Celle dove non ho falda
enddo
! Re-initialize qtot array and q counter
!write(*,*) a1dVarQoutSection; write(*,*) a1dVarQinDam
a2dVarQTot = 0.0
iTq = 0
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Discharge information time step
write(sVarQoutSectionFirst, sFMTVarQoutSection) a1dVarQoutSection(1)
write(sVarQoutSectionMax, sFMTVarQoutSection) maxval(a1dVarQoutSection)
call mprintf(.true., iINFO_Basic, ' Phys :: Convolution :: Discharge ::'// &
' Q SecFirst: '//sVarQoutSectionFirst//' [m^3/s]'// &
' Q SecMax: '//sVarQoutSectionMax//' [m^3/s]')
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Debug
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' ')
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQDisOut, oHMC_Vars(iID)%a2iMask, 'QDIS OUT END (Qout)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQVolOut, oHMC_Vars(iID)%a2iMask, 'QVOL OUT END (Qtmp)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQTot, oHMC_Vars(iID)%a2iMask, 'QTOT END (Q)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQout, oHMC_Vars(iID)%a2iMask, 'QOUT END (Qmap)') )
call mprintf(.true., iINFO_Extra, ' ========= DISCHARGE END =========== ')
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Updating model global variable(s)
oHMC_Vars(iID)%a2dQTot = a2dVarQTot
oHMC_Vars(iID)%a2dQout = a2dVarQout ! Distributed discharge
oHMC_Vars(iID)%a1dQoutSection = a1dVarQoutSection ! Section discharge OUT
! Info end
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: Discharge ... OK' )
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Convolution_Apps_Discharge_ChannelNetwork
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Subroutine for calculating discharge channel fraction type
subroutine HMC_Phys_Convolution_Apps_Discharge_ChannelFraction(iID, iRows, iCols, iNSection, &
dDtDataForcing, dDtAct, iNTime, iTq, dDtMax)
!------------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iID, iRows, iCols, iTAct, iFlagFlood
real(kind = 4) :: dDtDataForcing, dDtMax, dDtAct, dDtDischarge
integer(kind = 4), intent(inout) :: iNTime, iTq
integer(kind = 4) :: iI, iJ, iT, iS, iD
integer(kind = 4) :: iNSection
!integer(kind = 4) :: iFlagVarUc
real(kind = 4) :: dRate
!real(kind = 4), dimension (iRows, iCols) :: a2dVarHydroPrev, a2dVarHydro, a2dVarIntensity
real(kind = 4), dimension (iRows, iCols) :: a2dVarQDisOut, a2dVarQVolOut, a2dVarQTot, a2dVarQout
real(kind = 4), dimension (iRows, iCols) :: a2dVarQfloodCR, a2dVarQfloodCL, a2dVarQfloodIR, a2dVarQfloodIL
real(kind = 4), dimension (iNSection) :: a1dVarQoutSection
character(len = 10) :: sVarQoutSectionFirst, sVarQoutSectionMax
character(len = 10), parameter :: sFMTVarQoutSection = "(F10.2)"
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Variable(s) initialization
a2dVarQDisOut = 0.0; a2dVarQVolOut = 0.0; a2dVarQTot = 0.0; a2dVarQout = 0.0;
a2dVarQfloodCR = 0.0; a2dVarQfloodCL = 0.0; a2dVarQfloodIR = 0.0; a2dVarQfloodIL = 0.0;
a1dVarQoutSection = 0.0;
! Re-initializing Q distributed discharge matrix (in global memory)
oHMC_Vars(iID)%a2dQout = 0.0; ! Distributed discharge
oHMC_Vars(iID)%a1dQoutSection = 0.0; ! Section discharge OUT
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Flooding flag
iFlagFlood = oHMC_Namelist(iID)%iFlagFlood
! Temporal model step
iT = int(oHMC_Vars(iID)%iTime)
! Variable(s) time dependent from global declaration
a2dVarQDisOut = oHMC_Vars(iID)%a2dQC ! Q channel 2G (Forse non serve)
a2dVarQVolOut = oHMC_Vars(iID)%a2dQC ! Q channel 2G
a2dVarQTot = oHMC_Vars(iID)%a2dQTot
a2dVarQout = oHMC_Vars(iID)%a2dQout
a2dVarQfloodCR = oHMC_Vars(iID)%a2dQfloodCR !cumulative flooding right bank
a2dVarQfloodIR = oHMC_Vars(iID)%a2dQfloodIR !convolution step flooding right bank
a2dVarQfloodCL = oHMC_Vars(iID)%a2dQfloodCL !cumulative flooding left bank
a2dVarQfloodIL = oHMC_Vars(iID)%a2dQfloodIL !convolution step flooding left bank
a1dVarQoutSection = oHMC_Vars(iID)%a1dQoutSection
! Info start
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: Discharge ... ' )
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Debug
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' ========= DISCHARGE START =========== ')
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQDisOut, oHMC_Vars(iID)%a2iMask, 'QDIS OUT START (Qout)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQVolOut, oHMC_Vars(iID)%a2iMask, 'QVOL OUT START (Qtmp)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQTot, oHMC_Vars(iID)%a2iMask, 'QTOT START (Q)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQout, oHMC_Vars(iID)%a2iMask, 'QOUT START (Qmap)') )
call mprintf(.true., iINFO_Extra, checkvar(oHMC_Vars(iID)%a2dUc, oHMC_Vars(iID)%a2iMask, 'UC ') )
call mprintf(.true., iINFO_Extra, ' ')
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Integrating step (SurfaceFlow)
dDtDischarge = dDtAct
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Discharge in channel part of cells [m^3/s]
where( (oHMC_Vars(iID)%a2iChoice.le.1) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarQTot = a2dVarQTot + a2dVarQVolOut
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
!If Flooding is activated build the Average Qflooding maps
if (iFlagFlood.eq.1) then
where( (oHMC_Vars(iID)%a2iChoice.le.1) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarQfloodCR = a2dVarQfloodCR + a2dVarQfloodIR
a2dVarQfloodCL = a2dVarQfloodCL + a2dVarQfloodIL
endwhere
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Compute discharge from mm/s to m^3/s (constant is in 1/h)
if( (real(iTq)*dDtDischarge) .ge. (dDtDataForcing - dDtMax*1.001) ) then ! 1.001 for numerical approx
!------------------------------------------------------------------------------------------
! Compute distributed discharge (total and for each step)
where( oHMC_Vars(iID)%a2iMask.gt.0.0 )
a2dVarQTot = a2dVarQTot/(real(iTq))
a2dVarQout = a2dVarQTot
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
!If Flooding is activated build the Average Qflooding maps
if (iFlagFlood.eq.1) then
where( (oHMC_Vars(iID)%a2iChoice.le.1) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarQfloodCR = a2dVarQfloodCR / (real(iTq))
a2dVarQfloodCL = a2dVarQfloodCL / (real(iTq))
endwhere
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Calculating discharge in selected outlet(s)
a1dVarQoutSection = 0.0
do iS = 1, iNSection
! Extracting outlet(s) indexes
iI = 0; iJ = 0;
iI = oHMC_Vars(iID)%a2iXYSection(iS, 2); iJ = oHMC_Vars(iID)%a2iXYSection(iS, 1)
! Get section discharge
a1dVarQoutSection(iS) = a2dVarQTot(iI, iJ)
!a1dVarQoutSection(iS) = oHMC_Vars(iID)%a2dQC(iI, iJ)
! Subsurface flux
dRate = sin(oHMC_Vars(iID)%a2dBeta(iI, iJ))
if (dRate.gt.1) dRate = 0.99
if (dRate.lt.0) dRate = 0.1
if (oHMC_Vars(iID)%a2dWTableMax(iI,iJ).eq.0.0) dRate = 1.0 !Celle dove non ho falda
enddo
! Re-initialize qtot array and q counter
! write(*,*) a1dVarQoutSection; write(*,*) a1dVarQinDam
a2dVarQTot = 0.0
iTq = 0
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Discharge information time step
write(sVarQoutSectionFirst, sFMTVarQoutSection) a1dVarQoutSection(1)
write(sVarQoutSectionMax, sFMTVarQoutSection) maxval(a1dVarQoutSection)
call mprintf(.true., iINFO_Basic, ' Phys :: Convolution :: Discharge ::'// &
' Q SecFirst: '//sVarQoutSectionFirst//' [m^3/s]'// &
' Q SecMax: '//sVarQoutSectionMax//' [m^3/s]')
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Debug
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' ')
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQDisOut, oHMC_Vars(iID)%a2iMask, 'QDIS OUT END (Qout)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQVolOut, oHMC_Vars(iID)%a2iMask, 'QVOL OUT END (Qtmp)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQTot, oHMC_Vars(iID)%a2iMask, 'QTOT END (Q)') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarQout, oHMC_Vars(iID)%a2iMask, 'QOUT END (Qmap)') )
call mprintf(.true., iINFO_Extra, ' ========= DISCHARGE END =========== ')
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Updating model global variable(s)
oHMC_Vars(iID)%a2dQTot = a2dVarQTot
oHMC_Vars(iID)%a2dQout = a2dVarQout ! Distributed discharge
oHMC_Vars(iID)%a1dQoutSection = a1dVarQoutSection ! Section discharge OUT
oHMC_Vars(iID)%a2dQfloodCR = a2dVarQfloodCR !cumulative flooding right bank
oHMC_Vars(iID)%a2dQfloodCL = a2dVarQfloodCL !cumulative flooding left bank
! Info end
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: Discharge ... OK' )
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Convolution_Apps_Discharge_ChannelFraction
!------------------------------------------------------------------------------------------
end module HMC_Module_Phys_Convolution_Apps_Discharge
!------------------------------------------------------------------------------------------