-
Notifications
You must be signed in to change notification settings - Fork 2
/
wetdry_draft.f90
776 lines (728 loc) · 38.4 KB
/
wetdry_draft.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
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
module wetdry
real(sz) :: habsmin
real(sz) :: hoff
real(sz) :: h0 ! read h0 from fort.15
integer :: mnp ! number of nodes
integer :: mne ! number of cells
integer, allocatable :: nibcnt(:)
integer, allocatable :: noffold(:)
subroutine initializeWettingAndDrying(h0)
use sizes, only: mne, mnp
use global, only: h0, nodecode, nnodecode
implicit none
allocate (nibcnt(mnp))
allocate (noffold(mne))
nnodecode = 1
nodecode = 1
noffold(:) = 1
habsmin = 0.8d0*h0
hoff = 1.2d0*h0
!----------------------------------------------------------------------
end subroutine initializeWettingAndDrying
!----------------------------------------------------------------------
!----------------------------------------------------------------------
! S U B R O U T I N E
! C O M P U T E W E T T I N G A N D D R Y I N G
!----------------------------------------------------------------------
! Determines which nodes should be wet and which should be dry.
!----------------------------------------------------------------------
subroutine computeWettingAndDrying(it)
use sizes, only: sz, mne, mnp
use global, only: noff, nodecode, nnodecode, eta2, tk, nolifa, &
bsx1, bsy1, btime_end, C2DDI, C3D, g, h0, ifnlfa, nddt, &
nibnodecode, ilump, ncchange &
use mesh, only: ne, np, dp, mju, totalArea, nm, x, y, areas
use nodalattributes, only: BFCdLLimit, fgamma, ftheta, fric, &
manningsn, hbreak, ifhybf, ifnlbf, iflinbf, loadManningsN, &
loadZ0B_var, z0b_var
use global_3dvs, only: a, b, islip, kp, z0b, sigma, evtot, q
use subdomain, only: subdomainOn, enforceBN, enforceWDcb, enforceWDob
implicit none
integer, intent(in) :: it ! time step number
complex(sz) :: duds !jgf48.50 declare size SZ instead of plain COMPLEX
integer :: nc1, nc2, nc3
integer :: nm1, nm2, nm3
integer :: nm123
integer :: ncele
integer :: nctot
real(sz) :: kslip
real(sz) :: vel
real(sz) :: z0b1
real(sz) :: areaEle
real(sz) :: tkWet
real(sz) :: etaN1, etaN2, etaN3
real(sz) :: hTotN1, hTotN2, hTotN3
real(sz) :: deldist, deleta
real(sz) :: htot
real(sz) :: h1
integer :: nbnctot
integer :: i
integer :: ie
! WET...
! WET...WET/DRY - INITIALIZATIONS FOR WET/DRY LOOP
! WET...
DO I = 1, NP
NIBCNT(I) = 0
ENDDO
DO I = 1, NE
NOFFOLD(I) = NOFF(I)
NOFF(I) = 1
ENDDO
! WET...
! WET...WET/DRY - PART 1 - NODAL DRYING CRITERIA D1
! WET....Drying Criteria D1: this depends on NODECODE and updates NODECODE
! WET...
DO I = 1, NP
IF (NODECODE(I) .EQ. 1) THEN
HTOT = DP(I) + ETA2(I)
IF (HTOT .LE. H0) THEN
IF (HTOT .LT. HABSMIN) ETA2(I) = HABSMIN - DP(I)
NNODECODE(I) = 0
NODECODE(I) = 0
NCCHANGE = NCCHANGE + 1 !NCCHANGE=0 set near beginning of GWCE
! ENDIF
ENDIF
ENDIF
ENDDO
! WET...
! WET...END WET/DRY SECTION - PART 1
! WET...
! WET...
! WET...WET/DRY SECTION PART 2 - NODAL WETTING LOOPS W1 AND W2
! WET...
DO I = 1, NE
NM1 = NM(I, 1)
NM2 = NM(I, 2)
NM3 = NM(I, 3)
! WET...
! WET...Nodal Wetting Criteria W1: This depends on changes that occurred in D1
! WET...
NCTOT = NODECODE(NM1) + NODECODE(NM2) + NODECODE(NM3)
IF (NCTOT .EQ. 2) THEN
ETAN1 = ETA2(NM1)
ETAN2 = ETA2(NM2)
ETAN3 = ETA2(NM3)
HTOTN1 = DP(NM1) + ETA2(NM1)
HTOTN2 = DP(NM2) + ETA2(NM2)
HTOTN3 = DP(NM3) + ETA2(NM3)
IF ((NODECODE(NM1) .EQ. 1) .AND. (NODECODE(NM2) .EQ. 1)) THEN
IF ((HTOTN1 .GE. HOFF) .AND. (HTOTN2 .GE. HOFF)) THEN
NM123 = NM1
IF (ETA2(NM2) .GT. ETA2(NM1)) NM123 = NM2
DELDIST = SQRT((y(NM3) - y(NM123))**2.D0 &
+ (X(NM3) - X(NM123))**2.D0)
DELETA = ETA2(NM123) - ETA2(NM3)
! jgf50.60.18: Prevent numerical problems if DELETA is negative
IF (DELETA .lt. 0.d0) DELETA = 0.d0
H1 = ETA2(NM123) + DP(NM123)
! RJW merged from Casey 071219: Added the following logic for 3D friction.
! RJW modified the following for 3D friction
IF (C2DDI) THEN
! sb46.28sb02
!<< Convert Manning's N to Cd, if necessary.
IF (LoadManningsN) THEN
FRIC(NM123) = g*ManningsN(NM123)**2.d0 &
/((DP(NM123) + IFNLFA*ETA2(NM123)) &
**(1.d0/3.d0))
IF (FRIC(NM123) .LT. BFCdLLimit) THEN
FRIC(NM123) = BFCdLLimit
ENDIF
ENDIF
!>>
TKWET = FRIC(NM123)*(IFLINBF + (VELMIN/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)**(FGAMMA/FTHETA)))
IF (TKWET .LT. 0.0001d0) TKWET = 0.0001d0
VEL = G*(DELETA/DELDIST)/TKWET
ELSEIF (C3D) THEN
!C solve for the depth averaged velocity,U, from the relation :
!C tau=rho*g*(h+eta)*(deta/dx)=rho*Cd*|U|*U
!C U=sqrt(g*(h+eta)*(deta/dx)/Cd )
!C where: Cd=kappa^2/(ln(z+zo)/z0)^2 is the depth integrated drag coefficient
IF (LoadZ0B_var) THEN
Z0B1 = Z0B_var(NM123)
ELSEIF (LoadManningsN) THEN
Z0B1 = (DP(NM123) + IFNLFA*ETA2(NM123))*exp(-(1.0D0 + &
((0.41D0*(DP(NM123) + IFNLFA*ETA2(NM123))**(1.0D0/6.0D0))/ &
(ManningsN(NM123)*sqrt(g)))))
ELSE
Z0B1 = Z0B
ENDIF
VEL = sqrt(g*H1*(DELETA/DELDIST)) &
*((H1 + Z0B1)*LOG((H1 + Z0B1)/Z0B1) - H1)/(H1*0.41D0)
ENDIF
IF (VEL .GT. VELMIN) THEN
! .... third node met criteria and is also wet
NNODECODE(NM3) = 1
! RJW merged 08/26/20008 Casey 071219: Added the following logic to obtain the correct friction.
IF (C2DDI) THEN
! TK(NM123)=FRIC(NM123)*(IFLINBF+(VEL/H1)*
TK(NM123) = FRIC(NM123)*(IFLINBF + (VEL/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)** &
(FGAMMA/FTHETA)))
ELSEIF (C3D) THEN
IF (ISLIP .EQ. 0) THEN
DUDS = (Q(NM123, 2) - Q(NM123, 1)) &
/(SIGMA(2) - SIGMA(1))
BSX1(NM123) = EVTOT(1)*REAL(DUDS)
BSY1(NM123) = EVTOT(1)*AIMAG(DUDS)
BSX1(NM3) = EVTOT(1)*REAL(DUDS)
BSY1(NM3) = EVTOT(1)*AIMAG(DUDS)
ENDIF
IF (ISLIP .NE. 0) THEN
IF (ISLIP .EQ. 1) THEN
KSLIP = KP
ENDIF
IF (ISLIP .EQ. 2) THEN
KSLIP = (1.D0/ &
((1.D0/0.41D0)* &
LOG((ABS(((SIGMA(2) - SIGMA(1))/(A - B))* &
(DP(NM123) + IFNLFA*ETA2(NM123))) &
+ Z0B1) &
/(Z0B1))))**2.D0 &
*ABS(Q(NM123, 1))
IF (KSLIP .GT. 1.D0*ABS(Q(NM123, 1))) &
KSLIP = 1.D0*ABS(Q(NM123, 1))
IF (KSLIP .LT. 0.0025D0*ABS(Q(NM123, 1))) &
KSLIP = 0.0025D0*ABS(Q(NM123, 1))
ENDIF
BSX1(NM123) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM123) = KSLIP*AIMAG(Q(NM123, 1))
BSX1(NM3) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM3) = KSLIP*AIMAG(Q(NM123, 1))
ENDIF
ENDIF
ENDIF
ENDIF
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
ELSEIF ((NODECODE(NM2) .EQ. 1) .AND. (NODECODE(NM3) .EQ. 1)) &
THEN
IF ((HTOTN2 .GE. HOFF) .AND. (HTOTN3 .GE. HOFF)) THEN
NM123 = NM2
IF (ETA2(NM3) .GT. ETA2(NM2)) NM123 = NM3
DELDIST = SQRT((Y(NM1) - Y(NM123))**2.D0 &
+ (X(NM1) - X(NM123))**2.D0)
DELETA = ETA2(NM123) - ETA2(NM1)
! jgf50.60.18: Prevent numerical problems if DELETA is negative
IF (DELETA .lt. 0.d0) DELETA = 0.d0
H1 = ETA2(NM123) + DP(NM123)
! RJW merged 08/26/2008 Casey 071219: Added the following logic for 3D friction.
IF (C2DDI) THEN
! sb46.28sb02
!C<< Convert Manning's N to Cd, if necessary.
IF (LoadManningsN) THEN
FRIC(NM123) = g*ManningsN(NM123)**2.d0 &
/((DP(NM123) + IFNLFA*ETA2(NM123)) &
**(1.d0/3.d0))
IF (FRIC(NM123) .LT. BFCdLLimit) THEN
FRIC(NM123) = BFCdLLimit
ENDIF
ENDIF
!C>>
TKWET = FRIC(NM123)*(IFLINBF + (VELMIN/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)**(FGAMMA/FTHETA)))
IF (TKWET .LT. 0.0001d0) TKWET = 0.0001d0
VEL = G*(DELETA/DELDIST)/TKWET
ELSEIF (C3D) THEN
! solve for the depth averaged velocity,U, from the relation :
! tau=rho*g*(h+eta)*(deta/dx)=rho*Cd*|U|*U
! U=sqrt(g*(h+eta)*(deta/dx)/Cd )
! where: Cd=kappa^2/(ln(z+zo)/z0)^2 is the depth integrated drag coefficient
IF (LoadZ0B_var) THEN
Z0B1 = Z0B_var(NM123)
ELSEIF (LoadManningsN) THEN
Z0B1 = (DP(NM123) + IFNLFA*ETA2(NM123))*exp(-(1.0D0 + &
((0.41D0*(DP(NM123) + IFNLFA*ETA2(NM123))**(1.0D0/6.0D0))/ &
(ManningsN(NM123)*sqrt(g)))))
ELSE
Z0B1 = Z0B
ENDIF
VEL = sqrt(g*H1*(DELETA/DELDIST)) &
*((H1 + Z0B1)*LOG((H1 + Z0B1)/Z0B1) - H1)/(H1*0.41D0)
ENDIF
IF (VEL .GT. VELMIN) THEN
NNODECODE(NM1) = 1
! RJW merged 08/26/2008 Casey 071219: Added the following logic to obtain the correct friction.
IF (C2DDI) THEN
TK(NM123) = FRIC(NM123)*(IFLINBF + (VEL/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)** &
(FGAMMA/FTHETA)))
ELSEIF (C3D) THEN
IF (ISLIP .EQ. 0) THEN
DUDS = (Q(NM123, 2) - Q(NM123, 1)) &
/(SIGMA(2) - SIGMA(1))
BSX1(NM123) = EVTOT(1)*REAL(DUDS)
BSY1(NM123) = EVTOT(1)*AIMAG(DUDS)
BSX1(NM1) = EVTOT(1)*REAL(DUDS)
BSY1(NM1) = EVTOT(1)*AIMAG(DUDS)
ENDIF
IF (ISLIP .NE. 0) THEN
IF (ISLIP .EQ. 1) THEN
KSLIP = KP
ENDIF
IF (ISLIP .EQ. 2) THEN
KSLIP = (1.D0/ &
((1.D0/0.41D0)* &
LOG((ABS(((SIGMA(2) - SIGMA(1))/(A - B))* &
(DP(NM123) + IFNLFA*ETA2(NM123))) &
+ Z0B1) &
/(Z0B1))))**2.D0 &
*ABS(Q(NM123, 1))
IF (KSLIP .GT. 1.D0*ABS(Q(NM123, 1))) &
KSLIP = 1.D0*ABS(Q(NM123, 1))
IF (KSLIP .LT. 0.0025D0*ABS(Q(NM123, 1))) &
KSLIP = 0.0025D0*ABS(Q(NM123, 1))
ENDIF
BSX1(NM123) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM123) = KSLIP*AIMAG(Q(NM123, 1))
BSX1(NM1) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM1) = KSLIP*AIMAG(Q(NM123, 1))
ENDIF
ENDIF
ENDIF
ENDIF
! !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
ELSEIF ((NODECODE(NM3) .EQ. 1) .AND. (NODECODE(NM1) .EQ. 1)) &
THEN
IF ((HTOTN3 .GE. HOFF) .AND. (HTOTN1 .GE. HOFF)) THEN
NM123 = NM3
IF (ETA2(NM1) .GT. ETA2(NM3)) NM123 = NM1
DELDIST = SQRT((Y(NM2) - Y(NM123))**2.D0 &
+ (X(NM2) - X(NM123))**2.D0)
DELETA = ETA2(NM123) - ETA2(NM2)
! jgf50.60.18: Prevent numerical problems if DELETA is negative
IF (DELETA .lt. 0.d0) DELETA = 0.d0
H1 = ETA2(NM123) + DP(NM123)
! RJW merged 08/26/2008 Casey 071219: Added the following logic for 3D friction.
IF (C2DDI) THEN
! sb46.28sb02
!C<< Convert Manning's N to Cd, if necessary.
IF (LoadManningsN) THEN
FRIC(NM123) = g*ManningsN(NM123)**2.d0 &
/((DP(NM123) + IFNLFA*ETA2(NM123)) &
**(1.d0/3.d0))
IF (FRIC(NM123) .LT. BFCdLLimit) THEN
FRIC(NM123) = BFCdLLimit
ENDIF
ENDIF
!C>>
TKWET = FRIC(NM123)*(IFLINBF + (VELMIN/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)**(FGAMMA/FTHETA)))
IF (TKWET .LT. 0.0001d0) TKWET = 0.0001d0
VEL = G*(DELETA/DELDIST)/TKWET
ELSEIF (C3D) THEN
! solve for the depth averaged velocity,U, from the relation :
! tau=rho*g*(h+eta)*(deta/dx)=rho*Cd*|U|*U
! U=sqrt(g*(h+eta)*(deta/dx)/Cd )
! where: Cd=kappa^2/(ln(z+zo)/z0)^2 is the depth integrated drag coefficient
IF (LoadZ0B_var) THEN
Z0B1 = Z0B_var(NM123)
ELSEIF (LoadManningsN) THEN
Z0B1 = (DP(NM123) + IFNLFA*ETA2(NM123))*exp(-(1.0D0 + &
((0.41D0*(DP(NM123) + IFNLFA*ETA2(NM123))**(1.0D0/6.0D0))/ &
(ManningsN(NM123)*sqrt(g)))))
ELSE
Z0B1 = Z0B
ENDIF
VEL = sqrt(g*H1*(DELETA/DELDIST)) &
*((H1 + Z0B1)*LOG((H1 + Z0B1)/Z0B1) - H1)/(H1*0.41D0)
ENDIF
IF (VEL .GT. VELMIN) THEN
NNODECODE(NM2) = 1
! RJW merged 08/26/2008 Casey 071219: Added the following logic to obtain the correct friction.
IF (C2DDI) THEN
TK(NM123) = FRIC(NM123)*(IFLINBF + (VEL/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)** &
(FGAMMA/FTHETA)))
ELSEIF (C3D) THEN
IF (ISLIP .EQ. 0) THEN
DUDS = (Q(NM123, 2) - Q(NM123, 1)) &
/(SIGMA(2) - SIGMA(1))
BSX1(NM123) = EVTOT(1)*REAL(DUDS)
BSY1(NM123) = EVTOT(1)*AIMAG(DUDS)
BSX1(NM2) = EVTOT(1)*REAL(DUDS)
BSY1(NM2) = EVTOT(1)*AIMAG(DUDS)
ENDIF
IF (ISLIP .NE. 0) THEN
IF (ISLIP .EQ. 1) THEN
KSLIP = KP
ENDIF
IF (ISLIP .EQ. 2) THEN
KSLIP = (1.D0/ &
((1.D0/0.41D0)* &
LOG((ABS(((SIGMA(2) - SIGMA(1))/(A - B))* &
(DP(NM123) + IFNLFA*ETA2(NM123))) &
+ Z0B1) &
/(Z0B1))))**2.D0 &
*ABS(Q(NM123, 1))
IF (KSLIP .GT. 1.D0*ABS(Q(NM123, 1))) &
KSLIP = 1.D0*ABS(Q(NM123, 1))
IF (KSLIP .LT. 0.0025D0*ABS(Q(NM123, 1))) &
KSLIP = 0.0025D0*ABS(Q(NM123, 1))
ENDIF
BSX1(NM123) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM123) = KSLIP*AIMAG(Q(NM123, 1))
BSX1(NM2) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM2) = KSLIP*AIMAG(Q(NM123, 1))
ENDIF
ENDIF
ENDIF
ENDIF
ENDIF
ENDIF
! WET...
! WET...Nodal Wetting Criteria W2a
! WET...
NBNCTOT = NIBNODECODE(NM1) + NIBNODECODE(NM2) + NIBNODECODE(NM3)
NIBCNT(NM1) = NIBCNT(NM1) + NBNCTOT
NIBCNT(NM2) = NIBCNT(NM2) + NBNCTOT
NIBCNT(NM3) = NIBCNT(NM3) + NBNCTOT
ENDDO
if (subdomainOn .and. enforceBN .eq. 1) call enforceWDcb() ! NCSU Subdomain
if (subdomainOn .and. enforceBN .eq. 2) call enforceWDob() ! NCSU Subdomain
! wet...
! WET... ELEMENTAL WETTING CRITERIA WETBATHYCHANGE
! *******************************************************************************************
! tcm v50.66.01 -- This is an additional test for wetting only when time varying
! bathymetry is used and is only performed during the period of
! bathymetry evolution.
!
IF ((NDDT .NE. 0) .AND. (IT .LE. BTIME_END + 1)) THEN
DO I = 1, NE
NM1 = NM(I, 1)
NM2 = NM(I, 2)
NM3 = NM(I, 3)
NCTOT = NODECODE(NM1) + NODECODE(NM2) + NODECODE(NM3)
IF (NCTOT .lt. 3) THEN !If not wet from previous time step
NCTOT = NNODECODE(NM1) + NNODECODE(NM2) + NNODECODE(NM3)
if (NCTOT .lt. 3) then !if not alreay made wet for this time step
ETAN1 = ETA2(NM1)
ETAN2 = ETA2(NM2)
ETAN3 = ETA2(NM3)
HTOTN1 = DP(NM1) + ETA2(NM1)
HTOTN2 = DP(NM2) + ETA2(NM2)
HTOTN3 = DP(NM3) + ETA2(NM3)
! if all nodes have a depth greater than or equal to
! hoff = 1.2*H0, then make the element wet
IF ((HTOTN1 .GE. HOFF) .AND. (HTOTN2 .GE. HOFF) .AND. &
(HTOTN3 .GE. HOFF)) THEN
!THE ELEMENT SHOULD BE WET, SO WET THE DRY NODES
! Make Node 1 Wet and set parameters
IF (NNODECODE(NM1) .NE. 1) THEN !node 1
NNODECODE(NM1) = 1
NM123 = NM1
IF (C2DDI) THEN
!C<< Convert Manning's N to Cd, if necessary.
IF (LoadManningsN) THEN
FRIC(NM123) = g*ManningsN(NM123)**2.d0 &
/((DP(NM123) + IFNLFA*ETA2(NM123)) &
**(1.d0/3.d0))
IF (FRIC(NM123) .LT. BFCdLLimit) THEN
FRIC(NM123) = BFCdLLimit
ENDIF
ENDIF
ENDIF
VEL = VELMIN
H1 = HTOTN1
IF (C2DDI) THEN
TK(NM123) = FRIC(NM123)*(IFLINBF + (VEL/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)** &
(FGAMMA/FTHETA)))
ELSEIF (C3D) THEN
IF (ISLIP .EQ. 0) THEN
DUDS = (Q(NM123, 2) - Q(NM123, 1)) &
/(SIGMA(2) - SIGMA(1))
BSX1(NM123) = EVTOT(1)*REAL(DUDS)
BSY1(NM123) = EVTOT(1)*AIMAG(DUDS)
BSX1(NM3) = EVTOT(1)*REAL(DUDS)
BSY1(NM3) = EVTOT(1)*AIMAG(DUDS)
ENDIF
IF (ISLIP .NE. 0) THEN
IF (ISLIP .EQ. 1) THEN
KSLIP = KP
ENDIF
IF (ISLIP .EQ. 2) THEN
KSLIP = (1.D0/ &
((1.D0/0.41D0)* &
LOG((ABS(((SIGMA(2) - SIGMA(1))/(A - B))* &
(DP(NM123) + IFNLFA*ETA2(NM123))) &
+ Z0B1) &
/(Z0B1))))**2.D0 &
*ABS(Q(NM123, 1))
IF (KSLIP .GT. 1.D0*ABS(Q(NM123, 1))) &
KSLIP = 1.D0*ABS(Q(NM123, 1))
IF (KSLIP .LT. 0.0025D0*ABS(Q(NM123, 1))) &
KSLIP = 0.0025D0*ABS(Q(NM123, 1))
ENDIF
BSX1(NM123) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM123) = KSLIP*AIMAG(Q(NM123, 1))
BSX1(NM3) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM3) = KSLIP*AIMAG(Q(NM123, 1))
ENDIF
ENDIF
ENDIF !end node 1
! Make Node 2 Wet and set parameters
IF (NNODECODE(NM2) .NE. 1) THEN
NNODECODE(NM2) = 1
NM123 = NM2
IF (C2DDI) THEN
!C<< Convert Manning's N to Cd, if necessary.
IF (LoadManningsN) THEN
FRIC(NM123) = g*ManningsN(NM123)**2.d0 &
/((DP(NM123) + IFNLFA*ETA2(NM123)) &
**(1.d0/3.d0))
IF (FRIC(NM123) .LT. BFCdLLimit) THEN
FRIC(NM123) = BFCdLLimit
ENDIF
ENDIF
ENDIF
VEL = VELMIN
H1 = HTOTN2
IF (C2DDI) THEN
TK(NM123) = FRIC(NM123)*(IFLINBF + (VEL/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)** &
(FGAMMA/FTHETA)))
ELSEIF (C3D) THEN
IF (ISLIP .EQ. 0) THEN
DUDS = (Q(NM123, 2) - Q(NM123, 1)) &
/(SIGMA(2) - SIGMA(1))
BSX1(NM123) = EVTOT(1)*REAL(DUDS)
BSY1(NM123) = EVTOT(1)*AIMAG(DUDS)
BSX1(NM1) = EVTOT(1)*REAL(DUDS)
BSY1(NM1) = EVTOT(1)*AIMAG(DUDS)
ENDIF
IF (ISLIP .NE. 0) THEN
IF (ISLIP .EQ. 1) THEN
KSLIP = KP
ENDIF
IF (ISLIP .EQ. 2) THEN
KSLIP = (1.D0/ &
((1.D0/0.41D0)* &
LOG((ABS(((SIGMA(2) - SIGMA(1))/(A - B))* &
(DP(NM123) + IFNLFA*ETA2(NM123))) &
+ Z0B1) &
/(Z0B1))))**2.D0 &
*ABS(Q(NM123, 1))
IF (KSLIP .GT. 1.D0*ABS(Q(NM123, 1))) &
KSLIP = 1.D0*ABS(Q(NM123, 1))
IF (KSLIP .LT. 0.0025D0*ABS(Q(NM123, 1))) &
KSLIP = 0.0025D0*ABS(Q(NM123, 1))
ENDIF
BSX1(NM123) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM123) = KSLIP*AIMAG(Q(NM123, 1))
BSX1(NM1) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM1) = KSLIP*AIMAG(Q(NM123, 1))
ENDIF
ENDIF
ENDIF !node 2
! Make Node 3 Wet and set parameters
IF (NNODECODE(NM3) .NE. 1) THEN
NNODECODE(NM3) = 1
NM123 = NM3
IF (C2DDI) THEN
!C<< Convert Manning's N to Cd, if necessary.
IF (LoadManningsN) THEN
FRIC(NM123) = g*ManningsN(NM123)**2.d0 &
/((DP(NM123) + IFNLFA*ETA2(NM123)) &
**(1.d0/3.d0))
IF (FRIC(NM123) .LT. BFCdLLimit) THEN
FRIC(NM123) = BFCdLLimit
ENDIF
ENDIF
ENDIF
VEL = VELMIN
H1 = HTOTN3
IF (C2DDI) THEN
TK(NM123) = FRIC(NM123)*(IFLINBF + (VEL/H1)* &
(IFNLBF + IFHYBF* &
(1.D0 + (HBREAK/H1)**FTHETA)** &
(FGAMMA/FTHETA)))
ELSEIF (C3D) THEN
IF (ISLIP .EQ. 0) THEN
DUDS = (Q(NM123, 2) - Q(NM123, 1)) &
/(SIGMA(2) - SIGMA(1))
BSX1(NM123) = EVTOT(1)*REAL(DUDS)
BSY1(NM123) = EVTOT(1)*AIMAG(DUDS)
BSX1(NM2) = EVTOT(1)*REAL(DUDS)
BSY1(NM2) = EVTOT(1)*AIMAG(DUDS)
ENDIF
IF (ISLIP .NE. 0) THEN
IF (ISLIP .EQ. 1) THEN
KSLIP = KP
ENDIF
IF (ISLIP .EQ. 2) THEN
KSLIP = (1.D0/ &
((1.D0/0.41D0)* &
LOG((ABS(((SIGMA(2) - SIGMA(1))/(A - B))* &
(DP(NM123) + IFNLFA*ETA2(NM123))) &
+ Z0B1) &
/(Z0B1))))**2.D0 &
*ABS(Q(NM123, 1))
IF (KSLIP .GT. 1.D0*ABS(Q(NM123, 1))) &
KSLIP = 1.D0*ABS(Q(NM123, 1))
IF (KSLIP .LT. 0.0025D0*ABS(Q(NM123, 1))) &
KSLIP = 0.0025D0*ABS(Q(NM123, 1))
ENDIF
BSX1(NM123) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM123) = KSLIP*AIMAG(Q(NM123, 1))
BSX1(NM2) = KSLIP*REAL(Q(NM123, 1))
BSY1(NM2) = KSLIP*AIMAG(Q(NM123, 1))
ENDIF
ENDIF
ENDIF !node3
ENDIF !ALL DEPTHS GREATER THAN HOFF
ENDIF !IF NNODECODE SUM LESS THAN 3
ENDIF ! IF NODECODE SUM LESS THAN 3
ENDDO !LOOP OVER ELEMENTS
if (subdomainOn .and. enforceBN .eq. 1) call enforceWDcb() ! NCSU Subdomain
if (subdomainOn .and. enforceBN .eq. 2) call enforceWDob() ! NCSU Subdomain
ENDIF !IT TIME VARYING BATHYMETRY AND WITHIN CHANGE TIME
! ... END OF ADDITIONAL WETTING FOR TIME VARYING BATHYMETRY
! WET..
! WET... ELEMENTAL WETTING CRITERIA WETBATHYCHANGE
! *******************************************************************************************
! WET...
! WET...Nodal Wetting Criteria W2b
! WET...Check for adjacent nodes and force nodes wet when attached
! WET...to receiving barrier nodes
! WET...
DO I = 1, NP
IF ((NIBCNT(I) .GT. 0) .AND. (NNODECODE(I) .EQ. 0)) THEN
NNODECODE(I) = 1
ENDIF
ENDDO
! WET...
! WET...END WET/DRY SECTION - PART 2
! WET...
! WET...
! WET...START WET/DRY SECTION - PART 3
! WET...Elemental drying criteria DE1
! WET...This is an elemental check section designed to avoid artificial wetting of
! WET....of control sections
! WET...All elements where downhill flow originates from a barely wet node
! WET....into wet nodes are forced inactive; the only exception is receiving
! WET....overtopped barrier nodes
! WET...
DO I = 1, NE
NM1 = NM(I, 1)
NM2 = NM(I, 2)
NM3 = NM(I, 3)
NBNCTOT = NIBCNT(NM1)*NIBCNT(NM2)*NIBCNT(NM3)
IF (NBNCTOT .EQ. 0) THEN !No barrier/pipe receiving nodes in this elem
ETAN1 = ETA2(NM1)
ETAN2 = ETA2(NM2)
ETAN3 = ETA2(NM3)
HTOTN1 = DP(NM1) + ETA2(NM1)
HTOTN2 = DP(NM2) + ETA2(NM2)
HTOTN3 = DP(NM3) + ETA2(NM3)
! jgf52.08.08: Analyst can eliminate noff from
! consideration in fort.15 file.
if (noffActive .eqv. .true.) then
! ..ABC pattern
IF ((ETAN1 .GE. ETAN2) .AND. (ETAN2 .GT. ETAN3)) THEN
IF ((HTOTN1 .LT. HOFF) .OR. (HTOTN2 .LT. HOFF)) NOFF(I) = 0
ENDIF
IF ((ETAN2 .GE. ETAN3) .AND. (ETAN3 .GT. ETAN1)) THEN
IF ((HTOTN2 .LT. HOFF) .OR. (HTOTN3 .LT. HOFF)) NOFF(I) = 0
ENDIF
IF ((ETAN3 .GE. ETAN1) .AND. (ETAN1 .GT. ETAN2)) THEN
IF ((HTOTN3 .LT. HOFF) .OR. (HTOTN1 .LT. HOFF)) NOFF(I) = 0
ENDIF
! ..ACB pattern
IF ((ETAN1 .GE. ETAN3) .AND. (ETAN3 .GT. ETAN2)) THEN
IF ((HTOTN1 .LT. HOFF) .OR. (HTOTN3 .LT. HOFF)) NOFF(I) = 0
ENDIF
IF ((ETAN2 .GE. ETAN1) .AND. (ETAN1 .GT. ETAN3)) THEN
IF ((HTOTN2 .LT. HOFF) .OR. (HTOTN1 .LT. HOFF)) NOFF(I) = 0
ENDIF
IF ((ETAN3 .GE. ETAN2) .AND. (ETAN2 .GT. ETAN1)) THEN
IF ((HTOTN3 .LT. HOFF) .OR. (HTOTN2 .LT. HOFF)) NOFF(I) = 0
ENDIF
endif
ENDIF
ENDDO
! WET...
! WET...END WET/DRY SECTION - PART 3
! WET...
! WET...
! WET...START WET/DRY SECTION PART 4 - NODAL DRYING LOOP D2
! WET...Update number of active elements (MJU) and the total area (TotalArea) connected
! WET...to a node. If these are zero, the node is landlocked and should be dried.
! WET...These depend on NNODECODE which varies during the time step
! WET...
DO I = 1, NP
MJU(I) = 0
TotalArea(I) = 0.d0
ENDDO
DO IE = 1, NE
NM1 = NM(IE, 1)
NM2 = NM(IE, 2)
NM3 = NM(IE, 3)
NC1 = NNODECODE(NM1)
NC2 = NNODECODE(NM2)
NC3 = NNODECODE(NM3)
NCEle = NC1*NC2*NC3*NOFF(IE)
AreaEle = NCEle*Areas(IE)/2.d0
MJU(NM1) = MJU(NM1) + NCEle
MJU(NM2) = MJU(NM2) + NCEle
MJU(NM3) = MJU(NM3) + NCEle
TotalArea(NM1) = TotalArea(NM1) + AreaEle
TotalArea(NM2) = TotalArea(NM2) + AreaEle
TotalArea(NM3) = TotalArea(NM3) + AreaEle
ENDDO
! jjwnote - looks like this is used later in momentum equations
! jjwnote - this has implications on making this into a subroutine
DO I = 1, NP
IF ((NNODECODE(I) .EQ. 1) .AND. (MJU(I) .EQ. 0)) THEN
NNODECODE(I) = 0
ENDIF
IF (MJU(I) .EQ. 0) MJU(I) = 1 !Because MJU is also used to solve Mom Eq. !Eliminate this?
ENDDO
! WET...
! WET...END WET/DRY SECTION - PART 4
! WET...
! jjwnote - may have to pass TotalArea and mju as well
if (subdomainOn .and. enforceBN .eq. 1) call enforceWDcb() ! NCSU Subdomain
if (subdomainOn .and. enforceBN .eq. 2) call enforceWDob() ! NCSU Subdomain
! WET...
! WET...WET/DRY SECTION - PART 5 - RESET NODECODE USING NNODECODE
! WET...Check to see if any wetting occurred & update NODECODE
! WET...Note, NCCHANGE=0 set near the beginning of GWCE subroutine
! WET...
DO I = 1, NP
IF (NNODECODE(I) .NE. NODECODE(I)) THEN
NODECODE(I) = NNODECODE(I)
NCCHANGE = NCCHANGE + 1
ENDIF
ENDDO
! WET...
! WET...END WET/DRY SECTION - PART 5
! WET...
! WET...
! WET...WET/DRY SECTION - PART 6
! WET...Check to see if any NOFF changed requiring the matrix to be reset
! WET...Note, NCCHANGE=0 set near the beginning of GWCE subroutine
! WET...
DO I = 1, NE
IF (NOFF(I) .NE. NOFFOLD(I)) NCCHANGE = NCCHANGE + 1
ENDDO
! WET...
! WET... jgf45.06 If there has been any wetting or drying in any
! WET... of the subdomains, the NCCHANGE flag will be activated on all
! WET... of the subdomains, to prevent them from getting out of sync
! WET... with their MPI calls as some reset the GWCE and others do not.
! WET...
! WET...END WET/DRY SECTION - PART 6
! WET...
! ....
!-----------------------------------------------------------------------
end subroutine computeWettingAndDrying
!-----------------------------------------------------------------------
end module wetDry