-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathHMC_Module_Phys_Convolution_Apps_IntegrationStep.f90
524 lines (409 loc) · 27.3 KB
/
HMC_Module_Phys_Convolution_Apps_IntegrationStep.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
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
!------------------------------------------------------------------------------------
! File: HMC_Module_Phys_Convolution_Apps_IntegrationStep.f90
!
! Author(s): Fabio Delogu, Francesco Silvestro, Simone Gabellani
! Date: 20190410
!
! Convolution Apps IntegrationStep subroutine(s) for HMC model
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------
! Module Header
module HMC_Module_Phys_Convolution_Apps_IntegrationStep
!------------------------------------------------------------------------------------
! 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 to calculate dynamic integration step channel network
subroutine HMC_Phys_Convolution_Apps_IntegrationStep_ChannelNetwork(iID, iRows, iCols, dDtDataForcing, dDtIntegrAct)
!------------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iID, iRows, iCols
real(kind = 4) :: dDtDataForcing
real(kind = 4) :: dDtIntegr, dDtIntegrPStep, dDtPhysConv
real(kind = 4), intent(out) :: dDtIntegrAct
integer(kind = 4) :: iFlagVarDtPhysConv
real(kind = 4) :: dBc
real(kind = 4) :: dDtRef, dDtRefRatio, dDtIntegrMin
real(kind = 4) :: dDEMStepMean, dVarRainMax, dVarUDtMax
real(kind = 4) :: dVarHydroMax
real(kind = 4), dimension (iRows, iCols) :: a2dVarUc, a2dVarUh
real(kind = 4), dimension (iRows, iCols) :: a2dVarHydro, a2dVarRain, a2dVarRouting
real(kind = 4), dimension (iRows, iCols) :: a2dVarUcAct
real(kind = 4), dimension (iRows, iCols) :: a2dVarUDt
character(len = 10) :: sDtIntegrStep, sVarUDtMax, sVarRainMax
character(len = 10), parameter :: sFMTDtIntegrStep = "(I4)"
character(len = 10), parameter :: sFMTVarUDtMax = "(F8.5)"
character(len = 10), parameter :: sFMTVarRainMax = "(F8.5)"
integer(kind = 4) :: iDtPhysMethod
real(kind = 4), dimension(4) :: a1dDemStep, a1dIntegStep, a1dDtRef, a1dDtRefRatio
real(kind = 4) :: dDemDelta
integer(kind = 4) :: iIndexEnd, iIndexStart
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Variable(s) initialization
dDtIntegrAct = 0.0; dDtIntegrPStep = 0.0; dDtIntegr = 0.0;
dVarRainMax = 0.0; dDEMStepMean = 0.0;
dDtRef = 0.0; dDtRefRatio = 0.0; dDtIntegrMin = 0.0;
a2dVarHydro = 0.0;
a2dVarUc = 0.0; a2dVarUh = 0.0;
a2dVarUcAct = 0.0;
a2dVarUDt = 0.0;
iDtPhysMethod = 0;
a1dDemStep = -9999.0; a1dIntegStep = -9999.0; a1dDtRef = -9999.0; a1dDtRefRatio = -9999.0;
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Temporal integration at previous step
dDtIntegr = real(oHMC_Vars(iID)%iDtIntegr)
dDtIntegrPStep = real(oHMC_Vars(iID)%iDtIntegrPStep)
! Dynamic temporal integration step
iFlagVarDtPhysConv = oHMC_Namelist(iID)%iFlagVarDtPhysConv
dDtPhysConv = real(oHMC_Namelist(iID)%iDtPhysConv)
iDtPhysMethod = oHMC_Namelist(iID)%iDtPhysMethod
a1dDemStep = oHMC_Namelist(iID)%a1dDemStep
a1dIntegStep = oHMC_Namelist(iID)%a1dIntStep
a1dDtRef = oHMC_Namelist(iID)%a1dDtStep
a1dDtRefRatio = oHMC_Namelist(iID)%a1dDtRatioStep
! Exponent of dUcAct formula
dBc = oHMC_Namelist(iID)%dBc
! DEM mean step
dDEMStepMean = oHMC_Vars(iID)%dDEMStepMean
! Channel and hill surface velocity
a2dVarUc = oHMC_Vars(iID)%a2dUc
a2dVarUh = oHMC_Vars(iID)%a2dUh
! Dynamic variable(s)
a2dVarHydro = oHMC_Vars(iID)%a2dHydro
a2dVarRain = oHMC_Vars(iID)%a2dRain
a2dVarRouting = oHMC_Vars(iID)%a2dRouting
a2dVarUcAct = oHMC_Vars(iID)%a2dUcAct
a2dVarUDt = oHMC_Vars(iID)%a2dUDt
! Info start
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: IntegrationStep ... ' )
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Dynamic temporal integration step flag
if (iFlagVarDtPhysConv.eq.1 ) then
!------------------------------------------------------------------------------------------
! Actual temporal integration Step
dDtIntegrAct = dDtIntegr
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Definition of temporal integration parameter(s)
if (iDtPhysMethod .eq. 1) then
!------------------------------------------------------------------------------------------
! Get integration parameter(s) using scalar method
dDtRef = a1dDtRef(3); ! seconds
dDtIntegrMin = a1dIntegStep(3); ! seconds
dDtRefRatio = a1dDtRefRatio(3) ! seconds
!if ( (dDEMStepMean.ge.10) .and. (dDEMStepMean.lt.150) ) then ! NB: for vda integration step was too little
if ( (dDEMStepMean .ge. a1dDemStep(2)) .and. (dDEMStepMean .lt. a1dDemStep(3)) ) then
dDtRef = a1dDtRef(2); ! seconds
dDtIntegrMin = a1dIntegStep(2); ! seconds
dDtRefRatio = a1dDtRefRatio(2) ! seconds
elseif ( dDEMStepMean .lt. a1dDemStep(2) ) then
dDtRef = a1dDtRef(1); ! seconds
dDtIntegrMin = a1dIntegStep(1); ! seconds
dDtRefRatio = a1dDtRefRatio(1); ! seconds
elseif (dDEMStepMean .gt. a1dDemStep(4) ) then
dDtRef = a1dDtRef(4); ! seconds
dDtIntegrMin = a1dIntegStep(4); ! seconds
dDtRefRatio = a1dDtRefRatio(4); ! seconds
else
dDtRef = a1dDtRef(3); ! seconds
dDtIntegrMin = a1dIntegStep(3); ! seconds
dDtRefRatio = a1dDtRefRatio(3); ! seconds
endif
!------------------------------------------------------------------------------------------
elseif (iDtPhysMethod .eq. 2) then
!------------------------------------------------------------------------------------------
! Get integration parameter(s) using linear method
call getIntRange(a1dDemStep, a1dIntegStep, dDEMStepMean, dDemDelta, iIndexStart, iIndexEnd)
call getIntValue(a1dDemStep, a1dIntegStep, a1dDtRef, a1dDtRefRatio, &
dDEMStepMean, dDemDelta, iIndexStart, iIndexEnd, &
dDtIntegrMin, dDtRef, dDtRefRatio)
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Dynamic integration step evaluation
where( (oHMC_Vars(iID)%a2iChoice.eq.1) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarUcAct = a2dVarUc*(tan(oHMC_Vars(iID)%a2dBeta)**0.5)*a2dVarHydro**dBc
a2dVarUDt = a2dVarHydro*a2dVarUcAct/(1000*3600)*oHMC_Vars(iID)%a2dAreaCell !m^3/s
endwhere
where( (oHMC_Vars(iID)%a2iChoice.eq.0) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarUDt = a2dVarHydro*a2dVarUh/(1000*3600)*oHMC_Vars(iID)%a2dAreaCell !m^3/s
endwhere
! Checking waterlevel and updating
where( (a2dVarHydro.gt.0.0) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) )
a2dVarUDt = a2dVarUDt/(a2dVarHydro/1000*dDEMStepMean) !m/s
elsewhere
a2dVarUDt = 0.0
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Rain intensity max value
dVarRainMax = maxval(maxval(a2dVarRain,DIM = 1),DIM = 1)*3600.0/dDtDataForcing
! Velocity max value
dVarUDtMax = maxval(maxval(a2dVarUDt,DIM = 1),DIM = 1)
if (dVarUDtMax.le.0.1) dVarUDtMax = 0.1
! Hydro max value
dVarHydroMax = maxval(maxval(a2dVarHydro,DIM = 1),DIM = 1)
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! First estimation of integration step
dDtIntegrAct = dDEMStepMean/dVarUDtMax*0.6
! Checking Integration step value
if (dDtIntegrAct.gt.dDtDataForcing/dDtRefRatio) then
dDtIntegrAct = dDtDataForcing/dDtRefRatio
endif
if ( (dVarRainMax.gt.1.0) .and. (dDtIntegrAct.gt.dDtPhysConv) ) then
dDtIntegrAct = dDtPhysConv*exp(-dVarRainMax/3.0) + dDtPhysConv
endif
if (dDtIntegrAct.lt.dDtIntegrMin) then
dDtIntegrAct = dDtIntegrMin
endif
dDtIntegrAct = dDtDataForcing/int(dDtDataForcing/dDtIntegrAct)
dDtIntegrAct = int(dDtIntegrAct/dDtRef)*dDtRef
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Checking if dt integration is too big of previous dt integration
if(dDtIntegrAct.gt.dDtIntegrPStep + dDtRef) then
dDtIntegrAct = dDtIntegrPStep + dDtRef
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Updating routing variable; intensity must be the same between steps
a2dVarRouting = a2dVarRouting/dDtIntegrPStep*dDtIntegrAct
!------------------------------------------------------------------------------------------
else
!------------------------------------------------------------------------------------------
! Actual temporal integration Step
dDtIntegrAct = dDtIntegr
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Temporal integration step info
write(sDtIntegrStep, sFMTDtIntegrStep) int(dDtIntegrAct)
write(sVarUDtMax, sFMTVarUDtMax) dVarUDtMax
write(sVarRainMax, sFMTVarRainMax) dVarRainMax
call mprintf(.true., iINFO_Basic, ' Phys :: Convolution :: IntegrationStep :: '//sDtIntegrStep//' [s]')
call mprintf(.true., iINFO_Basic, ' Phys :: Convolution :: MaxValue :: '// &
' Uc: '//sVarUDtMax//' [m/s] '// &
' Rain: '//sVarRainMax//' [mm]')
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Passing variable(s) to global declaration
oHMC_Vars(iID)%iDtIntegr = int(dDtIntegrAct)
oHMC_Vars(iID)%iDtIntegrPStep = int(dDtIntegrAct)
oHMC_Vars(iID)%a2dRouting = a2dVarRouting
oHMC_Vars(iID)%a2dUcAct = a2dVarUcAct
oHMC_Vars(iID)%a2dUDt = a2dVarUDt
! Info end
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: IntegrationStep ... OK' )
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Convolution_Apps_IntegrationStep_ChannelNetwork
!------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Subroutine to calculate dynamic integration step channel fraction
subroutine HMC_Phys_Convolution_Apps_IntegrationStep_ChannelFraction(iID, iRows, iCols, dDtDataForcing, dDtIntegrAct)
!------------------------------------------------------------------------------------------
! Variable(s) declaration
integer(kind = 4) :: iID, iRows, iCols
real(kind = 4) :: dDtDataForcing
real(kind = 4) :: dDtIntegr, dDtIntegrPStep, dDtPhysConv
real(kind = 4), intent(out) :: dDtIntegrAct
integer(kind = 4) :: iFlagVarDtPhysConv
real(kind = 4) :: dBc
real(kind = 4) :: dDtRef, dDtRefRatio, dDtIntegrMin
real(kind = 4) :: dDEMStepMean, dVarRainMax, dVarUDtMax
real(kind = 4) :: dVarHydroMax
real(kind = 4), dimension (iRows, iCols) :: a2dVarUc, a2dVarUh
real(kind = 4), dimension (iRows, iCols) :: a2dVarHydro, a2dVarRain, a2dVarRouting
real(kind = 4), dimension (iRows, iCols) :: a2dVarUcAct
real(kind = 4), dimension (iRows, iCols) :: a2dVarUDt
character(len = 10) :: sDtIntegrStep, sVarUDtMax, sVarRainMax
character(len = 10), parameter :: sFMTDtIntegrStep = "(I4)"
character(len = 10), parameter :: sFMTVarUDtMax = "(F8.5)"
character(len = 10), parameter :: sFMTVarRainMax = "(F8.5)"
integer(kind = 4) :: iDtPhysMethod
real(kind = 4), dimension(4) :: a1dDemStep, a1dIntegStep, a1dDtRef, a1dDtRefRatio
real(kind = 4) :: dDemDelta
integer(kind = 4) :: iIndexEnd, iIndexStart
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Variable(s) initialization
dDtIntegrAct = 0.0; dDtIntegrPStep = 0.0; dDtIntegr = 0.0;
dVarRainMax = 0.0; dDEMStepMean = 0.0;
dDtRef = 0.0; dDtRefRatio = 0.0; dDtIntegrMin = 0.0;
a2dVarHydro = 0.0;
a2dVarUc = 0.0; a2dVarUh = 0.0;
a2dVarUcAct = 0.0;
a2dVarUDt = 0.0;
iDtPhysMethod = 0;
a1dDemStep = -9999.0; a1dIntegStep = -9999.0; a1dDtRef = -9999.0; a1dDtRefRatio = -9999.0;
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Temporal integration at previous step
dDtIntegr = real(oHMC_Vars(iID)%iDtIntegr)
dDtIntegrPStep = real(oHMC_Vars(iID)%iDtIntegrPStep)
! Dynamic temporal integration step
iFlagVarDtPhysConv = oHMC_Namelist(iID)%iFlagVarDtPhysConv
dDtPhysConv = real(oHMC_Namelist(iID)%iDtPhysConv)
iDtPhysMethod = oHMC_Namelist(iID)%iDtPhysMethod
a1dDemStep = oHMC_Namelist(iID)%a1dDemStep
a1dIntegStep = oHMC_Namelist(iID)%a1dIntStep
a1dDtRef = oHMC_Namelist(iID)%a1dDtStep
a1dDtRefRatio = oHMC_Namelist(iID)%a1dDtRatioStep
! Exponent of dUcAct formula
dBc = oHMC_Namelist(iID)%dBc
! DEM mean step
dDEMStepMean = oHMC_Vars(iID)%dDEMStepMean
! Channel and hill surface velocity
a2dVarUc = oHMC_Vars(iID)%a2dUc
a2dVarUh = oHMC_Vars(iID)%a2dUh
! Dynamic variable(s)
a2dVarRain = oHMC_Vars(iID)%a2dRain
a2dVarRouting = oHMC_Vars(iID)%a2dRouting
a2dVarUcAct = oHMC_Vars(iID)%a2dUcAct
a2dVarUDt = oHMC_Vars(iID)%a2dUDt
! Info start
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: IntegrationStep ... ' )
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Dynamic temporal integration step flag
if (iFlagVarDtPhysConv.eq.1 ) then
!------------------------------------------------------------------------------------------
! Actual temporal integration Step
dDtIntegrAct = dDtIntegr
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Definition of temporal integration parameter(s)
if (iDtPhysMethod .eq. 1) then
!------------------------------------------------------------------------------------------
! Get integration parameter(s) using scalar method
dDtRef = a1dDtRef(3); ! seconds
dDtIntegrMin = a1dIntegStep(3); ! seconds
dDtRefRatio = a1dDtRefRatio(3) ! seconds
!if ( (dDEMStepMean.ge.10) .and. (dDEMStepMean.lt.150) ) then ! NB: for vda integration step was too little
if ( (dDEMStepMean .ge. a1dDemStep(2)) .and. (dDEMStepMean .lt. a1dDemStep(3)) ) then
dDtRef = a1dDtRef(2); ! seconds
dDtIntegrMin = a1dIntegStep(2); ! seconds
dDtRefRatio = a1dDtRefRatio(2) ! seconds
elseif ( dDEMStepMean .lt. a1dDemStep(2) ) then
dDtRef = a1dDtRef(1); ! seconds
dDtIntegrMin = a1dIntegStep(1); ! seconds
dDtRefRatio = a1dDtRefRatio(1); ! seconds
elseif (dDEMStepMean .gt. a1dDemStep(4) ) then
dDtRef = a1dDtRef(4); ! seconds
dDtIntegrMin = a1dIntegStep(4); ! seconds
dDtRefRatio = a1dDtRefRatio(4); ! seconds
else
dDtRef = a1dDtRef(3); ! seconds
dDtIntegrMin = a1dIntegStep(3); ! seconds
dDtRefRatio = a1dDtRefRatio(3); ! seconds
endif
!------------------------------------------------------------------------------------------
elseif (iDtPhysMethod .eq. 2) then
!------------------------------------------------------------------------------------------
! Get integration parameter(s) using linear method
call getIntRange(a1dDemStep, a1dIntegStep, dDEMStepMean, dDemDelta, iIndexStart, iIndexEnd)
call getIntValue(a1dDemStep, a1dIntegStep, a1dDtRef, a1dDtRefRatio, &
dDEMStepMean, dDemDelta, iIndexStart, iIndexEnd, &
dDtIntegrMin, dDtRef, dDtRefRatio)
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Dynamic integration step evaluation
where( (oHMC_Vars(iID)%a2iChoice.eq.1) .and. (oHMC_Vars(iID)%a2iMask.gt.0.0) .and. (oHMC_Vars(iID)%a2dQC.gt.0.0) )
a2dVarUcAct = oHMC_Vars(iID)%a2dQC/(oHMC_Vars(iID)%a2dWidthC*oHMC_Vars(iID)%a2dHydroC) !m/s
a2dVarUDt = oHMC_Vars(iID)%a2dQC/(oHMC_Vars(iID)%a2dWidthC*oHMC_Vars(iID)%a2dHydroC) !m3/s
elsewhere
a2dVarUDt = 0.0
endwhere
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Rain intensity max value
dVarRainMax = maxval(maxval(a2dVarRain,DIM = 1),DIM = 1)*3600.0/dDtDataForcing
! Velocity max value
dVarUDtMax = maxval(maxval(a2dVarUDt,DIM = 1),DIM = 1)
if (dVarUDtMax.le.0.1) dVarUDtMax = 0.1
! Hydro max value
dVarHydroMax = maxval(maxval(oHMC_Vars(iID)%a2dHydroC,DIM = 1),DIM = 1)
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! First estimation of integration step
dDtIntegrAct = dDEMStepMean/dVarUDtMax*0.6
! Checking Integration step value
if (dDtIntegrAct.gt.dDtDataForcing/dDtRefRatio) then
dDtIntegrAct = dDtDataForcing/dDtRefRatio
endif
if ( (dVarRainMax.gt.1.0) .and. (dDtIntegrAct.gt.dDtPhysConv) ) then
dDtIntegrAct = dDtPhysConv*exp(-dVarRainMax/3.0) + dDtPhysConv
endif
if (dDtIntegrAct.lt.dDtIntegrMin) then
dDtIntegrAct = dDtIntegrMin
endif
dDtIntegrAct = dDtDataForcing/int(dDtDataForcing/dDtIntegrAct)
dDtIntegrAct = int(dDtIntegrAct/dDtRef)*dDtRef
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Checking if dt integration is too big of previous dt integration
if(dDtIntegrAct.gt.dDtIntegrPStep + dDtRef) then
dDtIntegrAct = dDtIntegrPStep + dDtRef
endif
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Updating routing variable; intensity must be the same between steps
a2dVarRouting = a2dVarRouting/dDtIntegrPStep*dDtIntegrAct
!------------------------------------------------------------------------------------------
else
!------------------------------------------------------------------------------------------
! Actual temporal integration Step
dDtIntegrAct = dDtIntegr
!------------------------------------------------------------------------------------------
endif
!------------------------------------------------------------------------------------------
! dDtIntegrAct = 30 !!!DA LEVARE
!------------------------------------------------------------------------------------------
! Temporal integration step info
write(sDtIntegrStep, sFMTDtIntegrStep) int(dDtIntegrAct)
write(sVarUDtMax, sFMTVarUDtMax) dVarUDtMax
write(sVarRainMax, sFMTVarRainMax) dVarRainMax
call mprintf(.true., iINFO_Basic, ' Phys :: Convolution :: IntegrationStep :: '//sDtIntegrStep//' [s]')
call mprintf(.true., iINFO_Basic, ' Phys :: Convolution :: MaxValue :: '// &
' Uc: '//sVarUDtMax//' [m/s] '// &
' Rain: '//sVarRainMax//' [mm]')
!------------------------------------------------------------------------------------------
!------------------------------------------------------------------------------------------
! Passing variable(s) to global declaration
oHMC_Vars(iID)%iDtIntegr = int(dDtIntegrAct)
oHMC_Vars(iID)%iDtIntegrPStep = int(dDtIntegrAct)
oHMC_Vars(iID)%a2dRouting = a2dVarRouting
oHMC_Vars(iID)%a2dUcAct = a2dVarUcAct
oHMC_Vars(iID)%a2dUDt = a2dVarUDt
! Info end
if (iDEBUG.gt.0) then
call mprintf(.true., iINFO_Extra, ' Phys :: Convolution :: IntegrationStep ... OK' )
endif
!------------------------------------------------------------------------------------------
end subroutine HMC_Phys_Convolution_Apps_IntegrationStep_ChannelFraction
!------------------------------------------------------------------------------------
end module HMC_Module_Phys_Convolution_Apps_IntegrationStep
!------------------------------------------------------------------------------------------