-
Notifications
You must be signed in to change notification settings - Fork 0
/
mass.f90
994 lines (825 loc) · 33.1 KB
/
mass.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
MODULE MASS
! Compute the mass equation differences
USE PRECISION_PARAMETERS
USE MESH_POINTERS
IMPLICIT NONE
PRIVATE
CHARACTER(255), PARAMETER :: massid='$Id: mass.f90 9967 2012-02-02 22:20:19Z randy.mcdermott $'
CHARACTER(255), PARAMETER :: massrev='$Revision: 9967 $'
CHARACTER(255), PARAMETER :: massdate='$Date: 2012-02-02 14:20:19 -0800 (Thu, 02 Feb 2012) $'
REAL(EB), POINTER, DIMENSION(:,:,:,:) :: ZZP
REAL(EB), POINTER, DIMENSION(:,:,:) :: UU,VV,WW,RHOP,DP
PUBLIC MASS_FINITE_DIFFERENCES,DENSITY,GET_REV_mass
CONTAINS
SUBROUTINE MASS_FINITE_DIFFERENCES(NM)
! Compute spatial differences for density equation
USE COMP_FUNCTIONS, ONLY: SECOND
USE GLOBAL_CONSTANTS, ONLY: N_TRACKED_SPECIES,NULL_BOUNDARY,OPEN_BOUNDARY,INTERPOLATED_BOUNDARY, &
PREDICTOR,CORRECTOR,EVACUATION_ONLY,SOLID_PHASE_ONLY,TUSED,DEBUG_OPENMP,SOLID_BOUNDARY, &
NO_MASS_FLUX,SPECIFIED_MASS_FLUX,HVAC_BOUNDARY,FLUX_LIMITER,ENTHALPY_TRANSPORT
INTEGER, INTENT(IN) :: NM
REAL(EB) :: TNOW,ZZZ(4),UN,RHO_D_DZDN
INTEGER :: I,J,K,N,II,JJ,KK,IIG,JJG,KKG,IW,IOR,SURF_INDEX
REAL(EB), POINTER, DIMENSION(:,:,:) :: FX=>NULL(),FY=>NULL(),FZ=>NULL()
REAL(EB), POINTER, DIMENSION(:,:,:,:) :: ZZP=>NULL()
TYPE(WALL_TYPE), POINTER :: WC=>NULL()
IF (EVACUATION_ONLY(NM) .OR. SOLID_PHASE_ONLY) RETURN
TNOW=SECOND()
CALL POINT_TO_MESH(NM)
IF (PREDICTOR) THEN
UU => U
VV => V
WW => W
DP => D
RHOP => RHO
IF (N_TRACKED_SPECIES > 0) ZZP => ZZ
ELSE
UU => US
VV => VS
WW => WS
DP => DS
RHOP => RHOS
IF (N_TRACKED_SPECIES > 0) ZZP => ZZS
ENDIF
FX=>WORK4
FY=>WORK5
FZ=>WORK6
FX=0._EB
FY=0._EB
FZ=0._EB
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(KBAR,JBAR,IBAR,KBM1,JBM1,IBM1,RHOP,FX,FY,FZ,UU,VV,WW,FLUX_LIMITER,R, &
!$OMP N_EXTERNAL_WALL_CELLS,N_INTERNAL_WALL_CELLS,WALL_INDEX,CELL_INDEX, &
!$OMP RHO_F,UVW_SAVE, &
!$OMP FRHO,SOLID,RDX,RDY,RDZ,RRN, &
!$OMP N_TRACKED_SPECIES,ZZP,ZZ_F)
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZZ)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBM1
ZZZ(1:4) = RHOP(I-1:I+2,J,K)
FX(I,J,K) = UU(I,J,K)*SCALAR_FACE_VALUE(UU(I,J,K),ZZZ,FLUX_LIMITER)*R(I)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZZ)
DO K=1,KBAR
DO J=1,JBM1
DO I=1,IBAR
ZZZ(1:4) = RHOP(I,J-1:J+2,K)
FY(I,J,K) = VV(I,J,K)*SCALAR_FACE_VALUE(VV(I,J,K),ZZZ,FLUX_LIMITER)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZZ)
DO K=1,KBM1
DO J=1,JBAR
DO I=1,IBAR
ZZZ(1:4) = RHOP(I,J,K-1:K+2)
FZ(I,J,K) = WW(I,J,K)*SCALAR_FACE_VALUE(WW(I,J,K),ZZZ,FLUX_LIMITER)
ENDDO
ENDDO
ENDDO
!$OMP END DO
!$OMP DO PRIVATE(IW,II,JJ,KK,IOR,SURF_INDEX,IIG,JJG,KKG,ZZZ,UN)
WLOOP_FL: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS
WC=>WALL(IW)
IF (WC%BOUNDARY_TYPE==NULL_BOUNDARY) CYCLE WLOOP_FL
II = WC%II
JJ = WC%JJ
KK = WC%KK
IOR = WC%IOR
SURF_INDEX = WC%SURF_INDEX
IIG = WC%IIG
JJG = WC%JJG
KKG = WC%KKG
! overwrite first off-wall advective flux if flow is away from the wall and if the face is not also a wall cell
IF (WC%BOUNDARY_TYPE/=INTERPOLATED_BOUNDARY .AND. WC%BOUNDARY_TYPE/=OPEN_BOUNDARY) THEN
OFF_WALL_SELECT_1: SELECT CASE(IOR)
CASE( 1) OFF_WALL_SELECT_1
! ghost FX/UU(II+1)
! /// II /// II+1 | II+2 | ...
! ^ WALL_INDEX(II+1,+1)
IF ((UU(II+1,JJ,KK)>0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II+1,JJ,KK),+1)>0)) THEN
ZZZ(1:3) = (/WC%RHO_F,RHOP(II+1:II+2,JJ,KK)/)
FX(II+1,JJ,KK) = UU(II+1,JJ,KK)*SCALAR_FACE_VALUE(UU(II+1,JJ,KK),ZZZ,FLUX_LIMITER)*R(II+1)
ENDIF
CASE(-1) OFF_WALL_SELECT_1
! FX/UU(II-2) ghost
! ... | II-2 | II-1 /// II ///
! ^ WALL_INDEX(II-1,-1)
IF ((UU(II-2,JJ,KK)<0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II-1,JJ,KK),-1)>0)) THEN
ZZZ(2:4) = (/RHOP(II-2:II-1,JJ,KK),WC%RHO_F/)
FX(II-2,JJ,KK) = UU(II-2,JJ,KK)*SCALAR_FACE_VALUE(UU(II-2,JJ,KK),ZZZ,FLUX_LIMITER)*R(II-2)
ENDIF
CASE( 2) OFF_WALL_SELECT_1
IF ((VV(II,JJ+1,KK)>0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ+1,KK),+2)>0)) THEN
ZZZ(1:3) = (/WC%RHO_F,RHOP(II,JJ+1:JJ+2,KK)/)
FY(II,JJ+1,KK) = VV(II,JJ+1,KK)*SCALAR_FACE_VALUE(VV(II,JJ+1,KK),ZZZ,FLUX_LIMITER)
ENDIF
CASE(-2) OFF_WALL_SELECT_1
IF ((VV(II,JJ-2,KK)<0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ-1,KK),-2)>0)) THEN
ZZZ(2:4) = (/RHOP(II,JJ-2:JJ-1,KK),WC%RHO_F/)
FY(II,JJ-2,KK) = VV(II,JJ-2,KK)*SCALAR_FACE_VALUE(VV(II,JJ-2,KK),ZZZ,FLUX_LIMITER)
ENDIF
CASE( 3) OFF_WALL_SELECT_1
IF ((WW(II,JJ,KK+1)>0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ,KK+1),+3)>0)) THEN
ZZZ(1:3) = (/WC%RHO_F,RHOP(II,JJ,KK+1:KK+2)/)
FZ(II,JJ,KK+1) = WW(II,JJ,KK+1)*SCALAR_FACE_VALUE(WW(II,JJ,KK+1),ZZZ,FLUX_LIMITER)
ENDIF
CASE(-3) OFF_WALL_SELECT_1
IF ((WW(II,JJ,KK-2)<0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ,KK-1),-3)>0)) THEN
ZZZ(2:4) = (/RHOP(II,JJ,KK-2:KK-1),WC%RHO_F/)
FZ(II,JJ,KK-2) = WW(II,JJ,KK-2)*SCALAR_FACE_VALUE(WW(II,JJ,KK-2),ZZZ,FLUX_LIMITER)
ENDIF
END SELECT OFF_WALL_SELECT_1
ENDIF
SELECT CASE(IOR)
CASE( 1)
UN = UU(II,JJ,KK)
CASE(-1)
UN = UU(II-1,JJ,KK)
CASE( 2)
UN = VV(II,JJ,KK)
CASE(-2)
UN = VV(II,JJ-1,KK)
CASE( 3)
UN = WW(II,JJ,KK)
CASE(-3)
UN = WW(II,JJ,KK-1)
END SELECT
! In case of interpolated boundary, use the original velocity, not the averaged value
IF (WC%BOUNDARY_TYPE==INTERPOLATED_BOUNDARY) UN = UVW_SAVE(IW)
IF (ENTHALPY_TRANSPORT) THEN
IF (SURFACE(SURF_INDEX)%SPECIES_BC_INDEX==SPECIFIED_MASS_FLUX) THEN
IF (PREDICTOR) UN = -SIGN(1._EB,REAL(IOR,EB))*WC%UWS
IF (CORRECTOR) UN = -SIGN(1._EB,REAL(IOR,EB))*WC%UW
ENDIF
ENDIF
! Compute flux on the face of the wall cell
SELECT CASE(IOR)
CASE( 1)
FX(II,JJ,KK) = UN*WC%RHO_F*R(II)
CASE(-1)
FX(II-1,JJ,KK) = UN*WC%RHO_F*R(II-1)
CASE( 2)
FY(II,JJ,KK) = UN*WC%RHO_F
CASE(-2)
FY(II,JJ-1,KK) = UN*WC%RHO_F
CASE( 3)
FZ(II,JJ,KK) = UN*WC%RHO_F
CASE(-3)
FZ(II,JJ,KK-1) = UN*WC%RHO_F
END SELECT
ENDDO WLOOP_FL
!$OMP END DO NOWAIT
!$OMP WORKSHARE
FRHO = 0._EB
!$OMP END WORKSHARE
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
FRHO(I,J,K) = (FX(I,J,K)-FX(I-1,J,K))*RDX(I)*RRN(I) &
+ (FY(I,J,K)-FY(I,J-1,K))*RDY(J) &
+ (FZ(I,J,K)-FZ(I,J,K-1))*RDZ(K)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
!$OMP END PARALLEL
! Compute the species equation differences
SPECIES_LOOP: DO N=1,N_TRACKED_SPECIES
FX=0._EB
FY=0._EB
FZ=0._EB
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(N,KBAR,JBAR,IBAR,KBM1,JBM1,IBM1,RHOP,ZZP,FX,FY,FZ,UU,VV,WW,FLUX_LIMITER,R, &
!$OMP N_EXTERNAL_WALL_CELLS,N_INTERNAL_WALL_CELLS, &
!$OMP WALL_INDEX,CELL_INDEX,RHO_F,ZZ_F,UVW_SAVE, &
!$OMP SURFACE,UWS,RHODW,RDN,MASSFLUX, &
!$OMP SOLID,DEL_RHO_D_DEL_Z,RDX,RDY,RDZ,RRN)
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZZ)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBM1
ZZZ(1:4) = RHOP(I-1:I+2,J,K)*ZZP(I-1:I+2,J,K,N)
FX(I,J,K) = UU(I,J,K)*SCALAR_FACE_VALUE(UU(I,J,K),ZZZ,FLUX_LIMITER)*R(I)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZZ)
DO K=1,KBAR
DO J=1,JBM1
DO I=1,IBAR
ZZZ(1:4) = RHOP(I,J-1:J+2,K)*ZZP(I,J-1:J+2,K,N)
FY(I,J,K) = VV(I,J,K)*SCALAR_FACE_VALUE(VV(I,J,K),ZZZ,FLUX_LIMITER)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZZ)
DO K=1,KBM1
DO J=1,JBAR
DO I=1,IBAR
ZZZ(1:4) = RHOP(I,J,K-1:K+2)*ZZP(I,J,K-1:K+2,N)
FZ(I,J,K) = WW(I,J,K)*SCALAR_FACE_VALUE(WW(I,J,K),ZZZ,FLUX_LIMITER)
ENDDO
ENDDO
ENDDO
!$OMP END DO
!$OMP DO SCHEDULE(STATIC) &
!$OMP PRIVATE(IW,II,JJ,KK,IOR,SURF_INDEX,IIG,JJG,KKG,ZZZ,UN,RHO_D_DZDN)
WLOOP2_FL: DO IW=1,N_EXTERNAL_WALL_CELLS+N_INTERNAL_WALL_CELLS
WC=>WALL(IW)
IF (WC%BOUNDARY_TYPE==NULL_BOUNDARY) CYCLE WLOOP2_FL
II = WC%II
JJ = WC%JJ
KK = WC%KK
IOR = WC%IOR
SURF_INDEX = WC%SURF_INDEX
IIG = WC%IIG
JJG = WC%JJG
KKG = WC%KKG
! overwrite first off-wall advective flux if flow is away from the wall and if the face is not also a wall cell
IF (WC%BOUNDARY_TYPE/=INTERPOLATED_BOUNDARY .AND. WC%BOUNDARY_TYPE/=OPEN_BOUNDARY) THEN
OFF_WALL_SELECT_2: SELECT CASE(IOR)
CASE( 1) OFF_WALL_SELECT_2
! ghost FX/UU(II+1)
! /// II /// II+1 | II+2 | ...
! ^ WALL_INDEX(II+1,+1)
IF ((UU(II+1,JJ,KK)>0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II+1,JJ,KK),+1)>0)) THEN
ZZZ(1:3) = (/WC%RHO_F,RHOP(II+1:II+2,JJ,KK)/)*(/WC%ZZ_F(N),ZZP(II+1:II+2,JJ,KK,N)/)
FX(II+1,JJ,KK) = UU(II+1,JJ,KK)*SCALAR_FACE_VALUE(UU(II+1,JJ,KK),ZZZ,FLUX_LIMITER)*R(II+1)
ENDIF
CASE(-1) OFF_WALL_SELECT_2
! FX/UU(II-2) ghost
! ... | II-2 | II-1 /// II ///
! ^ WALL_INDEX(II-1,-1)
IF ((UU(II-2,JJ,KK)<0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II-1,JJ,KK),-1)>0)) THEN
ZZZ(2:4) = (/RHOP(II-2:II-1,JJ,KK),WC%RHO_F/)*(/ZZP(II-2:II-1,JJ,KK,N),WC%ZZ_F(N)/)
FX(II-2,JJ,KK) = UU(II-2,JJ,KK)*SCALAR_FACE_VALUE(UU(II-2,JJ,KK),ZZZ,FLUX_LIMITER)*R(II-2)
ENDIF
CASE( 2) OFF_WALL_SELECT_2
IF ((VV(II,JJ+1,KK)>0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ+1,KK),+2)>0)) THEN
ZZZ(1:3) = (/WC%RHO_F,RHOP(II,JJ+1:JJ+2,KK)/)*(/WC%ZZ_F(N),ZZP(II,JJ+1:JJ+2,KK,N)/)
FY(II,JJ+1,KK) = VV(II,JJ+1,KK)*SCALAR_FACE_VALUE(VV(II,JJ+1,KK),ZZZ,FLUX_LIMITER)
ENDIF
CASE(-2) OFF_WALL_SELECT_2
IF ((VV(II,JJ-2,KK)<0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ-1,KK),-2)>0)) THEN
ZZZ(2:4) = (/RHOP(II,JJ-2:JJ-1,KK),WC%RHO_F/)*(/ZZP(II,JJ-2:JJ-1,KK,N),WC%ZZ_F(N)/)
FY(II,JJ-2,KK) = VV(II,JJ-2,KK)*SCALAR_FACE_VALUE(VV(II,JJ-2,KK),ZZZ,FLUX_LIMITER)
ENDIF
CASE( 3) OFF_WALL_SELECT_2
IF ((WW(II,JJ,KK+1)>0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ,KK+1),+3)>0)) THEN
ZZZ(1:3) = (/WC%RHO_F,RHOP(II,JJ,KK+1:KK+2)/)*(/WC%ZZ_F(N),ZZP(II,JJ,KK+1:KK+2,N)/)
FZ(II,JJ,KK+1) = WW(II,JJ,KK+1)*SCALAR_FACE_VALUE(WW(II,JJ,KK+1),ZZZ,FLUX_LIMITER)
ENDIF
CASE(-3) OFF_WALL_SELECT_2
IF ((WW(II,JJ,KK-2)<0._EB) .AND. .NOT.(WALL_INDEX(CELL_INDEX(II,JJ,KK-1),-3)>0)) THEN
ZZZ(2:4) = (/RHOP(II,JJ,KK-2:KK-1),WC%RHO_F/)*(/ZZP(II,JJ,KK-2:KK-1,N),WC%ZZ_F(N)/)
FZ(II,JJ,KK-2) = WW(II,JJ,KK-2)*SCALAR_FACE_VALUE(WW(II,JJ,KK-2),ZZZ,FLUX_LIMITER)
ENDIF
END SELECT OFF_WALL_SELECT_2
ENDIF
! Get the normal components of velocity at the wall
SELECT CASE(IOR)
CASE( 1)
UN = UU(II,JJ,KK)
CASE(-1)
UN = UU(II-1,JJ,KK)
CASE( 2)
UN = VV(II,JJ,KK)
CASE(-2)
UN = VV(II,JJ-1,KK)
CASE( 3)
UN = WW(II,JJ,KK)
CASE(-3)
UN = WW(II,JJ,KK-1)
END SELECT
! At interpolated boundaries, use the actual normal components of velocity, not the ones that have been forced to match.
! This line is temporarily commented out because it may not be necessary.
!! IF (WC%BOUNDARY_TYPE==INTERPOLATED_BOUNDARY) UN = UVW_SAVE(IW)
! Special handling of normal component of velocity at surfaces with a specified mass flux.
IF ((SURFACE(SURF_INDEX)%SPECIES_BC_INDEX==SPECIFIED_MASS_FLUX .OR. &
(SURFACE(SURF_INDEX)%SPECIES_BC_INDEX==HVAC_BOUNDARY .OR. &
ANY(SURFACE(SURF_INDEX)%LEAK_PATH>0._EB)) .AND. WC%UWS<0._EB) .AND. WC%ZZ_F(N)>0._EB) THEN
! recreate diffusive flux from divg b/c UWP based on old RHODW
RHO_D_DZDN = 2._EB*WC%RHODW(N)*(ZZP(IIG,JJG,KKG,N)-WC%ZZ_F(N))*WC%RDN
UN = SIGN(1._EB,REAL(IOR,EB))*(WC%MASSFLUX(N) + RHO_D_DZDN)/(WC%RHO_F*WC%ZZ_F(N))
ENDIF
! Compute species mass flux on the face of the wall cell
SELECT CASE(IOR)
CASE( 1)
FX(II,JJ,KK) = UN*WC%RHO_F*WC%ZZ_F(N)*R(II)
CASE(-1)
FX(II-1,JJ,KK) = UN*WC%RHO_F*WC%ZZ_F(N)*R(II-1)
CASE( 2)
FY(II,JJ,KK) = UN*WC%RHO_F*WC%ZZ_F(N)
CASE(-2)
FY(II,JJ-1,KK) = UN*WC%RHO_F*WC%ZZ_F(N)
CASE( 3)
FZ(II,JJ,KK) = UN*WC%RHO_F*WC%ZZ_F(N)
CASE(-3)
FZ(II,JJ,KK-1) = UN*WC%RHO_F*WC%ZZ_F(N)
END SELECT
ENDDO WLOOP2_FL
!$OMP END DO
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
DEL_RHO_D_DEL_Z(I,J,K,N) = -DEL_RHO_D_DEL_Z(I,J,K,N) &
+ (FX(I,J,K)-FX(I-1,J,K))*RDX(I)*RRN(I) &
+ (FY(I,J,K)-FY(I,J-1,K))*RDY(J) &
+ (FZ(I,J,K)-FZ(I,J,K-1))*RDZ(K)
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
!$OMP END PARALLEL
ENDDO SPECIES_LOOP
TUSED(3,NM)=TUSED(3,NM)+SECOND()-TNOW
END SUBROUTINE MASS_FINITE_DIFFERENCES
SUBROUTINE DENSITY(NM)
! Update the density and species mass fractions
USE COMP_FUNCTIONS, ONLY: SECOND
USE PHYSICAL_FUNCTIONS, ONLY : GET_SPECIFIC_GAS_CONSTANT,GET_SENSIBLE_ENTHALPY
USE GLOBAL_CONSTANTS, ONLY: N_TRACKED_SPECIES,TMPMAX,TMPMIN,EVACUATION_ONLY, &
PREDICTOR,CORRECTOR,CHANGE_TIME_STEP,TMPA,N_ZONE, &
GAS_SPECIES, R0,SOLID_PHASE_ONLY,TUSED, &
DEBUG_OPENMP,CLIP_MASS_FRACTION,ENTHALPY_TRANSPORT
REAL(EB) :: DTRATIO,OMDTRATIO,TNOW,ZZ_GET(0:N_TRACKED_SPECIES),H_S
INTEGER :: I,J,K,N
INTEGER, INTENT(IN) :: NM
IF (EVACUATION_ONLY(NM)) RETURN
IF (SOLID_PHASE_ONLY) RETURN
TNOW=SECOND()
CALL POINT_TO_MESH(NM)
PREDICTOR_STEP: SELECT CASE (PREDICTOR)
CASE(.TRUE.) PREDICTOR_STEP
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(CHANGE_TIME_STEP,NM,N_TRACKED_SPECIES,KBAR,JBAR,IBAR,SOLID,CELL_INDEX,ZZS,DT,DEL_RHO_D_DEL_Z, &
!$OMP RHO,ZZ, &
!$OMP DTRATIO,DT_PREV,OMDTRATIO,RHOS, &
!$OMP FRHO,CLIP_MASS_FRACTION,N_ZONE,PBAR_S,PBAR,D_PBAR_DT,KBP1,JBP1,IBP1,RSUM,TMP,PRESSURE_ZONE, &
!$OMP TMPMIN,TMPMAX)
IF (.NOT.CHANGE_TIME_STEP(NM)) THEN
! NOTE: This IF statement is required because the source terms for species are zeroed out at
! the beginning of DIVERGENCE_PART_1, but the array also stores the divergence of the advective
! flux which is computed once in MASS_FINITE_DIFFERNENCES above, outside the CHANGE_TIME_STEP loop.
! DIVERGENCE_PART_1 is inside the loop. The source terms are then applied to the next substep in
! MASS_FINITE_DIFFERENCES.
! Store enthalpy for time derivative in divergence
IF (ENTHALPY_TRANSPORT) THEN
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZ_GET)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (N_TRACKED_SPECIES>0) ZZ_GET(1:N_TRACKED_SPECIES) = ZZ(I,J,K,1:N_TRACKED_SPECIES)
CALL GET_SENSIBLE_ENTHALPY(ZZ_GET,H_S,TMP(I,J,K))
RHO_H_S_OVER_PBAR(I,J,K) = RHO(I,J,K)*H_S/PBAR(K,PRESSURE_ZONE(I,J,K))
ENDDO
ENDDO
ENDDO
!$OMP END DO
ENDIF
!$OMP DO COLLAPSE(4) SCHEDULE(DYNAMIC) PRIVATE(N,K,J,I)
DO N=1,N_TRACKED_SPECIES
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
ZZS(I,J,K,N) = RHO(I,J,K)*ZZ(I,J,K,N) - DT*DEL_RHO_D_DEL_Z(I,J,K,N)
ENDDO
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
ELSE
!$OMP SINGLE
DTRATIO = DT/DT_PREV
OMDTRATIO = 1._EB - DTRATIO
!$OMP END SINGLE
!$OMP DO COLLAPSE(4) SCHEDULE(DYNAMIC) PRIVATE(N,K,J,I)
DO N=1,N_TRACKED_SPECIES
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
ZZS(I,J,K,N) = OMDTRATIO*RHO(I,J,K)*ZZ(I,J,K,N) + DTRATIO*RHOS(I,J,K)*ZZS(I,J,K,N)
ENDDO
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
ENDIF
! Predict the density at the next time step (RHOS or RHO^*)
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
RHOS(I,J,K) = RHO(I,J,K)-DT*FRHO(I,J,K)
ENDDO
ENDDO
ENDDO
!$OMP END DO
! Correct densities above or below clip limits
!$OMP SINGLE
CALL CHECK_DENSITY
!$OMP END SINGLE
! Extract mass fraction from RHO * ZZ
!$OMP DO COLLAPSE(4) SCHEDULE(DYNAMIC) PRIVATE(N,K,J,I)
DO N=1,N_TRACKED_SPECIES
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
ZZS(I,J,K,N) = ZZS(I,J,K,N)/RHOS(I,J,K)
ENDDO
ENDDO
ENDDO
ENDDO
!$OMP END DO
! Correct mass fractions above or below clip limits
IF (CLIP_MASS_FRACTION) THEN
!$OMP WORKSHARE
ZZS(1:IBAR,1:JBAR,1:KBAR,1:N_TRACKED_SPECIES) = MAX(0._EB,MIN(1._EB,ZZS(1:IBAR,1:JBAR,1:KBAR,1:N_TRACKED_SPECIES)))
!$OMP END WORKSHARE
ELSE
!$OMP SINGLE
CALL CHECK_MASS_FRACTION
!$OMP END SINGLE
ENDIF
! Predict background pressure at next time step
!$OMP DO PRIVATE(I)
DO I=1,N_ZONE
PBAR_S(:,I) = PBAR(:,I) + D_PBAR_DT(I)*DT
ENDDO
!$OMP END DO NOWAIT
! Compute molecular weight term RSUM=R0*SUM(Y_i/M_i)
IF (N_TRACKED_SPECIES>0) THEN
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZ_GET)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
ZZ_GET(1:N_TRACKED_SPECIES) = ZZS(I,J,K,1:N_TRACKED_SPECIES)
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM(I,J,K))
ENDDO
ENDDO
ENDDO
!$OMP END DO
ENDIF
! Extract predicted temperature at next time step from Equation of State
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
TMP(I,J,K) = PBAR_S(K,PRESSURE_ZONE(I,J,K))/(RSUM(I,J,K)*RHOS(I,J,K))
ENDDO
ENDDO
ENDDO
!$OMP END DO
!$OMP WORKSHARE
TMP = MAX(TMPMIN,MIN(TMPMAX,TMP))
!$OMP END WORKSHARE NOWAIT
!$OMP END PARALLEL
! The CORRECTOR step
CASE(.FALSE.) PREDICTOR_STEP
! Store enthalpy for time derivative in divergence
IF (ENTHALPY_TRANSPORT) THEN
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZ_GET)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (N_TRACKED_SPECIES>0) ZZ_GET(1:N_TRACKED_SPECIES) = ZZS(I,J,K,1:N_TRACKED_SPECIES)
CALL GET_SENSIBLE_ENTHALPY(ZZ_GET,H_S,TMP(I,J,K))
RHO_H_S_OVER_PBAR(I,J,K) = 0.5_EB*( RHO_H_S_OVER_PBAR(I,J,K) + RHOS(I,J,K)*H_S/PBAR_S(K,PRESSURE_ZONE(I,J,K)) )
ENDDO
ENDDO
ENDDO
!$OMP END DO
ENDIF
! Correct species mass fraction at next time step (ZZ here actually means ZZ*RHO)
!$OMP PARALLEL DEFAULT(NONE) &
!$OMP SHARED(N_TRACKED_SPECIES,KBAR,JBAR,IBAR,SOLID,CELL_INDEX,ZZ,RHO,RHOS,ZZS,DT,DEL_RHO_D_DEL_Z, &
!$OMP FRHO,CLIP_MASS_FRACTION,N_ZONE,PBAR,PBAR_S,D_PBAR_S_DT,KBP1,JBP1,IBP1,RSUM,TMP,PRESSURE_ZONE, &
!$OMP TMPMIN,TMPMAX)
!$OMP DO COLLAPSE(4) SCHEDULE(DYNAMIC) PRIVATE(N,K,J,I)
DO N=1,N_TRACKED_SPECIES
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
ZZ(I,J,K,N) = .5_EB*(RHO(I,J,K)*ZZ(I,J,K,N) + RHOS(I,J,K)*ZZS(I,J,K,N) - DT*DEL_RHO_D_DEL_Z(I,J,K,N))
ENDDO
ENDDO
ENDDO
ENDDO
!$OMP END DO NOWAIT
! Correct density at next time step
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
RHO(I,J,K) = .5_EB*(RHO(I,J,K)+RHOS(I,J,K)-DT*FRHO(I,J,K))
ENDDO
ENDDO
ENDDO
!$OMP END DO
! Correct densities above or below clip limits
!$OMP SINGLE
CALL CHECK_DENSITY
!$OMP END SINGLE
! Extract Y_n from rho*Y_n
!$OMP DO COLLAPSE(4) SCHEDULE(DYNAMIC) PRIVATE(N,K,J,I)
DO N=1,N_TRACKED_SPECIES
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (SOLID(CELL_INDEX(I,J,K))) CYCLE
ZZ(I,J,K,N) = ZZ(I,J,K,N)/RHO(I,J,K)
ENDDO
ENDDO
ENDDO
ENDDO
!$OMP END DO
! Correct mass fractions above or below clip limits
IF (CLIP_MASS_FRACTION) THEN
!$OMP WORKSHARE
ZZ(1:IBAR,1:JBAR,1:KBAR,1:N_TRACKED_SPECIES) = MAX(0._EB,MIN(1._EB,ZZ(1:IBAR,1:JBAR,1:KBAR,1:N_TRACKED_SPECIES)))
!$OMP END WORKSHARE NOWAIT
ELSE
!$OMP SINGLE
CALL CHECK_MASS_FRACTION
!$OMP END SINGLE
ENDIF
! Correct background pressure
!$OMP DO PRIVATE(I)
DO I=1,N_ZONE
PBAR(:,I) = .5_EB*(PBAR(:,I) + PBAR_S(:,I) + D_PBAR_S_DT(I)*DT)
ENDDO
!$OMP END DO NOWAIT
! Compute molecular weight term RSUM=R0*SUM(Y_i/M_i)
IF (N_TRACKED_SPECIES>0) THEN
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I,ZZ_GET)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
ZZ_GET(1:N_TRACKED_SPECIES) = ZZ(I,J,K,1:N_TRACKED_SPECIES)
CALL GET_SPECIFIC_GAS_CONSTANT(ZZ_GET,RSUM(I,J,K))
ENDDO
ENDDO
ENDDO
!$OMP END DO
ENDIF
! Extract predicted temperature at next time step from Equation of State
!$OMP DO COLLAPSE(3) SCHEDULE(STATIC) PRIVATE(K,J,I)
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
TMP(I,J,K) = PBAR(K,PRESSURE_ZONE(I,J,K))/(RSUM(I,J,K)*RHO(I,J,K))
ENDDO
ENDDO
ENDDO
!$OMP END DO
!$OMP WORKSHARE
TMP = MAX(TMPMIN,MIN(TMPMAX,TMP))
!$OMP END WORKSHARE NOWAIT
!$OMP END PARALLEL
END SELECT PREDICTOR_STEP
TUSED(3,NM)=TUSED(3,NM)+SECOND()-TNOW
END SUBROUTINE DENSITY
SUBROUTINE CHECK_DENSITY
! Redistribute mass from cells below or above the density cut-off limits
! Do not apply OpenMP to this routine
USE GLOBAL_CONSTANTS, ONLY : PREDICTOR,RHOMIN,RHOMAX
REAL(EB) :: MASS_N(-3:3),CONST,MASS_C,RHO_CUT,VC(-3:3),SIGN_FACTOR,SUM_MASS_N
INTEGER :: IC,I,J,K
REAL(EB), POINTER, DIMENSION(:,:,:) :: DELTA_RHO=>NULL()
DELTA_RHO => WORK2
DELTA_RHO = 0._EB
IF (PREDICTOR) THEN
RHOP=>RHOS
ELSE
RHOP=>RHO
ENDIF
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (RHOP(I,J,K)>=RHOMIN .AND. RHOP(I,J,K)<=RHOMAX) CYCLE
IC = CELL_INDEX(I,J,K)
IF (SOLID(IC)) CYCLE
IF (RHOP(I,J,K)<RHOMIN) THEN
RHO_CUT = RHOMIN
SIGN_FACTOR = 1._EB
ELSE
RHO_CUT = RHOMAX
SIGN_FACTOR = -1._EB
ENDIF
MASS_N = 0._EB
VC( 0) = DX(I) *DY(J) *DZ(K)
VC(-1) = DX(I-1)*DY(J) *DZ(K)
VC( 1) = DX(I+1)*DY(J) *DZ(K)
VC(-2) = DX(I) *DY(J-1)*DZ(K)
VC( 2) = DX(I) *DY(J+1)*DZ(K)
VC(-3) = DX(I) *DY(J) *DZ(K-1)
VC( 3) = DX(I) *DY(J) *DZ(K+1)
MASS_C = ABS(RHO_CUT-RHOP(I,J,K))*VC(0)
IF (WALL_INDEX(IC,-1)==0) MASS_N(-1) = ABS(MIN(RHOMAX,MAX(RHOMIN,RHOP(I-1,J,K)))-RHO_CUT)*VC(-1)
IF (WALL_INDEX(IC, 1)==0) MASS_N( 1) = ABS(MIN(RHOMAX,MAX(RHOMIN,RHOP(I+1,J,K)))-RHO_CUT)*VC( 1)
IF (WALL_INDEX(IC,-2)==0) MASS_N(-2) = ABS(MIN(RHOMAX,MAX(RHOMIN,RHOP(I,J-1,K)))-RHO_CUT)*VC(-2)
IF (WALL_INDEX(IC, 2)==0) MASS_N( 2) = ABS(MIN(RHOMAX,MAX(RHOMIN,RHOP(I,J+1,K)))-RHO_CUT)*VC( 2)
IF (WALL_INDEX(IC,-3)==0) MASS_N(-3) = ABS(MIN(RHOMAX,MAX(RHOMIN,RHOP(I,J,K-1)))-RHO_CUT)*VC(-3)
IF (WALL_INDEX(IC, 3)==0) MASS_N( 3) = ABS(MIN(RHOMAX,MAX(RHOMIN,RHOP(I,J,K+1)))-RHO_CUT)*VC( 3)
SUM_MASS_N = SUM(MASS_N)
IF (SUM_MASS_N<=ZERO_P) CYCLE
CONST = SIGN_FACTOR*MIN(1._EB,MASS_C/SUM_MASS_N)
DELTA_RHO(I,J,K) = DELTA_RHO(I,J,K) + CONST*SUM_MASS_N/VC( 0)
DELTA_RHO(I-1,J,K) = DELTA_RHO(I-1,J,K) - CONST*MASS_N(-1)/VC(-1)
DELTA_RHO(I+1,J,K) = DELTA_RHO(I+1,J,K) - CONST*MASS_N( 1)/VC( 1)
DELTA_RHO(I,J-1,K) = DELTA_RHO(I,J-1,K) - CONST*MASS_N(-2)/VC(-2)
DELTA_RHO(I,J+1,K) = DELTA_RHO(I,J+1,K) - CONST*MASS_N( 2)/VC( 2)
DELTA_RHO(I,J,K-1) = DELTA_RHO(I,J,K-1) - CONST*MASS_N(-3)/VC(-3)
DELTA_RHO(I,J,K+1) = DELTA_RHO(I,J,K+1) - CONST*MASS_N( 3)/VC( 3)
ENDDO
ENDDO
ENDDO
RHOP(1:IBAR,1:JBAR,1:KBAR) = MIN(RHOMAX,MAX(RHOMIN,RHOP(1:IBAR,1:JBAR,1:KBAR)+DELTA_RHO(1:IBAR,1:JBAR,1:KBAR)))
END SUBROUTINE CHECK_DENSITY
SUBROUTINE CHECK_MASS_FRACTION
! Redistribute species mass from cells below or above the cut-off limits
! Do not apply OpenMP to this routine
USE GLOBAL_CONSTANTS, ONLY : PREDICTOR,N_TRACKED_SPECIES
REAL(EB) :: SUM,CONST,MASS_C,MASS_N(-3:3),ZZ_CUT,SIGN_FACTOR,SUM_MASS_N,VC(-3:3)
INTEGER :: IC,N,I,J,K
REAL(EB), POINTER, DIMENSION(:,:,:) :: DELTA_ZZ=>NULL()
DELTA_ZZ => WORK1
IF (PREDICTOR) THEN
RHOP => RHOS
ZZP => ZZS
ELSE
RHOP => RHO
ZZP => ZZ
ENDIF
SPECIES_LOOP: DO N=1,N_TRACKED_SPECIES
DELTA_ZZ = 0._EB
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
IF (ZZP(I,J,K,N)>=0._EB .AND. ZZP(I,J,K,N)<=1._EB) CYCLE
IC = CELL_INDEX(I,J,K)
IF (SOLID(IC)) CYCLE
IF (ZZP(I,J,K,N)<0._EB) THEN
ZZ_CUT = 0._EB
SIGN_FACTOR = 1._EB
ELSE
ZZ_CUT = 1._EB
SIGN_FACTOR = -1._EB
ENDIF
MASS_N = 0._EB
VC( 0) = DX(I) *DY(J) *DZ(K)
VC(-1) = DX(I-1)*DY(J) *DZ(K)
VC( 1) = DX(I+1)*DY(J) *DZ(K)
VC(-2) = DX(I) *DY(J-1)*DZ(K)
VC( 2) = DX(I) *DY(J+1)*DZ(K)
VC(-3) = DX(I) *DY(J) *DZ(K-1)
VC( 3) = DX(I) *DY(J) *DZ(K+1)
IF (WALL_INDEX(IC,-1)==0) MASS_N(-1) = RHOP(I-1,J,K)*ABS(MIN(1._EB,MAX(0._EB,ZZP(I-1,J,K,N)))-ZZ_CUT)*VC(-1)
IF (WALL_INDEX(IC, 1)==0) MASS_N( 1) = RHOP(I+1,J,K)*ABS(MIN(1._EB,MAX(0._EB,ZZP(I+1,J,K,N)))-ZZ_CUT)*VC( 1)
IF (WALL_INDEX(IC,-2)==0) MASS_N(-2) = RHOP(I,J-1,K)*ABS(MIN(1._EB,MAX(0._EB,ZZP(I,J-1,K,N)))-ZZ_CUT)*VC(-2)
IF (WALL_INDEX(IC, 2)==0) MASS_N( 2) = RHOP(I,J+1,K)*ABS(MIN(1._EB,MAX(0._EB,ZZP(I,J+1,K,N)))-ZZ_CUT)*VC( 2)
IF (WALL_INDEX(IC,-3)==0) MASS_N(-3) = RHOP(I,J,K-1)*ABS(MIN(1._EB,MAX(0._EB,ZZP(I,J,K-1,N)))-ZZ_CUT)*VC(-3)
IF (WALL_INDEX(IC, 3)==0) MASS_N( 3) = RHOP(I,J,K+1)*ABS(MIN(1._EB,MAX(0._EB,ZZP(I,J,K+1,N)))-ZZ_CUT)*VC( 3)
SUM_MASS_N = SUM(MASS_N)
IF (SUM_MASS_N<=ZERO_P) CYCLE
MASS_C = RHOP(I,J,K)*ABS(ZZP(I,J,K,N)-ZZ_CUT)*VC(0)
CONST = SIGN_FACTOR*MIN(1._EB,MASS_C/SUM_MASS_N)
DELTA_ZZ(I,J,K) = DELTA_ZZ(I,J,K) + CONST*SUM_MASS_N/(RHOP(I,J,K) *VC( 0))
DELTA_ZZ(I-1,J,K) = DELTA_ZZ(I-1,J,K) - CONST*MASS_N(-1)/(RHOP(I-1,J,K)*VC(-1))
DELTA_ZZ(I+1,J,K) = DELTA_ZZ(I+1,J,K) - CONST*MASS_N( 1)/(RHOP(I+1,J,K)*VC( 1))
DELTA_ZZ(I,J-1,K) = DELTA_ZZ(I,J-1,K) - CONST*MASS_N(-2)/(RHOP(I,J-1,K)*VC(-2))
DELTA_ZZ(I,J+1,K) = DELTA_ZZ(I,J+1,K) - CONST*MASS_N( 2)/(RHOP(I,J+1,K)*VC( 2))
DELTA_ZZ(I,J,K-1) = DELTA_ZZ(I,J,K-1) - CONST*MASS_N(-3)/(RHOP(I,J,K-1)*VC(-3))
DELTA_ZZ(I,J,K+1) = DELTA_ZZ(I,J,K+1) - CONST*MASS_N( 3)/(RHOP(I,J,K+1)*VC( 3))
ENDDO
ENDDO
ENDDO
ZZP(1:IBAR,1:JBAR,1:KBAR,N) = ZZP(1:IBAR,1:JBAR,1:KBAR,N) + DELTA_ZZ(1:IBAR,1:JBAR,1:KBAR)
ENDDO SPECIES_LOOP
END SUBROUTINE CHECK_MASS_FRACTION
REAL(EB) FUNCTION SCALAR_FACE_VALUE(A,U,LIMITER)
REAL(EB), INTENT(IN) :: A,U(4)
INTEGER, INTENT(IN) :: LIMITER
REAL(EB) :: R,B,DU_UP,DU_LOC,V(5)
! This function computes the scalar value on a face.
! The scalar is denoted U, and the velocity is denoted A.
! The divergence (computed elsewhere) uses a central difference across
! the cell subject to a flux LIMITER. The flux LIMITER choices are:
!
! LIMITER = 0 implements central differencing
! LIMITER = 1 implements first-order upwinding (monotone)
! LIMITER = 2 implements the SUPERBEE (SB) LIMITER of Roe
! LIMITER = 3 implements the MINMOD LIMITER
! LIMITER = 4 implements the CHARM LIMITER
!
! location of face
!
! f
! | o | o | o | o |
! A
! U(1) U(2) U(3) U(4)
WIND_DIRECTION_IF: IF (A>0._EB) THEN
! the flow is left to right
DU_UP = U(2)-U(1)
DU_LOC = U(3)-U(2)
R = 0._EB
B = 0._EB
SELECT CASE(LIMITER)
CASE(0) ! central differencing
SCALAR_FACE_VALUE = 0.5_EB*(U(2)+U(3))
CASE(1) ! first-order upwinding
SCALAR_FACE_VALUE = U(2)
CASE(2) ! SUPERBEE, Roe (1986)
IF (ABS(DU_LOC)>ZERO_P) R = DU_UP/DU_LOC
B = MAX(0._EB,MIN(2._EB*R,1._EB),MIN(R,2._EB))
SCALAR_FACE_VALUE = U(2) + 0.5_EB*B*DU_LOC
CASE(3) ! MINMOD
IF (ABS(DU_LOC)>ZERO_P) R = DU_UP/DU_LOC
B = MAX(0._EB,MIN(1._EB,R))
SCALAR_FACE_VALUE = U(2) + 0.5_EB*B*DU_LOC
CASE(4) ! CHARM
IF (ABS(DU_UP)>ZERO_P) R = DU_LOC/DU_UP
IF (R>0._EB) B = R*(3._EB*R+1._EB)/((R+1._EB)**2)
SCALAR_FACE_VALUE = U(2) + 0.5_EB*B*DU_UP
CASE(5) ! MP5, Suresh and Huynh (1997)
V = (/2._EB*U(1)-U(2),U(1:4)/)
SCALAR_FACE_VALUE = MP5(V)
END SELECT
ELSE WIND_DIRECTION_IF
! the flow is right to left
DU_UP = U(4)-U(3)
DU_LOC = U(3)-U(2)
R = 0._EB
B = 0._EB
SELECT CASE(LIMITER)
CASE(0) ! central differencing
SCALAR_FACE_VALUE = 0.5_EB*(U(2)+U(3))
CASE(1) ! first-order upwinding
SCALAR_FACE_VALUE = U(3)
CASE(2) ! SUPERBEE, Roe (1986)
IF (ABS(DU_LOC)>ZERO_P) R = DU_UP/DU_LOC
B = MAX(0._EB,MIN(2._EB*R,1._EB),MIN(R,2._EB))
SCALAR_FACE_VALUE = U(3) - 0.5_EB*B*DU_LOC
CASE(3) ! MINMOD
IF (ABS(DU_LOC)>ZERO_P) R = DU_UP/DU_LOC
B = MAX(0._EB,MIN(1._EB,R))
SCALAR_FACE_VALUE = U(3) - 0.5_EB*B*DU_LOC
CASE(4) ! CHARM
IF (ABS(DU_UP)>ZERO_P) R = DU_LOC/DU_UP
IF (R>0._EB) B = R*(3._EB*R+1._EB)/((R+1._EB)**2)
SCALAR_FACE_VALUE = U(3) - 0.5_EB*B*DU_UP
CASE(5) ! MP5, Suresh and Huynh (1997)
V = (/2._EB*U(4)-U(3),U(4),U(3),U(2),U(1)/)
SCALAR_FACE_VALUE = MP5(V)
END SELECT
ENDIF WIND_DIRECTION_IF
END FUNCTION SCALAR_FACE_VALUE
REAL(EB) FUNCTION MP5(V)
USE MATH_FUNCTIONS, ONLY: MINMOD2,MINMOD4
REAL(EB), INTENT(IN) :: V(-2:2)
REAL(EB), PARAMETER :: B1 = 0.016666666666667_EB, B2 = 1.333333333333_EB, ALPHA=4._EB, EPSM=1.E-10_EB
REAL(EB) :: VOR,VMP,DJM1,DJ,DJP1,DM4JPH,DM4JMH,VUL,VAV,VMD,VLC,VMIN,VMAX
! Monotonicity preserving 5th-order scheme (MP5) of Suresh and Huynh, JCP 136, 83-99 (1997)
VOR = B1*(2._EB*V(-2)-13._EB*V(-1)+47._EB*V(0)+27._EB*V(1)-3._EB*V(2))
VMP = V(0) + MINMOD2(V(1)-V(0),ALPHA*(V(0)-V(-1)))
IF ((VOR-V(0))*(VOR-VMP)<EPSM) THEN
MP5=VOR
ELSE
DJM1 = V(-2)-2._EB*V(-1)+V(0)
DJ = V(-1)-2._EB*V(0) +V(1)
DJP1 = V(0) -2._EB*V(1) +V(2)
DM4JPH = MINMOD4(4._EB*DJ-DJP1,4._EB*DJP1-DJ,DJ,DJP1)
DM4JMH = MINMOD4(4._EB*DJ-DJM1,4._EB*DJM1-DJ,DJ,DJM1)
VUL = V(0) + ALPHA*(V(0)-V(-1))
VAV = 0.5_EB*(V(0)+V(1))
VMD = VAV - 0.5_EB*DM4JPH
VLC = V(0) + 0.5_EB*(V(0)-V(-1)) + B2*DM4JMH
VMIN = MAX(MIN(V(0),V(1),VMD),MIN(V(0),VUL,VLC))
VMAX = MIN(MAX(V(0),V(1),VMD),MAX(V(0),VUL,VLC))
MP5 = VOR + MINMOD2(VMIN-VOR,VMAX-VOR)
ENDIF
END FUNCTION MP5
SUBROUTINE GET_REV_mass(MODULE_REV,MODULE_DATE)
INTEGER,INTENT(INOUT) :: MODULE_REV
CHARACTER(255),INTENT(INOUT) :: MODULE_DATE
WRITE(MODULE_DATE,'(A)') massrev(INDEX(massrev,':')+1:LEN_TRIM(massrev)-2)
READ (MODULE_DATE,'(I5)') MODULE_REV
WRITE(MODULE_DATE,'(A)') massdate
END SUBROUTINE GET_REV_mass
END MODULE MASS