/
velo.f90
3571 lines (3006 loc) · 117 KB
/
velo.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
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
MODULE VELO
! Module computes the velocity flux terms, baroclinic torque correction terms, and performs the CFL Check
USE PRECISION_PARAMETERS
USE GLOBAL_CONSTANTS
USE MESH_POINTERS
USE COMP_FUNCTIONS, ONLY: SECOND
IMPLICIT NONE
PRIVATE
PUBLIC COMPUTE_VELOCITY_FLUX,VELOCITY_PREDICTOR,VELOCITY_CORRECTOR,NO_FLUX,BAROCLINIC_CORRECTION, &
MATCH_VELOCITY,MATCH_VELOCITY_FLUX,VELOCITY_BC,COMPUTE_VISCOSITY,VISCOSITY_BC
PRIVATE VELOCITY_FLUX,VELOCITY_FLUX_CYLINDRICAL
CONTAINS
SUBROUTINE COMPUTE_VELOCITY_FLUX(T,DT,NM,FUNCTION_CODE)
REAL(EB), INTENT(IN) :: T,DT
REAL(EB) :: TNOW
INTEGER, INTENT(IN) :: NM,FUNCTION_CODE
IF (SOLID_PHASE_ONLY .OR. FREEZE_VELOCITY) RETURN
TNOW = SECOND()
SELECT CASE(FUNCTION_CODE)
CASE(1)
IF (PREDICTOR .OR. COMPUTE_VISCOSITY_TWICE) CALL COMPUTE_VISCOSITY(T,NM)
CASE(2)
BAROCLINIC_TERMS_ATTACHED = .FALSE.
IF (PREDICTOR .OR. COMPUTE_VISCOSITY_TWICE) CALL VISCOSITY_BC(NM)
IF (.NOT.CYLINDRICAL) CALL VELOCITY_FLUX(T,DT,NM)
IF ( CYLINDRICAL) CALL VELOCITY_FLUX_CYLINDRICAL(T,NM)
END SELECT
T_USED(4) = T_USED(4) + SECOND() - TNOW
END SUBROUTINE COMPUTE_VELOCITY_FLUX
SUBROUTINE COMPUTE_VISCOSITY(T,NM)
USE PHYSICAL_FUNCTIONS, ONLY: GET_VISCOSITY,LES_FILTER_WIDTH_FUNCTION,GET_POTENTIAL_TEMPERATURE
USE TURBULENCE, ONLY: VARDEN_DYNSMAG,TEST_FILTER,FILL_EDGES,WALL_MODEL,RNG_EDDY_VISCOSITY,WALE_VISCOSITY
USE MATH_FUNCTIONS, ONLY:EVALUATE_RAMP
REAL(EB), INTENT(IN) :: T
INTEGER, INTENT(IN) :: NM
REAL(EB) :: ZZ_GET(1:N_TRACKED_SPECIES),NU_EDDY,DELTA,KSGS,U2,V2,W2,AA,A_IJ(3,3),BB,B_IJ(3,3),&
DUDX,DUDY,DUDZ,DVDX,DVDY,DVDZ,DWDX,DWDY,DWDZ,MU_EFF,SLIP_COEF,VEL_GAS,VEL_T,RAMP_T,TSI,&
VDF,LS,THETA_0,THETA_1,THETA_2,DTDZBAR,WGT
REAL(EB), PARAMETER :: RAPLUS=1._EB/26._EB, C_LS=0.76_EB
INTEGER :: I,J,K,IIG,JJG,KKG,II,JJ,KK,IW,IOR
REAL(EB), POINTER, DIMENSION(:,:,:) :: RHOP=>NULL(),UP=>NULL(),VP=>NULL(),WP=>NULL(), &
UP_HAT=>NULL(),VP_HAT=>NULL(),WP_HAT=>NULL(), &
UU=>NULL(),VV=>NULL(),WW=>NULL(),DTDZ=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:,:) :: ZZP=>NULL()
INTEGER, POINTER, DIMENSION(:,:,:) :: CELL_COUNTER=>NULL()
TYPE(WALL_TYPE), POINTER :: WC=>NULL()
TYPE(SURFACE_TYPE), POINTER :: SF=>NULL()
IF (EVACUATION_ONLY(NM)) RETURN ! No need to update viscosity, use initial one
CALL POINT_TO_MESH(NM)
IF (PREDICTOR) THEN
RHOP => RHO
UU => U
VV => V
WW => W
ZZP => ZZ
ELSE
RHOP => RHOS
UU => US
VV => VS
WW => WS
ZZP => ZZS
ENDIF
! Compute viscosity for DNS using primitive species
!$OMP PARALLEL DO FIRSTPRIVATE(ZZ_GET) SCHEDULE(guided)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
ZZ_GET(1:N_TRACKED_SPECIES) = ZZP(I,J,K,1:N_TRACKED_SPECIES)
CALL GET_VISCOSITY(ZZ_GET,MU_DNS(I,J,K),TMP(I,J,K))
ENDDO
ENDDO
ENDDO
!$OMP END PARALLEL DO
CALL COMPUTE_STRAIN_RATE(NM)
SELECT_TURB: SELECT CASE (TURB_MODEL)
CASE (NO_TURB_MODEL)
MU = MU_DNS
CASE (CONSMAG,DYNSMAG) SELECT_TURB ! Smagorinsky (1963) eddy viscosity
IF (PREDICTOR .AND. TURB_MODEL==DYNSMAG) CALL VARDEN_DYNSMAG(NM) ! dynamic procedure, Moin et al. (1991)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
MU(I,J,K) = MU_DNS(I,J,K) + RHOP(I,J,K)*CSD2(I,J,K)*STRAIN_RATE(I,J,K)
ENDDO
ENDDO
ENDDO
CASE (DEARDORFF) SELECT_TURB ! Deardorff (1980) eddy viscosity model (current default)
! Velocities relative to the p-cell center
UP => WORK1
VP => WORK2
WP => WORK3
UP=0._EB
VP=0._EB
WP=0._EB
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
UP(I,J,K) = 0.5_EB*(UU(I,J,K) + UU(I-1,J,K))
VP(I,J,K) = 0.5_EB*(VV(I,J,K) + VV(I,J-1,K))
WP(I,J,K) = 0.5_EB*(WW(I,J,K) + WW(I,J,K-1))
ENDDO
ENDDO
ENDDO
! fill mesh boundary ghost cells
DO IW=1,N_EXTERNAL_WALL_CELLS
WC=>WALL(IW)
SELECT CASE(WC%BOUNDARY_TYPE)
CASE(INTERPOLATED_BOUNDARY)
II = WC%ONE_D%II
JJ = WC%ONE_D%JJ
KK = WC%ONE_D%KK
UP(II,JJ,KK) = U_GHOST(IW)
VP(II,JJ,KK) = V_GHOST(IW)
WP(II,JJ,KK) = W_GHOST(IW)
CASE(OPEN_BOUNDARY,MIRROR_BOUNDARY)
II = WC%ONE_D%II
JJ = WC%ONE_D%JJ
KK = WC%ONE_D%KK
IIG = WC%ONE_D%IIG
JJG = WC%ONE_D%JJG
KKG = WC%ONE_D%KKG
UP(II,JJ,KK) = UP(IIG,JJG,KKG)
VP(II,JJ,KK) = VP(IIG,JJG,KKG)
WP(II,JJ,KK) = WP(IIG,JJG,KKG)
END SELECT
ENDDO
! fill edge and corner ghost cells
CALL FILL_EDGES(UP)
CALL FILL_EDGES(VP)
CALL FILL_EDGES(WP)
UP_HAT => WORK4
VP_HAT => WORK5
WP_HAT => WORK6
UP_HAT=0._EB
VP_HAT=0._EB
WP_HAT=0._EB
CALL TEST_FILTER(UP_HAT,UP)
CALL TEST_FILTER(VP_HAT,VP)
CALL TEST_FILTER(WP_HAT,WP)
POTENTIAL_TEMPERATURE_IF: IF (.NOT.POTENTIAL_TEMPERATURE_CORRECTION) THEN
!$OMP PARALLEL DO PRIVATE(DELTA, KSGS, NU_EDDY) SCHEDULE(static)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
DELTA = LES_FILTER_WIDTH_FUNCTION(DX(I),DY(J),DZ(K))
KSGS = 0.5_EB*( (UP(I,J,K)-UP_HAT(I,J,K))**2 + (VP(I,J,K)-VP_HAT(I,J,K))**2 + (WP(I,J,K)-WP_HAT(I,J,K))**2 )
NU_EDDY = C_DEARDORFF*DELTA*SQRT(KSGS)
MU(I,J,K) = MU_DNS(I,J,K) + RHOP(I,J,K)*NU_EDDY
ENDDO
ENDDO
ENDDO
!$OMP END PARALLEL DO
ELSE POTENTIAL_TEMPERATURE_IF
DTDZ => WORK7
DTDZ = 0._EB
DO K=0,KBAR
DO J=0,JBAR
DO I=0,IBAR
THETA_1 = GET_POTENTIAL_TEMPERATURE(TMP(I,J,K),ZC(K))
THETA_2 = GET_POTENTIAL_TEMPERATURE(TMP(I,J,K+1),ZC(K+1))
DTDZ(I,J,K) = (THETA_2-THETA_1)*RDZN(K)
ENDDO
ENDDO
ENDDO
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
DELTA = LES_FILTER_WIDTH_FUNCTION(DX(I),DY(J),DZ(K))
LS = DELTA
KSGS = 0.5_EB*( (UP(I,J,K)-UP_HAT(I,J,K))**2 + (VP(I,J,K)-VP_HAT(I,J,K))**2 + (WP(I,J,K)-WP_HAT(I,J,K))**2 )
DTDZBAR = 0.5_EB*(DTDZ(I,J,K)+DTDZ(I,J,K+1))
IF (DTDZBAR>0._EB) THEN
THETA_0 = GET_POTENTIAL_TEMPERATURE(TMP_0(K),ZC(K))
LS = C_LS*SQRT(KSGS)/SQRT(ABS(GVEC(3))/THETA_0*DTDZBAR) ! von Schoenberg Eq. (3.19)
ENDIF
NU_EDDY = C_DEARDORFF*MIN(LS,DELTA)*SQRT(KSGS)
MU(I,J,K) = MU_DNS(I,J,K) + RHOP(I,J,K)*NU_EDDY
ENDDO
ENDDO
ENDDO
ENDIF POTENTIAL_TEMPERATURE_IF
CASE (VREMAN) SELECT_TURB ! Vreman (2004) eddy viscosity model (experimental)
! A. W. Vreman. An eddy-viscosity subgrid-scale model for turbulent shear flow: Algebraic theory and applications.
! Phys. Fluids, 16(10):3670-3681, 2004.
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
DUDX = RDX(I)*(UU(I,J,K)-UU(I-1,J,K))
DVDY = RDY(J)*(VV(I,J,K)-VV(I,J-1,K))
DWDZ = RDZ(K)*(WW(I,J,K)-WW(I,J,K-1))
DUDY = 0.25_EB*RDY(J)*(UU(I,J+1,K)-UU(I,J-1,K)+UU(I-1,J+1,K)-UU(I-1,J-1,K))
DUDZ = 0.25_EB*RDZ(K)*(UU(I,J,K+1)-UU(I,J,K-1)+UU(I-1,J,K+1)-UU(I-1,J,K-1))
DVDX = 0.25_EB*RDX(I)*(VV(I+1,J,K)-VV(I-1,J,K)+VV(I+1,J-1,K)-VV(I-1,J-1,K))
DVDZ = 0.25_EB*RDZ(K)*(VV(I,J,K+1)-VV(I,J,K-1)+VV(I,J-1,K+1)-VV(I,J-1,K-1))
DWDX = 0.25_EB*RDX(I)*(WW(I+1,J,K)-WW(I-1,J,K)+WW(I+1,J,K-1)-WW(I-1,J,K-1))
DWDY = 0.25_EB*RDY(J)*(WW(I,J+1,K)-WW(I,J-1,K)+WW(I,J+1,K-1)-WW(I,J-1,K-1))
! Vreman, Eq. (6)
A_IJ(1,1)=DUDX; A_IJ(2,1)=DUDY; A_IJ(3,1)=DUDZ
A_IJ(1,2)=DVDX; A_IJ(2,2)=DVDY; A_IJ(3,2)=DVDZ
A_IJ(1,3)=DWDX; A_IJ(2,3)=DWDY; A_IJ(3,3)=DWDZ
AA=0._EB
DO JJ=1,3
DO II=1,3
AA = AA + A_IJ(II,JJ)*A_IJ(II,JJ)
ENDDO
ENDDO
! Vreman, Eq. (7)
B_IJ(1,1)=(DX(I)*A_IJ(1,1))**2 + (DY(J)*A_IJ(2,1))**2 + (DZ(K)*A_IJ(3,1))**2
B_IJ(2,2)=(DX(I)*A_IJ(1,2))**2 + (DY(J)*A_IJ(2,2))**2 + (DZ(K)*A_IJ(3,2))**2
B_IJ(3,3)=(DX(I)*A_IJ(1,3))**2 + (DY(J)*A_IJ(2,3))**2 + (DZ(K)*A_IJ(3,3))**2
B_IJ(1,2)=DX(I)**2*A_IJ(1,1)*A_IJ(1,2) + DY(J)**2*A_IJ(2,1)*A_IJ(2,2) + DZ(K)**2*A_IJ(3,1)*A_IJ(3,2)
B_IJ(1,3)=DX(I)**2*A_IJ(1,1)*A_IJ(1,3) + DY(J)**2*A_IJ(2,1)*A_IJ(2,3) + DZ(K)**2*A_IJ(3,1)*A_IJ(3,3)
B_IJ(2,3)=DX(I)**2*A_IJ(1,2)*A_IJ(1,3) + DY(J)**2*A_IJ(2,2)*A_IJ(2,3) + DZ(K)**2*A_IJ(3,2)*A_IJ(3,3)
BB = B_IJ(1,1)*B_IJ(2,2) - B_IJ(1,2)**2 &
+ B_IJ(1,1)*B_IJ(3,3) - B_IJ(1,3)**2 &
+ B_IJ(2,2)*B_IJ(3,3) - B_IJ(2,3)**2 ! Vreman, Eq. (8)
IF (ABS(AA)>TWO_EPSILON_EB .AND. BB>TWO_EPSILON_EB) THEN
NU_EDDY = C_VREMAN*SQRT(BB/AA) ! Vreman, Eq. (5)
ELSE
NU_EDDY=0._EB
ENDIF
MU(I,J,K) = MU_DNS(I,J,K) + RHOP(I,J,K)*NU_EDDY
ENDDO
ENDDO
ENDDO
CASE (RNG) SELECT_TURB
! A. Yakhot, S. A. Orszag, V. Yakhot, and M. Israeli. Renormalization Group Formulation of Large-Eddy Simulation.
! Journal of Scientific Computing, 1(1):1-51, 1989.
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
DELTA = LES_FILTER_WIDTH_FUNCTION(DX(I),DY(J),DZ(K))
CALL RNG_EDDY_VISCOSITY(MU_EFF,MU_DNS(I,J,K),RHOP(I,J,K),STRAIN_RATE(I,J,K),DELTA)
MU(I,J,K) = MU_EFF
ENDDO
ENDDO
ENDDO
CASE (WALE) SELECT_TURB
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
DELTA = LES_FILTER_WIDTH_FUNCTION(DX(I),DY(J),DZ(K))
! compute velocity gradient tensor
DUDX = RDX(I)*(UU(I,J,K)-UU(I-1,J,K))
DVDY = RDY(J)*(VV(I,J,K)-VV(I,J-1,K))
DWDZ = RDZ(K)*(WW(I,J,K)-WW(I,J,K-1))
DUDY = 0.25_EB*RDY(J)*(UU(I,J+1,K)-UU(I,J-1,K)+UU(I-1,J+1,K)-UU(I-1,J-1,K))
DUDZ = 0.25_EB*RDZ(K)*(UU(I,J,K+1)-UU(I,J,K-1)+UU(I-1,J,K+1)-UU(I-1,J,K-1))
DVDX = 0.25_EB*RDX(I)*(VV(I+1,J,K)-VV(I-1,J,K)+VV(I+1,J-1,K)-VV(I-1,J-1,K))
DVDZ = 0.25_EB*RDZ(K)*(VV(I,J,K+1)-VV(I,J,K-1)+VV(I,J-1,K+1)-VV(I,J-1,K-1))
DWDX = 0.25_EB*RDX(I)*(WW(I+1,J,K)-WW(I-1,J,K)+WW(I+1,J,K-1)-WW(I-1,J,K-1))
DWDY = 0.25_EB*RDY(J)*(WW(I,J+1,K)-WW(I,J-1,K)+WW(I,J+1,K-1)-WW(I,J-1,K-1))
A_IJ(1,1)=DUDX; A_IJ(1,2)=DUDY; A_IJ(1,3)=DUDZ
A_IJ(2,1)=DVDX; A_IJ(2,2)=DVDY; A_IJ(2,3)=DVDZ
A_IJ(3,1)=DWDX; A_IJ(3,2)=DWDY; A_IJ(3,3)=DWDZ
CALL WALE_VISCOSITY(NU_EDDY,A_IJ,STRAIN_RATE(I,J,K),DELTA)
MU(I,J,K) = MU_DNS(I,J,K) + RHOP(I,J,K)*NU_EDDY
ENDDO
ENDDO
ENDDO
END SELECT SELECT_TURB
! Compute resolved kinetic energy per unit mass
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
U2 = 0.25_EB*(UU(I-1,J,K)+UU(I,J,K))**2
V2 = 0.25_EB*(VV(I,J-1,K)+VV(I,J,K))**2
W2 = 0.25_EB*(WW(I,J,K-1)+WW(I,J,K))**2
KRES(I,J,K) = 0.5_EB*(U2+V2+W2)
ENDDO
ENDDO
ENDDO
! Mirror viscosity into solids and exterior boundary cells
CELL_COUNTER => IWORK1 ; CELL_COUNTER = 0
WALL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS
WC=>WALL(IW)
IF (WC%BOUNDARY_TYPE==NULL_BOUNDARY) CYCLE WALL_LOOP
II = WC%ONE_D%II
JJ = WC%ONE_D%JJ
KK = WC%ONE_D%KK
IOR = WC%ONE_D%IOR
IIG = WC%ONE_D%IIG
JJG = WC%ONE_D%JJG
KKG = WC%ONE_D%KKG
SF=>SURFACE(WC%SURF_INDEX)
SELECT CASE(WC%BOUNDARY_TYPE)
CASE(SOLID_BOUNDARY)
IF (ABS(SF%T_IGN-T_BEGIN)<=SPACING(SF%T_IGN) .AND. SF%RAMP_INDEX(TIME_VELO)>=1) THEN
TSI = T
ELSE
TSI = T-SF%T_IGN
ENDIF
RAMP_T = EVALUATE_RAMP(TSI,SF%TAU(TIME_VELO),SF%RAMP_INDEX(TIME_VELO))
VEL_T = RAMP_T*SQRT(SF%VEL_T(1)**2 + SF%VEL_T(2)**2)
SELECT CASE(ABS(IOR))
CASE(1)
VEL_GAS = SQRT( 0.25_EB*( (VV(IIG,JJG,KKG)+VV(IIG,JJG-1,KKG))**2 + (WW(IIG,JJG,KKG)+WW(IIG,JJG,KKG-1))**2 ) )
CASE(2)
VEL_GAS = SQRT( 0.25_EB*( (UU(IIG,JJG,KKG)+UU(IIG-1,JJG,KKG))**2 + (WW(IIG,JJG,KKG)+WW(IIG,JJG,KKG-1))**2 ) )
CASE(3)
VEL_GAS = SQRT( 0.25_EB*( (UU(IIG,JJG,KKG)+UU(IIG-1,JJG,KKG))**2 + (VV(IIG,JJG,KKG)+VV(IIG,JJG-1,KKG))**2 ) )
END SELECT
CALL WALL_MODEL(SLIP_COEF,WC%U_TAU,WC%Y_PLUS,VEL_GAS-VEL_T,&
MU_DNS(IIG,JJG,KKG)/RHO(IIG,JJG,KKG),1._EB/WC%ONE_D%RDN,SURFACE(WC%SURF_INDEX)%ROUGHNESS)
IF (LES) THEN
DELTA = LES_FILTER_WIDTH_FUNCTION(DX(IIG),DY(JJG),DZ(KKG))
SELECT CASE(NEAR_WALL_TURB_MODEL)
CASE DEFAULT ! Constant Smagorinsky with Van Driest damping
VDF = 1._EB-EXP(-WC%Y_PLUS*RAPLUS)
NU_EDDY = (VDF*C_SMAGORINSKY*DELTA)**2*STRAIN_RATE(IIG,JJG,KKG)
CASE(WALE)
! compute velocity gradient tensor
DUDX = RDX(IIG)*(UU(IIG,JJG,KKG)-UU(IIG-1,JJG,KKG))
DVDY = RDY(JJG)*(VV(IIG,JJG,KKG)-VV(IIG,JJG-1,KKG))
DWDZ = RDZ(KKG)*(WW(IIG,JJG,KKG)-WW(IIG,JJG,KKG-1))
DUDY = 0.25_EB*RDY(JJG)*(UU(IIG,JJG+1,KKG)-UU(IIG,JJG-1,KKG)+UU(IIG-1,JJG+1,KKG)-UU(IIG-1,JJG-1,KKG))
DUDZ = 0.25_EB*RDZ(KKG)*(UU(IIG,JJG,KKG+1)-UU(IIG,JJG,KKG-1)+UU(IIG-1,JJG,KKG+1)-UU(IIG-1,JJG,KKG-1))
DVDX = 0.25_EB*RDX(IIG)*(VV(IIG+1,JJG,KKG)-VV(IIG-1,JJG,KKG)+VV(IIG+1,JJG-1,KKG)-VV(IIG-1,JJG-1,KKG))
DVDZ = 0.25_EB*RDZ(KKG)*(VV(IIG,JJG,KKG+1)-VV(IIG,JJG,KKG-1)+VV(IIG,JJG-1,KKG+1)-VV(IIG,JJG-1,KKG-1))
DWDX = 0.25_EB*RDX(IIG)*(WW(IIG+1,JJG,KKG)-WW(IIG-1,JJG,KKG)+WW(IIG+1,JJG,KKG-1)-WW(IIG-1,JJG,KKG-1))
DWDY = 0.25_EB*RDY(JJG)*(WW(IIG,JJG+1,KKG)-WW(IIG,JJG-1,KKG)+WW(IIG,JJG+1,KKG-1)-WW(IIG,JJG-1,KKG-1))
A_IJ(1,1)=DUDX; A_IJ(1,2)=DUDY; A_IJ(1,3)=DUDZ
A_IJ(2,1)=DVDX; A_IJ(2,2)=DVDY; A_IJ(2,3)=DVDZ
A_IJ(3,1)=DWDX; A_IJ(3,2)=DWDY; A_IJ(3,3)=DWDZ
CALL WALE_VISCOSITY(NU_EDDY,A_IJ,STRAIN_RATE(IIG,JJG,KKG),DELTA)
END SELECT
IF (CELL_COUNTER(IIG,JJG,KKG)==0) MU(IIG,JJG,KKG) = 0._EB
CELL_COUNTER(IIG,JJG,KKG) = CELL_COUNTER(IIG,JJG,KKG) + 1
WGT = 1._EB/REAL(CELL_COUNTER(IIG,JJG,KKG),EB)
MU(IIG,JJG,KKG) = (1._EB-WGT)*MU(IIG,JJG,KKG) + WGT*(MU_DNS(IIG,JJG,KKG) + RHOP(IIG,JJG,KKG)*NU_EDDY)
ELSE
MU(IIG,JJG,KKG) = MU_DNS(IIG,JJG,KKG)
ENDIF
IF (SOLID(CELL_INDEX(II,JJ,KK))) MU(II,JJ,KK) = MU(IIG,JJG,KKG)
CASE(OPEN_BOUNDARY,MIRROR_BOUNDARY)
MU(II,JJ,KK) = MU(IIG,JJG,KKG)
KRES(II,JJ,KK) = KRES(IIG,JJG,KKG)
END SELECT
ENDDO WALL_LOOP
MU( 0,0:JBP1, 0) = MU( 1,0:JBP1,1)
MU(IBP1,0:JBP1, 0) = MU(IBAR,0:JBP1,1)
MU(IBP1,0:JBP1,KBP1) = MU(IBAR,0:JBP1,KBAR)
MU( 0,0:JBP1,KBP1) = MU( 1,0:JBP1,KBAR)
MU(0:IBP1, 0, 0) = MU(0:IBP1, 1,1)
MU(0:IBP1,JBP1,0) = MU(0:IBP1,JBAR,1)
MU(0:IBP1,JBP1,KBP1) = MU(0:IBP1,JBAR,KBAR)
MU(0:IBP1,0,KBP1) = MU(0:IBP1, 1,KBAR)
MU(0, 0,0:KBP1) = MU( 1, 1,0:KBP1)
MU(IBP1,0,0:KBP1) = MU(IBAR, 1,0:KBP1)
MU(IBP1,JBP1,0:KBP1) = MU(IBAR,JBAR,0:KBP1)
MU(0,JBP1,0:KBP1) = MU( 1,JBAR,0:KBP1)
END SUBROUTINE COMPUTE_VISCOSITY
SUBROUTINE COMPUTE_STRAIN_RATE(NM)
INTEGER, INTENT(IN) :: NM
REAL(EB) :: DUDX,DUDY,DUDZ,DVDX,DVDY,DVDZ,DWDX,DWDY,DWDZ,S11,S22,S33,S12,S13,S23,ONTHDIV
INTEGER :: I,J,K,IOR,IIG,JJG,KKG,IW,SURF_INDEX
REAL(EB), POINTER, DIMENSION(:,:,:) :: UU=>NULL(),VV=>NULL(),WW=>NULL()
TYPE(WALL_TYPE), POINTER :: WC=>NULL()
CALL POINT_TO_MESH(NM)
IF (PREDICTOR) THEN
UU => U
VV => V
WW => W
ELSE
UU => US
VV => VS
WW => WS
ENDIF
SELECT CASE (TURB_MODEL)
CASE DEFAULT
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
DUDX = RDX(I)*(UU(I,J,K)-UU(I-1,J,K))
DVDY = RDY(J)*(VV(I,J,K)-VV(I,J-1,K))
DWDZ = RDZ(K)*(WW(I,J,K)-WW(I,J,K-1))
DUDY = 0.25_EB*RDY(J)*(UU(I,J+1,K)-UU(I,J-1,K)+UU(I-1,J+1,K)-UU(I-1,J-1,K))
DUDZ = 0.25_EB*RDZ(K)*(UU(I,J,K+1)-UU(I,J,K-1)+UU(I-1,J,K+1)-UU(I-1,J,K-1))
DVDX = 0.25_EB*RDX(I)*(VV(I+1,J,K)-VV(I-1,J,K)+VV(I+1,J-1,K)-VV(I-1,J-1,K))
DVDZ = 0.25_EB*RDZ(K)*(VV(I,J,K+1)-VV(I,J,K-1)+VV(I,J-1,K+1)-VV(I,J-1,K-1))
DWDX = 0.25_EB*RDX(I)*(WW(I+1,J,K)-WW(I-1,J,K)+WW(I+1,J,K-1)-WW(I-1,J,K-1))
DWDY = 0.25_EB*RDY(J)*(WW(I,J+1,K)-WW(I,J-1,K)+WW(I,J+1,K-1)-WW(I,J-1,K-1))
ONTHDIV = ONTH*(DUDX+DVDY+DWDZ)
S11 = DUDX - ONTHDIV
S22 = DVDY - ONTHDIV
S33 = DWDZ - ONTHDIV
S12 = 0.5_EB*(DUDY+DVDX)
S13 = 0.5_EB*(DUDZ+DWDX)
S23 = 0.5_EB*(DVDZ+DWDY)
STRAIN_RATE(I,J,K) = SQRT(2._EB*(S11**2 + S22**2 + S33**2 + 2._EB*(S12**2 + S13**2 + S23**2)))
ENDDO
ENDDO
ENDDO
CASE (DEARDORFF)
! Here we omit the 3D loop, we only need the wall cell values of STRAIN_RATE
END SELECT
WALL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS
WC=>WALL(IW)
IF (WC%BOUNDARY_TYPE/=SOLID_BOUNDARY) CYCLE WALL_LOOP
SURF_INDEX = WC%SURF_INDEX
IIG = WC%ONE_D%IIG
JJG = WC%ONE_D%JJG
KKG = WC%ONE_D%KKG
IOR = WC%ONE_D%IOR
! Handle the case where OBST lives on an external boundary
IF (IW>N_EXTERNAL_WALL_CELLS) THEN
SELECT CASE(IOR)
CASE( 1); IF (IIG>IBAR) CYCLE WALL_LOOP
CASE(-1); IF (IIG<1) CYCLE WALL_LOOP
CASE( 2); IF (JJG>JBAR) CYCLE WALL_LOOP
CASE(-2); IF (JJG<1) CYCLE WALL_LOOP
CASE( 3); IF (KKG>KBAR) CYCLE WALL_LOOP
CASE(-3); IF (KKG<1) CYCLE WALL_LOOP
END SELECT
ENDIF
DUDX = RDX(IIG)*(UU(IIG,JJG,KKG)-UU(IIG-1,JJG,KKG))
DVDY = RDY(JJG)*(VV(IIG,JJG,KKG)-VV(IIG,JJG-1,KKG))
DWDZ = RDZ(KKG)*(WW(IIG,JJG,KKG)-WW(IIG,JJG,KKG-1))
ONTHDIV = ONTH*(DUDX+DVDY+DWDZ)
S11 = DUDX - ONTHDIV
S22 = DVDY - ONTHDIV
S33 = DWDZ - ONTHDIV
DUDY = 0.25_EB*RDY(JJG)*(UU(IIG,JJG+1,KKG)-UU(IIG,JJG-1,KKG)+UU(IIG-1,JJG+1,KKG)-UU(IIG-1,JJG-1,KKG))
DUDZ = 0.25_EB*RDZ(KKG)*(UU(IIG,JJG,KKG+1)-UU(IIG,JJG,KKG-1)+UU(IIG-1,JJG,KKG+1)-UU(IIG-1,JJG,KKG-1))
DVDX = 0.25_EB*RDX(IIG)*(VV(IIG+1,JJG,KKG)-VV(IIG-1,JJG,KKG)+VV(IIG+1,JJG-1,KKG)-VV(IIG-1,JJG-1,KKG))
DVDZ = 0.25_EB*RDZ(KKG)*(VV(IIG,JJG,KKG+1)-VV(IIG,JJG,KKG-1)+VV(IIG,JJG-1,KKG+1)-VV(IIG,JJG-1,KKG-1))
DWDX = 0.25_EB*RDX(IIG)*(WW(IIG+1,JJG,KKG)-WW(IIG-1,JJG,KKG)+WW(IIG+1,JJG,KKG-1)-WW(IIG-1,JJG,KKG-1))
DWDY = 0.25_EB*RDY(JJG)*(WW(IIG,JJG+1,KKG)-WW(IIG,JJG-1,KKG)+WW(IIG,JJG+1,KKG-1)-WW(IIG,JJG-1,KKG-1))
S12 = 0.5_EB*(DUDY+DVDX)
S13 = 0.5_EB*(DUDZ+DWDX)
S23 = 0.5_EB*(DVDZ+DWDY)
STRAIN_RATE(IIG,JJG,KKG) = SQRT(2._EB*(S11**2 + S22**2 + S33**2 + 2._EB*(S12**2 + S13**2 + S23**2)))
ENDDO WALL_LOOP
END SUBROUTINE COMPUTE_STRAIN_RATE
SUBROUTINE VISCOSITY_BC(NM)
! Specify ghost cell values of the viscosity array MU
INTEGER, INTENT(IN) :: NM
REAL(EB) :: MU_OTHER,DP_OTHER,KRES_OTHER
INTEGER :: II,JJ,KK,IW,IIO,JJO,KKO,NOM,N_INT_CELLS
TYPE(WALL_TYPE),POINTER :: WC=>NULL()
TYPE(EXTERNAL_WALL_TYPE),POINTER :: EWC=>NULL()
CALL POINT_TO_MESH(NM)
! Mirror viscosity into solids and exterior boundary cells
WALL_LOOP: DO IW=1,N_EXTERNAL_WALL_CELLS
WC =>WALL(IW)
EWC=>EXTERNAL_WALL(IW)
IF (EWC%NOM==0) CYCLE WALL_LOOP
II = WC%ONE_D%II
JJ = WC%ONE_D%JJ
KK = WC%ONE_D%KK
NOM = EWC%NOM
MU_OTHER = 0._EB
DP_OTHER = 0._EB
KRES_OTHER = 0._EB
DO KKO=EWC%KKO_MIN,EWC%KKO_MAX
DO JJO=EWC%JJO_MIN,EWC%JJO_MAX
DO IIO=EWC%IIO_MIN,EWC%IIO_MAX
MU_OTHER = MU_OTHER + OMESH(NOM)%MU(IIO,JJO,KKO)
KRES_OTHER = KRES_OTHER + OMESH(NOM)%KRES(IIO,JJO,KKO)
IF (PREDICTOR) THEN
DP_OTHER = DP_OTHER + OMESH(NOM)%D(IIO,JJO,KKO)
ELSE
DP_OTHER = DP_OTHER + OMESH(NOM)%DS(IIO,JJO,KKO)
ENDIF
ENDDO
ENDDO
ENDDO
N_INT_CELLS = (EWC%IIO_MAX-EWC%IIO_MIN+1) * (EWC%JJO_MAX-EWC%JJO_MIN+1) * (EWC%KKO_MAX-EWC%KKO_MIN+1)
MU_OTHER = MU_OTHER/REAL(N_INT_CELLS,EB)
KRES_OTHER = KRES_OTHER/REAL(N_INT_CELLS,EB)
DP_OTHER = DP_OTHER/REAL(N_INT_CELLS,EB)
MU(II,JJ,KK) = MU_OTHER
KRES(II,JJ,KK) = KRES_OTHER
IF (PREDICTOR) THEN
D(II,JJ,KK) = DP_OTHER
ELSE
DS(II,JJ,KK) = DP_OTHER
ENDIF
ENDDO WALL_LOOP
END SUBROUTINE VISCOSITY_BC
SUBROUTINE VELOCITY_FLUX(T,DT,NM)
! Compute convective and diffusive terms of the momentum equations
USE MATH_FUNCTIONS, ONLY: EVALUATE_RAMP
USE COMPLEX_GEOMETRY, ONLY: CCIBM_VELOCITY_FLUX
INTEGER, INTENT(IN) :: NM
REAL(EB), INTENT(IN) :: T,DT
REAL(EB) :: MUX,MUY,MUZ,UP,UM,VP,VM,WP,WM,VTRM,OMXP,OMXM,OMYP,OMYM,OMZP,OMZM,TXYP,TXYM,TXZP,TXZM,TYZP,TYZM, &
DTXYDY,DTXZDZ,DTYZDZ,DTXYDX,DTXZDX,DTYZDY, &
DUDX,DVDY,DWDZ,DUDY,DUDZ,DVDX,DVDZ,DWDX,DWDY, &
VOMZ,WOMY,UOMY,VOMX,UOMZ,WOMX, &
RRHO,GX(0:IBAR_MAX),GY(0:IBAR_MAX),GZ(0:IBAR_MAX),TXXP,TXXM,TYYP,TYYM,TZZP,TZZM,DTXXDX,DTYYDY,DTZZDZ, &
DUMMY=0._EB
REAL(EB) :: VEG_UMAG
INTEGER :: I,J,K,IEXP,IEXM,IEYP,IEYM,IEZP,IEZM,IC,IC1,IC2
REAL(EB), POINTER, DIMENSION(:,:,:) :: TXY=>NULL(),TXZ=>NULL(),TYZ=>NULL(),OMX=>NULL(),OMY=>NULL(),OMZ=>NULL(), &
UU=>NULL(),VV=>NULL(),WW=>NULL(),RHOP=>NULL(),DP=>NULL()
CALL POINT_TO_MESH(NM)
IF (PREDICTOR) THEN
UU => U
VV => V
WW => W
DP => D
RHOP => RHO
ELSE
UU => US
VV => VS
WW => WS
DP => DS
RHOP => RHOS
ENDIF
TXY => WORK1
TXZ => WORK2
TYZ => WORK3
OMX => WORK4
OMY => WORK5
OMZ => WORK6
! Compute vorticity and stress tensor components
!$OMP PARALLEL DO PRIVATE(DUDY, DVDX, DUDZ, DWDX, DVDZ, DWDY, &
!$OMP& MUX, MUY, MUZ) SCHEDULE(STATIC)
DO K=0,KBAR
DO J=0,JBAR
DO I=0,IBAR
DUDY = RDYN(J)*(UU(I,J+1,K)-UU(I,J,K))
DVDX = RDXN(I)*(VV(I+1,J,K)-VV(I,J,K))
DUDZ = RDZN(K)*(UU(I,J,K+1)-UU(I,J,K))
DWDX = RDXN(I)*(WW(I+1,J,K)-WW(I,J,K))
DVDZ = RDZN(K)*(VV(I,J,K+1)-VV(I,J,K))
DWDY = RDYN(J)*(WW(I,J+1,K)-WW(I,J,K))
OMX(I,J,K) = DWDY - DVDZ
OMY(I,J,K) = DUDZ - DWDX
OMZ(I,J,K) = DVDX - DUDY
MUX = 0.25_EB*(MU(I,J+1,K)+MU(I,J,K)+MU(I,J,K+1)+MU(I,J+1,K+1))
MUY = 0.25_EB*(MU(I+1,J,K)+MU(I,J,K)+MU(I,J,K+1)+MU(I+1,J,K+1))
MUZ = 0.25_EB*(MU(I+1,J,K)+MU(I,J,K)+MU(I,J+1,K)+MU(I+1,J+1,K))
TXY(I,J,K) = MUZ*(DVDX + DUDY)
TXZ(I,J,K) = MUY*(DUDZ + DWDX)
TYZ(I,J,K) = MUX*(DVDZ + DWDY)
ENDDO
ENDDO
ENDDO
!$OMP END PARALLEL DO
! Wannier Flow (Stokes flow) test case
IF (PERIODIC_TEST==5) THEN
OMX=0._EB
OMY=0._EB
OMZ=0._EB
OME_E = -1.E6_EB
ENDIF
! Compute gravity components
IF (.NOT.SPATIAL_GRAVITY_VARIATION) THEN
GX(0:IBAR) = EVALUATE_RAMP(T,DUMMY,I_RAMP_GX)*GVEC(1)
GY(0:IBAR) = EVALUATE_RAMP(T,DUMMY,I_RAMP_GY)*GVEC(2)
GZ(0:IBAR) = EVALUATE_RAMP(T,DUMMY,I_RAMP_GZ)*GVEC(3)
ELSE
DO I=0,IBAR
GX(I) = EVALUATE_RAMP(X(I),DUMMY,I_RAMP_GX)*GVEC(1)
GY(I) = EVALUATE_RAMP(X(I),DUMMY,I_RAMP_GY)*GVEC(2)
GZ(I) = EVALUATE_RAMP(X(I),DUMMY,I_RAMP_GZ)*GVEC(3)
ENDDO
ENDIF
! Compute x-direction flux term FVX
!$OMP PARALLEL PRIVATE(WP, WM, VP, VM, UP, UM, &
!$OMP& OMXP, OMXM, OMYP, OMYM, OMZP, OMZM, &
!$OMP& TXZP, TXZM, TXYP, TXYM, TYZP, TYZM, &
!$OMP& IC, IEXP, IEXM, IEYP, IEYM, IEZP, IEZM, &
!$OMP& RRHO, DUDX, DVDY, DWDZ, VTRM)
!$OMP DO SCHEDULE(static) &
!$OMP& PRIVATE(WOMY, VOMZ, TXXP, TXXM, DTXXDX, DTXYDY, DTXZDZ)
DO K=1,KBAR
DO J=1,JBAR
DO I=0,IBAR
WP = WW(I,J,K) + WW(I+1,J,K)
WM = WW(I,J,K-1) + WW(I+1,J,K-1)
VP = VV(I,J,K) + VV(I+1,J,K)
VM = VV(I,J-1,K) + VV(I+1,J-1,K)
OMYP = OMY(I,J,K)
OMYM = OMY(I,J,K-1)
OMZP = OMZ(I,J,K)
OMZM = OMZ(I,J-1,K)
TXZP = TXZ(I,J,K)
TXZM = TXZ(I,J,K-1)
TXYP = TXY(I,J,K)
TXYM = TXY(I,J-1,K)
IC = CELL_INDEX(I,J,K)
IEYP = EDGE_INDEX(8,IC)
IEYM = EDGE_INDEX(6,IC)
IEZP = EDGE_INDEX(12,IC)
IEZM = EDGE_INDEX(10,IC)
IF (OME_E(-1,IEYP)>-1.E5_EB) THEN
OMYP = OME_E(-1,IEYP)
TXZP = TAU_E(-1,IEYP)
ENDIF
IF (OME_E( 1,IEYM)>-1.E5_EB) THEN
OMYM = OME_E( 1,IEYM)
TXZM = TAU_E( 1,IEYM)
ENDIF
IF (OME_E(-2,IEZP)>-1.E5_EB) THEN
OMZP = OME_E(-2,IEZP)
TXYP = TAU_E(-2,IEZP)
ENDIF
IF (OME_E( 2,IEZM)>-1.E5_EB) THEN
OMZM = OME_E( 2,IEZM)
TXYM = TAU_E( 2,IEZM)
ENDIF
WOMY = WP*OMYP + WM*OMYM
VOMZ = VP*OMZP + VM*OMZM
RRHO = 2._EB/(RHOP(I,J,K)+RHOP(I+1,J,K))
DVDY = (VV(I+1,J,K)-VV(I+1,J-1,K))*RDY(J)
DWDZ = (WW(I+1,J,K)-WW(I+1,J,K-1))*RDZ(K)
TXXP = MU(I+1,J,K)*( FOTH*DP(I+1,J,K) - 2._EB*(DVDY+DWDZ) )
DVDY = (VV(I,J,K)-VV(I,J-1,K))*RDY(J)
DWDZ = (WW(I,J,K)-WW(I,J,K-1))*RDZ(K)
TXXM = MU(I,J,K) *( FOTH*DP(I,J,K) - 2._EB*(DVDY+DWDZ) )
DTXXDX= RDXN(I)*(TXXP-TXXM)
DTXYDY= RDY(J) *(TXYP-TXYM)
DTXZDZ= RDZ(K) *(TXZP-TXZM)
VTRM = DTXXDX + DTXYDY + DTXZDZ
FVX(I,J,K) = 0.25_EB*(WOMY - VOMZ) - GX(I) + RRHO*(GX(I)*RHO_0(K) - VTRM)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
! Compute y-direction flux term FVY
!$OMP DO SCHEDULE(static) &
!$OMP& PRIVATE(WOMX, UOMZ, TYYP, TYYM, DTXYDX, DTYYDY, DTYZDZ)
DO K=1,KBAR
DO J=0,JBAR
DO I=1,IBAR
UP = UU(I,J,K) + UU(I,J+1,K)
UM = UU(I-1,J,K) + UU(I-1,J+1,K)
WP = WW(I,J,K) + WW(I,J+1,K)
WM = WW(I,J,K-1) + WW(I,J+1,K-1)
OMXP = OMX(I,J,K)
OMXM = OMX(I,J,K-1)
OMZP = OMZ(I,J,K)
OMZM = OMZ(I-1,J,K)
TYZP = TYZ(I,J,K)
TYZM = TYZ(I,J,K-1)
TXYP = TXY(I,J,K)
TXYM = TXY(I-1,J,K)
IC = CELL_INDEX(I,J,K)
IEXP = EDGE_INDEX(4,IC)
IEXM = EDGE_INDEX(2,IC)
IEZP = EDGE_INDEX(12,IC)
IEZM = EDGE_INDEX(11,IC)
IF (OME_E(-2,IEXP)>-1.E5_EB) THEN
OMXP = OME_E(-2,IEXP)
TYZP = TAU_E(-2,IEXP)
ENDIF
IF (OME_E( 2,IEXM)>-1.E5_EB) THEN
OMXM = OME_E( 2,IEXM)
TYZM = TAU_E( 2,IEXM)
ENDIF
IF (OME_E(-1,IEZP)>-1.E5_EB) THEN
OMZP = OME_E(-1,IEZP)
TXYP = TAU_E(-1,IEZP)
ENDIF
IF (OME_E( 1,IEZM)>-1.E5_EB) THEN
OMZM = OME_E( 1,IEZM)
TXYM = TAU_E( 1,IEZM)
ENDIF
WOMX = WP*OMXP + WM*OMXM
UOMZ = UP*OMZP + UM*OMZM
RRHO = 2._EB/(RHOP(I,J,K)+RHOP(I,J+1,K))
DUDX = (UU(I,J+1,K)-UU(I-1,J+1,K))*RDX(I)
DWDZ = (WW(I,J+1,K)-WW(I,J+1,K-1))*RDZ(K)
TYYP = MU(I,J+1,K)*( FOTH*DP(I,J+1,K) - 2._EB*(DUDX+DWDZ) )
DUDX = (UU(I,J,K)-UU(I-1,J,K))*RDX(I)
DWDZ = (WW(I,J,K)-WW(I,J,K-1))*RDZ(K)
TYYM = MU(I,J,K) *( FOTH*DP(I,J,K) - 2._EB*(DUDX+DWDZ) )
DTXYDX= RDX(I) *(TXYP-TXYM)
DTYYDY= RDYN(J)*(TYYP-TYYM)
DTYZDZ= RDZ(K) *(TYZP-TYZM)
VTRM = DTXYDX + DTYYDY + DTYZDZ
FVY(I,J,K) = 0.25_EB*(UOMZ - WOMX) - GY(I) + RRHO*(GY(I)*RHO_0(K) - VTRM)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
! Compute z-direction flux term FVZ
!$OMP DO SCHEDULE(static) &
!$OMP& PRIVATE(UOMY, VOMX, TZZP, TZZM, DTXZDX, DTYZDY, DTZZDZ)
DO K=0,KBAR
DO J=1,JBAR
DO I=1,IBAR
UP = UU(I,J,K) + UU(I,J,K+1)
UM = UU(I-1,J,K) + UU(I-1,J,K+1)
VP = VV(I,J,K) + VV(I,J,K+1)
VM = VV(I,J-1,K) + VV(I,J-1,K+1)
OMYP = OMY(I,J,K)
OMYM = OMY(I-1,J,K)
OMXP = OMX(I,J,K)
OMXM = OMX(I,J-1,K)
TXZP = TXZ(I,J,K)
TXZM = TXZ(I-1,J,K)
TYZP = TYZ(I,J,K)
TYZM = TYZ(I,J-1,K)
IC = CELL_INDEX(I,J,K)
IEXP = EDGE_INDEX(4,IC)
IEXM = EDGE_INDEX(3,IC)
IEYP = EDGE_INDEX(8,IC)
IEYM = EDGE_INDEX(7,IC)
IF (OME_E(-1,IEXP)>-1.E5_EB) THEN
OMXP = OME_E(-1,IEXP)
TYZP = TAU_E(-1,IEXP)
ENDIF
IF (OME_E( 1,IEXM)>-1.E5_EB) THEN
OMXM = OME_E( 1,IEXM)
TYZM = TAU_E( 1,IEXM)
ENDIF
IF (OME_E(-2,IEYP)>-1.E5_EB) THEN
OMYP = OME_E(-2,IEYP)
TXZP = TAU_E(-2,IEYP)
ENDIF
IF (OME_E( 2,IEYM)>-1.E5_EB) THEN
OMYM = OME_E( 2,IEYM)
TXZM = TAU_E( 2,IEYM)
ENDIF
UOMY = UP*OMYP + UM*OMYM
VOMX = VP*OMXP + VM*OMXM
RRHO = 2._EB/(RHOP(I,J,K)+RHOP(I,J,K+1))
DUDX = (UU(I,J,K+1)-UU(I-1,J,K+1))*RDX(I)
DVDY = (VV(I,J,K+1)-VV(I,J-1,K+1))*RDY(J)
TZZP = MU(I,J,K+1)*( FOTH*DP(I,J,K+1) - 2._EB*(DUDX+DVDY) )
DUDX = (UU(I,J,K)-UU(I-1,J,K))*RDX(I)
DVDY = (VV(I,J,K)-VV(I,J-1,K))*RDY(J)
TZZM = MU(I,J,K) *( FOTH*DP(I,J,K) - 2._EB*(DUDX+DVDY) )
DTXZDX= RDX(I) *(TXZP-TXZM)
DTYZDY= RDY(J) *(TYZP-TYZM)
DTZZDZ= RDZN(K)*(TZZP-TZZM)
VTRM = DTXZDX + DTYZDY + DTZZDZ
FVZ(I,J,K) = 0.25_EB*(VOMX - UOMY) - GZ(I) + RRHO*(GZ(I)*0.5_EB*(RHO_0(K)+RHO_0(K+1)) - VTRM)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
!$OMP END PARALLEL
IF (EVACUATION_ONLY(NM)) THEN
FVZ = 0._EB
RETURN
END IF
! Additional force terms
IF (ANY(MEAN_FORCING) .AND. .NOT.NEW_MOMENTUM_NUDGING) CALL MOMENTUM_NUDGING ! Mean forcing
IF (ANY(MEAN_FORCING) .AND. NEW_MOMENTUM_NUDGING) CALL MOMENTUM_NUDGING_2 ! Test new mean forcing
IF (ANY(ABS(FVEC)>TWO_EPSILON_EB)) CALL DIRECT_FORCE ! Direct force
IF (ANY(ABS(OVEC)>TWO_EPSILON_EB)) CALL CORIOLIS_FORCE ! Coriolis force
IF (WFDS_BNDRYFUEL) CALL VEGETATION_DRAG ! Surface vegetation drag
IF (PATCH_VELOCITY) CALL PATCH_VELOCITY_FLUX(DT,NM) ! Specified patch velocity
IF (CC_IBM) CALL CCIBM_VELOCITY_FLUX(DT,NM) ! Direct-forcing Immersed Boundary Method
IF (PERIODIC_TEST==7) CALL MMS_VELOCITY_FLUX(NM,T) ! Source term in manufactured solution
CONTAINS
SUBROUTINE MOMENTUM_NUDGING
! Add a force vector to the momentum equation that moves the flow field towards the direction of the mean flow.
REAL(EB) :: UBAR,VBAR,WBAR,INTEGRAL,SUM_VOLUME,VC,UMEAN,VMEAN,WMEAN,DU_FORCING,DV_FORCING,DW_FORCING,DT_LOC
DT_LOC = MAX(DT,DT_MEAN_FORCING)
MEAN_FORCING_X: IF (MEAN_FORCING(1)) THEN
SELECT_RAMP_U: SELECT CASE(I_RAMP_U0_Z)
CASE(0) SELECT_RAMP_U
INTEGRAL = 0._EB
SUM_VOLUME = 0._EB
DO K=1,KBAR
DO J=1,JBAR
DO I=0,IBAR
IC1 = CELL_INDEX(I,J,K)
IC2 = CELL_INDEX(I+1,J,K)
IF (SOLID(IC1)) CYCLE
IF (SOLID(IC2)) CYCLE
IF (.NOT.MEAN_FORCING_CELL(I,J,K) ) CYCLE
IF (.NOT.MEAN_FORCING_CELL(I+1,J,K)) CYCLE
VC = DXN(I)*DY(J)*DZ(K)
INTEGRAL = INTEGRAL + UU(I,J,K)*VC
SUM_VOLUME = SUM_VOLUME + VC
ENDDO
ENDDO
ENDDO
IF (SUM_VOLUME>TWO_EPSILON_EB) THEN
UMEAN = INTEGRAL/SUM_VOLUME
ELSE
UMEAN = 0._EB
ENDIF
UBAR = U0*EVALUATE_RAMP(T,DUMMY,I_RAMP_U0_T)
DU_FORCING = (UBAR-UMEAN)/DT_LOC
DO K=1,KBAR
DO J=1,JBAR
DO I=0,IBAR
IF (.NOT.MEAN_FORCING_CELL(I,J,K) ) CYCLE
IF (.NOT.MEAN_FORCING_CELL(I+1,J,K)) CYCLE
FVX(I,J,K) = FVX(I,J,K) - DU_FORCING
ENDDO
ENDDO
ENDDO
CASE(1:) SELECT_RAMP_U
K_LOOP_U: DO K=1,KBAR
INTEGRAL = 0._EB
SUM_VOLUME = 0._EB
DO J=1,JBAR
DO I=0,IBAR
IC1 = CELL_INDEX(I,J,K)
IC2 = CELL_INDEX(I+1,J,K)
IF (SOLID(IC1)) CYCLE
IF (SOLID(IC2)) CYCLE
VC = DXN(I)*DY(J)*DZ(K)
INTEGRAL = INTEGRAL + UU(I,J,K)*VC
SUM_VOLUME = SUM_VOLUME + VC
ENDDO
ENDDO
IF (SUM_VOLUME>TWO_EPSILON_EB) THEN
UMEAN = INTEGRAL/SUM_VOLUME
ELSE
! this can happen if all cells in a given row, k, are solid
UMEAN = 0._EB
ENDIF
UBAR = U0*EVALUATE_RAMP(T,DUMMY,I_RAMP_U0_T)*EVALUATE_RAMP(ZC(K),DUMMY,I_RAMP_U0_Z)
DU_FORCING = (UBAR-UMEAN)/DT_LOC
FVX(:,:,K) = FVX(:,:,K) - DU_FORCING
ENDDO K_LOOP_U
END SELECT SELECT_RAMP_U
ENDIF MEAN_FORCING_X
MEAN_FORCING_Y: IF (MEAN_FORCING(2)) THEN
SELECT_RAMP_V: SELECT CASE(I_RAMP_V0_Z)
CASE(0) SELECT_RAMP_V
INTEGRAL = 0._EB
SUM_VOLUME = 0._EB
DO K=1,KBAR
DO J=0,JBAR
DO I=1,IBAR
IC1 = CELL_INDEX(I,J,K)
IC2 = CELL_INDEX(I,J+1,K)
IF (SOLID(IC1)) CYCLE
IF (SOLID(IC2)) CYCLE
IF (.NOT.MEAN_FORCING_CELL(I,J,K) ) CYCLE
IF (.NOT.MEAN_FORCING_CELL(I,J+1,K)) CYCLE
VC = DX(I)*DYN(J)*DZ(K)
INTEGRAL = INTEGRAL + VV(I,J,K)*VC
SUM_VOLUME = SUM_VOLUME + VC
ENDDO
ENDDO
ENDDO
IF (SUM_VOLUME>TWO_EPSILON_EB) THEN
VMEAN = INTEGRAL/SUM_VOLUME
ELSE
VMEAN = 0._EB
ENDIF
VBAR = V0*EVALUATE_RAMP(T,DUMMY,I_RAMP_V0_T)
DV_FORCING = (VBAR-VMEAN)/DT_LOC
DO K=1,KBAR
DO J=0,JBAR
DO I=1,IBAR
IF (.NOT.MEAN_FORCING_CELL(I,J,K) ) CYCLE
IF (.NOT.MEAN_FORCING_CELL(I,J+1,K)) CYCLE
FVY(I,J,K) = FVY(I,J,K) - DV_FORCING
ENDDO
ENDDO
ENDDO
CASE(1:) SELECT_RAMP_V
K_LOOP_V: DO K=1,KBAR
INTEGRAL = 0._EB
SUM_VOLUME = 0._EB
DO J=0,JBAR
DO I=1,IBAR
IC1 = CELL_INDEX(I,J,K)