-
Notifications
You must be signed in to change notification settings - Fork 4
Expand file tree
/
Copy pathHMC_Module_Phys_Convolution_Apps_SubFlow.f90
More file actions
376 lines (281 loc) · 17 KB
/
HMC_Module_Phys_Convolution_Apps_SubFlow.f90
File metadata and controls
376 lines (281 loc) · 17 KB
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
!------------------------------------------------------------------------------------
! File: HMC_Module_Phys_Convolution_Apps_SubFlow.f90
!
! Author(s): Fabio Delogu, Francesco Silvestro, Simone Gabellani
! Date: 20260320
!
! Convolution Apps SubFlow subroutine(s) for HMC model
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Module Header
module HMC_Module_Phys_Convolution_Apps_SubFlow
!------------------------------------------------------------------------------------
! External module(s)
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 hypodermic flow
subroutine HMC_Phys_Convolution_Apps_SubFlow(iID, iRows, iCols, dDtDataForcing, dDtAct, iNDam)
!------------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iID, iRows, iCols
integer(kind = 4) :: iNDam
real(kind = 4) :: dDtDataForcing, dDtAct
integer(kind = 4) :: iI, iII, iIII, iJ, iJJ, iJJJ, iD
!integer(kind = 4) :: iVarPNT
integer(kind = 4) :: iFlagFlowDeep
real(kind = 4) :: dDtSubflow
real(kind = 4) :: dRate, dRateMin
real(kind = 4), dimension (iRows, iCols) :: a2dVarVTot, a2dVarVTotStep, a2dVarVLoss
real(kind = 4), dimension (iRows, iCols) :: a2dVarFlowExf
!real(kind = 4), dimension (iRows, iCols) :: a2dVarVSub
real(kind = 4), dimension (iNDam) :: a1dVarVDam
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Initialize variable(s)
dDtSubflow = 0.0; dRate = 0.0; dRateMin = 0.0; iFlagFlowDeep = 0;
a2dVarVTot = 0.0; a2dVarVTotStep = 0.0; a2dVarVLoss = 0.0;
a2dVarFlowExf = 0.0
a1dVarVDam = 0.0
oHMC_Vars(iID)%a2dFlowExf = 0.0
!iVarPNT = 0; a2dVarVSub = 0.0;
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Integrating step (Subflow)
dDtSubflow = dDtAct
! Hypodermic flow minimum rate
dRateMin = oHMC_Namelist(iID)%dRateMin
iFlagFlowDeep = oHMC_Namelist(iID)%iFlagFlowDeep
! Dam(s) data from global declaration
a1dVarVDam = oHMC_Vars(iID)%a1dVDam
! Variable(s) from global declaration
a2dVarVTot = oHMC_Vars(iID)%a2dVTot
!a2dVarVSub = oHMC_Vars(iID)%a2dVSub
a2dVarVLoss = oHMC_Vars(iID)%a2dVLoss
a2dVarFlowExf = oHMC_Vars(iID)%a2dFlowExf
! Info start
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: SubFlow ... ' )
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Debug
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' ========= SUBFLOW START ========= ')
call mprintf(.true., iINFO_Extra, checkvar(a2dVarVTot, oHMC_Vars(iID)%a2iMask, 'VTOT START ') )
!call mprintf(.true., iINFO_Extra, checkvar(a2dVarVSub, oHMC_Vars(iID)%a2iMask, 'VSUB START ') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarVLoss, oHMC_Vars(iID)%a2iMask, 'VLOSS START ') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarFlowExf, oHMC_Vars(iID)%a2iMask, 'FLOWEXF START ') )
call mprintf(.true., iINFO_Extra, '')
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Conditions on total and loss volume
where( (oHMC_Vars(iID)%a2iMask.gt.0.0).and.(a2dVarVTot.lt.0.0) ) a2dVarVTot = 0.0
where( (oHMC_Vars(iID)%a2iMask.gt.0.0).and.(a2dVarVLoss.lt.0.0) ) a2dVarVLoss = 0.0
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Check dam availability
if (iNDam .gt. 0) then
! Set VTot equal zero into Dam(s) cell(s)
do iD = 1,iNDam
iI = 0; iJ = 0;
iI = oHMC_Vars(iID)%a2iXYDam(iD, 2); iJ = oHMC_Vars(iID)%a2iXYDam(iD, 1)
!a1dVarVDam(iD) = a1dVarVDam(iD) + a2dVarVTot(iI,iJ)*(oHMC_Vars(iID)%a2dAreaCell(iI,iJ))/1000
a2dVarVTot(iI,iJ) = 0.0
enddo
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Cycles on each pixels (using grid or index methods)
! NOTE:
! Calcola il volume di uscita dalla cella nei due casi: a2dV > o < di a2dS;
! Qsup � la portata che esce dalla parte superiore della cella e si aggiunge al deflusso superficiale
! il contatore punta alla cella successiva (controllare se vale per l'ultima cella)
! Evaluate flag to use index or grid routing
if (oHMC_Namelist(iID)%iFlagRoutingType .eq. 2 .and. &
oHMC_Vars(iID)%bRoutingIndex) then
call compute_routing_index( &
iID, iRows, iCols, iFlagFlowDeep, dRateMin, &
a2dVarVTotStep, a2dVarVLoss )
else
call compute_routing_grid( &
iID, iRows, iCols, iFlagFlowDeep, dRateMin, &
a2dVarVTotStep, a2dVarVLoss )
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Updating total volume
a2dVarVTot = a2dVarVTot + a2dVarVTotStep
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Checking total volume and calculating exfiltration volume
where ( (oHMC_Vars(iID)%a2iMask.gt.0.0) .and. (a2dVarVTot.gt.oHMC_Vars(iID)%a2dS) )
! Calculating esfiltration flow [m/seconds]
a2dVarFlowExf = (a2dVarVTot - oHMC_Vars(iID)%a2dS)/(1000.0*dDtSubflow) !in m/sec
! Updating total volume information
a2dVarVTot = 1.0*oHMC_Vars(iID)%a2dS
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Hypodermic flow to fill cell(s)-lake of the dam
if (iNDam .gt. 0) then
do iD = 1,iNdam
iI = 0; iJ = 0;
iI = oHMC_Vars(iID)%a2iXYDam(iD,2); iJ = oHMC_Vars(iID)%a2iXYDam(iD,1)
! Amount of upstream volume at dam section. V set zero at the subroutine begin
a1dVarVDam(iD) = a1dVarVDam(iD) + a2dVarVTot(iI,iJ)*(oHMC_Vars(iID)%a2dAreaCell(iI,iJ))/1000
enddo
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Debug
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, '')
call mprintf(.true., iINFO_Extra, checkvar(a2dVarVTot, oHMC_Vars(iID)%a2iMask, 'VTOT END ') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarVTotStep, oHMC_Vars(iID)%a2iMask, 'VTOTSTEP START ') )
!call mprintf(.true., iINFO_Extra, checkvar(a2dVarVSub, oHMC_Vars(iID)%a2iMask, 'VSUB END ') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarVLoss, oHMC_Vars(iID)%a2iMask, 'VLOSS END ') )
call mprintf(.true., iINFO_Extra, checkvar(a2dVarFlowExf, oHMC_Vars(iID)%a2iMask, 'FLOWEXF END ') )
call mprintf(.true., iINFO_Extra, ' ========= SUBFLOW END ========= ')
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Updating variable(s)
oHMC_Vars(iID)%a2dVTot = a2dVarVTot
oHMC_Vars(iID)%a2dVLoss = a2dVarVLoss
oHMC_Vars(iID)%a2dFlowExf = a2dVarFlowExf
oHMC_Vars(iID)%a1dVDam = a1dVarVDam
! Info end
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: SubFlow ... OK' )
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Convolution_Apps_SubFlow
!------------------------------------------------------------------------------------------
subroutine compute_routing_index( &
iID, iRows, iCols, iFlagFlowDeep, dRateMin, a2dVarVTotStep, a2dVarVLoss )
!------------------------------------------------------------------------------------------
implicit none
integer(kind = 4), intent(in) :: iID
integer(kind = 4), intent(in) :: iRows, iCols
integer(kind = 4), intent(in) :: iFlagFlowDeep
real(kind = 4), intent(in) :: dRateMin
real(kind = 4), dimension(iRows, iCols), intent(inout) :: a2dVarVTotStep
real(kind = 4), dimension(iRows, iCols), intent(inout) :: a2dVarVLoss
integer(kind = 4) :: iCell
integer(kind = 4) :: iI, iJ
integer(kind = 4) :: iIII, iJJJ
integer(kind = 4) :: iActiveCells
real(kind = 4) :: dRate, dRateRescaling
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! No rescaling if dRateRescaling>=0.99; otherwise, dRateRescaling will be the new maximum
dRateRescaling = oHMC_Namelist(iID)%dRateRescaling
if (dRateRescaling .ge. 0.99) dRateRescaling = 0.99
if (oHMC_Namelist(iID)%dRateRescaling .ge. 0.99) &
oHMC_Namelist(iID)%dRateRescaling = 0.99
! iterate over active cells
iActiveCells = oHMC_Vars(iID)%iRoutingCellsActive
do iCell = 1, iActiveCells
iI = oHMC_Vars(iID)%a1iRoutingSrcI(iCell)
iJ = oHMC_Vars(iID)%a1iRoutingSrcJ(iCell)
iIII = oHMC_Vars(iID)%a1iRoutingDstI(iCell)
iJJJ = oHMC_Vars(iID)%a1iRoutingDstJ(iCell)
! Safety check: keep behavior explicit and robust
if (iIII.ge.1 .and. iIII.le.iRows .and. iJJJ.ge.1 .and. iJJJ.le.iCols) then
! Calculating VTot and VLoss using flowdeep condition
if (iFlagFlowDeep .eq. 0) then
! VTot (VLoss == 0)
a2dVarVTotStep(iIII,iJJJ) = a2dVarVTotStep(iIII,iJJJ) + &
oHMC_Vars(iID)%a2dVSub(iI,iJ)
else
! Rate definition
dRate = sin(oHMC_Vars(iID)%a2dBeta(iI,iJ))
! Checking rate value
if (dRate .gt. 0.99) dRate = 0.99
if (dRate .lt. dRateMin) dRate = dRateMin
! Rescaling dRate according to dRateRescaling set in the namelist
! (dRate=fraction to hypodermic flow; (1-dRate)=fraction to percolation)
dRate = dRateMin + (dRate - dRateMin)/(0.99 - dRateMin) * &
(dRateRescaling - dRateMin)
! VTot
a2dVarVTotStep(iIII,iJJJ) = a2dVarVTotStep(iIII,iJJJ) + &
oHMC_Vars(iID)%a2dVSub(iI,iJ)*dRate
! VLoss
a2dVarVLoss(iIII,iJJJ) = a2dVarVLoss(iIII,iJJJ) + &
oHMC_Vars(iID)%a2dVSub(iI,iJ)*(1.0 - dRate)
endif
endif
enddo
!------------------------------------------------------------------------------------------
end subroutine compute_routing_index
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Subroutine for routing hypodermic flow using full-grid iteration
subroutine compute_routing_grid( &
iID, iRows, iCols, iFlagFlowDeep, dRateMin, a2dVarVTotStep, a2dVarVLoss )
!------------------------------------------------------------------------------------------
! variable(s)
integer(kind = 4), intent(in) :: iID
integer(kind = 4), intent(in) :: iRows, iCols
integer(kind = 4), intent(in) :: iFlagFlowDeep
real(kind = 4), intent(in) :: dRateMin
real(kind = 4), dimension(iRows, iCols), intent(inout) :: a2dVarVTotStep
real(kind = 4), dimension(iRows, iCols), intent(inout) :: a2dVarVLoss
integer(kind = 4) :: iI, iII, iIII, iJ, iJJ, iJJJ
real(kind = 4) :: dRate, dRateRescaling
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! No rescaling if dRateRescaling>=0.99; otherwise, dRateRescaling will be the new maximum
dRateRescaling = oHMC_Namelist(iID)%dRateRescaling
if (dRateRescaling .ge. 0.99) dRateRescaling = 0.99
iI = 0; iJ = 0;
do iJ = 1, iCols
do iI = 1, iRows
! DEM condition
if (oHMC_Vars(iID)%a2iMask(iI,iJ) .gt. 0.0) then
! Defining flow directions
iII = int((int(oHMC_Vars(iID)%a2iPNT(iI,iJ)) - 1)/3) - 1
iJJ = int(oHMC_Vars(iID)%a2iPNT(iI,iJ)) - 5 - 3*iII
iIII = iI + iII
iJJJ = iJ + iJJ
! Check destination cell inside the domain
if (iIII .ge. 1 .and. iIII .le. iRows .and. &
iJJJ .ge. 1 .and. iJJJ .le. iCols) then
! Calculating VTot and VLoss using flowdeep condition
if (iFlagFlowDeep .eq. 0) then
! VTot (VLoss == 0)
a2dVarVTotStep(iIII,iJJJ) = a2dVarVTotStep(iIII,iJJJ) + &
oHMC_Vars(iID)%a2dVSub(iI,iJ)
else
! Rate definition
dRate = sin(oHMC_Vars(iID)%a2dBeta(iI,iJ))
! Checking rate value
if (dRate .gt. 0.99) dRate = 0.99
if (dRate .lt. dRateMin) dRate = dRateMin
! Rescaling dRate according to dRateRescaling set in the namelist
! (dRate = fraction to hypodermic flow; (1-dRate) = fraction to percolation)
dRate = dRateMin + (dRate - dRateMin)/(0.99 - dRateMin) * &
(dRateRescaling - dRateMin)
! VTot
a2dVarVTotStep(iIII,iJJJ) = a2dVarVTotStep(iIII,iJJJ) + &
oHMC_Vars(iID)%a2dVSub(iI,iJ)*dRate
! VLoss
a2dVarVLoss(iIII,iJJJ) = a2dVarVLoss(iIII,iJJJ) + &
oHMC_Vars(iID)%a2dVSub(iI,iJ)*(1.0 - dRate)
endif
endif
endif
enddo
enddo
!------------------------------------------------------------------------------------------
end subroutine compute_routing_grid
!------------------------------------------------------------------------------------------
end module HMC_Module_Phys_Convolution_Apps_SubFlow
!------------------------------------------------------------------------------------------