-
Notifications
You must be signed in to change notification settings - Fork 0
/
turb.f90
2523 lines (2231 loc) · 71.6 KB
/
turb.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
! This module is useful for verification tests and development of
! turbulence models.
MODULE TURBULENCE
USE PRECISION_PARAMETERS
USE GLOBAL_CONSTANTS
USE MESH_POINTERS
USE MESH_VARIABLES
USE COMP_FUNCTIONS
IMPLICIT NONE
CHARACTER(255), PARAMETER :: turbid='$Id: turb.f90 9962 2012-02-02 15:53:17Z randy.mcdermott $'
CHARACTER(255), PARAMETER :: turbrev='$Revision: 9962 $'
CHARACTER(255), PARAMETER :: turbdate='$Date: 2012-02-02 07:53:17 -0800 (Thu, 02 Feb 2012) $'
PRIVATE
PUBLIC :: NS_ANALYTICAL_SOLUTION, INIT_TURB_ARRAYS, VARDEN_DYNSMAG, &
GET_REV_turb, WERNER_WENGLE_WALL_MODEL, COMPRESSION_WAVE, VELTAN2D,VELTAN3D,STRATIFIED_MIXING_LAYER, &
SURFACE_HEAT_FLUX_MODEL, SYNTHETIC_TURBULENCE, SYNTHETIC_EDDY_SETUP, TEST_FILTER, EX2G3D, TENSOR_DIFFUSIVITY_MODEL
CONTAINS
SUBROUTINE INIT_TURB_ARRAYS(NM)
USE MEMORY_FUNCTIONS, ONLY: ChkMemErr
INTEGER, INTENT(IN) :: NM
INTEGER :: IZERO
TYPE (MESH_TYPE), POINTER :: M
CALL POINT_TO_MESH(NM)
M => MESHES(NM)
!IF (PERIODIC_TEST==2 .OR. TURB_MODEL==DYNSMAG) THEN
! ! real work arrays
! ALLOCATE(M%TURB_WORK1(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
! CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK1',IZERO)
! M%TURB_WORK1 = 0._EB
! ALLOCATE(M%TURB_WORK2(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
! CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK2',IZERO)
! M%TURB_WORK2 = 0._EB
! ALLOCATE(M%TURB_WORK3(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
! CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK3',IZERO)
! M%TURB_WORK3 = 0._EB
! ALLOCATE(M%TURB_WORK4(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
! CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK4',IZERO)
! M%TURB_WORK4 = 0._EB
!ENDIF
IF (TURB_MODEL==DYNSMAG) THEN
! real work arrays
ALLOCATE(M%TURB_WORK5(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK5',IZERO)
M%TURB_WORK5 = 0._EB
ALLOCATE(M%TURB_WORK6(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK6',IZERO)
M%TURB_WORK6 = 0._EB
ALLOCATE(M%TURB_WORK7(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK7',IZERO)
M%TURB_WORK7 = 0._EB
ALLOCATE(M%TURB_WORK8(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK8',IZERO)
M%TURB_WORK8 = 0._EB
ALLOCATE(M%TURB_WORK9(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK9',IZERO)
M%TURB_WORK9 = 0._EB
ALLOCATE(M%TURB_WORK10(0:IBP1,0:JBP1,0:KBP1),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK10',IZERO)
M%TURB_WORK10 = 0._EB
ENDIF
! 1D working arrays
IF (LES) THEN
ALLOCATE(M%TURB_WORK11(0:MAX(IBP1,JBP1,KBP1)),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK11',IZERO)
M%TURB_WORK11 = 0._EB
ALLOCATE(M%TURB_WORK12(0:MAX(IBP1,JBP1,KBP1)),STAT=IZERO)
CALL ChkMemErr('INIT_TURB_ARRAYS','TURB_WORK12',IZERO)
M%TURB_WORK12 = 0._EB
ENDIF
END SUBROUTINE INIT_TURB_ARRAYS
SUBROUTINE NS_ANALYTICAL_SOLUTION(NM)
! Initialize flow variables with an analytical solution of the governing equations
INTEGER, INTENT(IN) :: NM
INTEGER :: I,J,K
REAL(EB) :: UU,WW
CALL POINT_TO_MESH(NM)
DO K=1,KBAR
DO J=1,JBAR
DO I=0,IBAR
U(I,J,K) = 1._EB - 2._EB*COS(X(I))*SIN(ZC(K))
ENDDO
ENDDO
ENDDO
DO K=0,KBAR
DO J=1,JBAR
DO I=1,IBAR
W(I,J,K) = 1._EB + 2._EB*SIN(XC(I))*COS(Z(K))
ENDDO
ENDDO
ENDDO
DO K=0,KBP1
DO J=0,JBP1
DO I=0,IBP1
UU = 1._EB - 2._EB*COS(XC(I))*SIN(ZC(K))
WW = 1._EB + 2._EB*SIN(XC(I))*COS(ZC(K))
H(I,J,K) = -( COS(2._EB*XC(I)) + COS(2._EB*ZC(K)) ) + 0.5_EB*(UU**2+WW**2)
ENDDO
ENDDO
ENDDO
END SUBROUTINE NS_ANALYTICAL_SOLUTION
SUBROUTINE COMPRESSION_WAVE(NM,T,ITEST)
INTEGER, INTENT(IN) :: NM,ITEST
REAL(EB), INTENT(IN) :: T
INTEGER :: I,J,K
CALL POINT_TO_MESH(NM)
SELECT CASE(ITEST)
CASE(3) ! stationary compression wave
DO K=1,KBAR
DO J=1,JBAR
DO I=0,IBAR
U(I,J,K) = 2._EB + SIN(X(I))
US(I,J,K) = U(I,J,K)
ENDDO
ENDDO
ENDDO
DO K=1,KBAR
DO J=0,JBAR
DO I=1,IBAR
V(I,J,K) = 0._EB
VS(I,J,K) = 0._EB
ENDDO
ENDDO
ENDDO
DO K=0,KBAR
DO J=1,JBAR
DO I=1,IBAR
W(I,J,K) = 3._EB + SIN(Z(K))
WS(I,J,K) = W(I,J,K)
ENDDO
ENDDO
ENDDO
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
D(I,J,K) = (U(I,J,K)-U(I-1,J,K))*RDX(I) + (W(I,J,K)-W(I,J,K-1))*RDZ(K)
DS(I,J,K) = D(I,J,K)
ENDDO
ENDDO
ENDDO
CASE(4) ! pulsating dilation
PREDICTOR_IF: IF (PREDICTOR) THEN
DO K=1,KBAR
DO J=1,JBAR
DO I=0,IBAR
US(I,J,K) = SIN(X(I))*COS(T)
ENDDO
ENDDO
ENDDO
DO K=1,KBAR
DO J=0,JBAR
DO I=1,IBAR
VS(I,J,K) = 0._EB
ENDDO
ENDDO
ENDDO
DO K=0,KBAR
DO J=1,JBAR
DO I=1,IBAR
WS(I,J,K) = SIN(Z(K))*COS(T)
ENDDO
ENDDO
ENDDO
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
DS(I,J,K) = (US(I,J,K)-US(I-1,J,K))*RDX(I) + (WS(I,J,K)-WS(I,J,K-1))*RDZ(K)
ENDDO
ENDDO
ENDDO
ELSE PREDICTOR_IF
DO K=1,KBAR
DO J=1,JBAR
DO I=0,IBAR
U(I,J,K) = SIN(X(I))*COS(T)
ENDDO
ENDDO
ENDDO
DO K=1,KBAR
DO J=0,JBAR
DO I=1,IBAR
V(I,J,K) = 0._EB
ENDDO
ENDDO
ENDDO
DO K=0,KBAR
DO J=1,JBAR
DO I=1,IBAR
W(I,J,K) = SIN(Z(K))*COS(T)
ENDDO
ENDDO
ENDDO
DO K=1,KBAR
DO J=1,JBAR
DO I=1,IBAR
D(I,J,K) = (U(I,J,K)-U(I-1,J,K))*RDX(I) + (W(I,J,K)-W(I,J,K-1))*RDZ(K)
ENDDO
ENDDO
ENDDO
ENDIF PREDICTOR_IF
END SELECT
END SUBROUTINE COMPRESSION_WAVE
SUBROUTINE STRATIFIED_MIXING_LAYER(NM)
USE PHYSICAL_FUNCTIONS, ONLY: GET_SPECIFIC_GAS_CONSTANT
INTEGER, INTENT(IN) :: NM
INTEGER :: I,J,K
REAL(EB) :: ZZ_GET(0:N_TRACKED_SPECIES)
CALL POINT_TO_MESH(NM)
DO K=0,KBP1
U(:,:,K) = 0.5_EB*MIXING_LAYER_U0* TANH(2._EB*ZC(K)/MIXING_LAYER_H0)
RHO(:,:,K) = RHOA*(1._EB - 0.5_EB*MIXING_LAYER_THETA0*TANH(2._EB*ZC(K)/MIXING_LAYER_H0))
ENDDO
! initialize temperature
DO K=0,KBP1
DO J=0,JBP1
DO I=0,IBP1
IF (N_TRACKED_SPECIES>0) 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))
TMP(I,J,K) = PBAR(K,PRESSURE_ZONE(I,J,K))/(RSUM(I,J,K)*RHO(I,J,K))
ENDDO
ENDDO
ENDDO
END SUBROUTINE STRATIFIED_MIXING_LAYER
SUBROUTINE VARDEN_DYNSMAG(NM)
INTEGER, INTENT(IN) :: NM
REAL(EB) :: TEMP_TERM,DUDY,DUDZ,DVDX,DVDZ,DWDX,DWDY,ONTHDIV
INTEGER :: I,J,K
REAL(EB), POINTER, DIMENSION(:,:,:) :: UU,VV,WW,UP,VP,WP,RHOP,RHOPHAT
REAL(EB), POINTER, DIMENSION(:,:,:) :: S11,S22,S33,S12,S13,S23,SS
REAL(EB), POINTER, DIMENSION(:,:,:) :: SHAT11,SHAT22,SHAT33,SHAT12,SHAT13,SHAT23,SSHAT
REAL(EB), POINTER, DIMENSION(:,:,:) :: BETA11,BETA22,BETA33,BETA12,BETA13,BETA23
REAL(EB), POINTER, DIMENSION(:,:,:) :: BETAHAT11,BETAHAT22,BETAHAT33,BETAHAT12,BETAHAT13,BETAHAT23
REAL(EB), POINTER, DIMENSION(:,:,:) :: M11,M22,M33,M12,M13,M23,MM,MMHAT
REAL(EB), POINTER, DIMENSION(:,:,:) :: L11,L22,L33,L12,L13,L23,ML,MLHAT
REAL(EB), PARAMETER :: ALPHA = 6.0_EB ! See Lund, 1997 CTR briefs.
! References:
!
! M. Germano, U. Piomelli, P. Moin, and W. Cabot. A dynamic subgrid-scale eddy viscosity model.
! Phys. Fluids A, 3(7):1760-1765, 1991.
!
! M. Pino Martin, U. Piomelli, and G. Candler. Subgrid-scale models for compressible large-eddy
! simulation. Theoret. Comput. Fluid Dynamics, 13:361-376, 2000.
!
! P. Moin, K. Squires, W. Cabot, and S. Lee. A dynamic subgrid-scale model for compressible
! turbulence and scalar transport. Phys. Fluids A, 3(11):2746-2757, 1991.
!
! T. S. Lund. On the use of discrete filters for large eddy simulation. Center for Turbulence
! Research Annual Research Briefs, 1997.
!
! R. McDermott. Variable density formulation of the dynamic Smagorinsky model.
! http://randy.mcdermott.googlepages.com/dynsmag_comp.pdf
! *****************************************************************************
! CAUTION WHEN MODIFYING: The order in which the tensor components are computed
! is important because we overwrite pointers several times to conserve memory.
! *****************************************************************************
IF (EVACUATION_ONLY(NM)) RETURN
CALL POINT_TO_MESH(NM)
IF (PREDICTOR) THEN
UU=>U
VV=>V
WW=>W
RHOP=>RHO
ELSE
UU=>US
VV=>VS
WW=>WS
RHOP=>RHOS
ENDIF
UP => TURB_WORK1
VP => TURB_WORK2
WP => TURB_WORK3
S11 => WORK1
S22 => WORK2
S33 => WORK3
S12 => WORK4
S13 => WORK5
S23 => WORK6
SS => WORK7
DO K = 1,KBAR
DO J = 1,JBAR
DO I = 1,IBAR
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))
S11(I,J,K) = RDX(I)*(UU(I,J,K)-UU(I-1,J,K))
S22(I,J,K) = RDY(J)*(VV(I,J,K)-VV(I,J-1,K))
S33(I,J,K) = RDZ(K)*(WW(I,J,K)-WW(I,J,K-1))
ONTHDIV = ONTH*(S11(I,J,K)+S22(I,J,K)+S33(I,J,K))
S11(I,J,K) = S11(I,J,K)-ONTHDIV
S22(I,J,K) = S22(I,J,K)-ONTHDIV
S33(I,J,K) = S33(I,J,K)-ONTHDIV
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))
S12(I,J,K) = 0.5_EB*(DUDY+DVDX)
S13(I,J,K) = 0.5_EB*(DUDZ+DWDX)
S23(I,J,K) = 0.5_EB*(DVDZ+DWDY)
! calculate magnitude of the grid strain rate
SS(I,J,K) = RHOP(I,J,K)*STRAIN_RATE(I,J,K)
ENDDO
ENDDO
ENDDO
! second-order extrapolation to ghost cells
CALL EX2G3D(UP,-1.E10_EB,1.E10_EB)
CALL EX2G3D(VP,-1.E10_EB,1.E10_EB)
CALL EX2G3D(WP,-1.E10_EB,1.E10_EB)
CALL EX2G3D(S11,-1.E10_EB,1.E10_EB)
CALL EX2G3D(S22,-1.E10_EB,1.E10_EB)
CALL EX2G3D(S33,-1.E10_EB,1.E10_EB)
CALL EX2G3D(S12,-1.E10_EB,1.E10_EB)
CALL EX2G3D(S13,-1.E10_EB,1.E10_EB)
CALL EX2G3D(S23,-1.E10_EB,1.E10_EB)
CALL EX2G3D(SS,0._EB,1.E10_EB)
! test filter the strain rate
SHAT11 => TURB_WORK4
SHAT22 => TURB_WORK5
SHAT33 => TURB_WORK6
SHAT12 => TURB_WORK7
SHAT13 => TURB_WORK8
SHAT23 => TURB_WORK9
SSHAT => TURB_WORK10
CALL TEST_FILTER(SHAT11,S11,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(SHAT22,S22,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(SHAT33,S33,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(SHAT12,S12,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(SHAT13,S13,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(SHAT23,S23,-1.E10_EB,1.E10_EB)
! calculate magnitude of test filtered strain rate
DO K = 1,KBAR
DO J = 1,JBAR
DO I = 1,IBAR
SSHAT(I,J,K) = SQRT(2._EB*(SHAT11(I,J,K)*SHAT11(I,J,K) + &
SHAT22(I,J,K)*SHAT22(I,J,K) + &
SHAT33(I,J,K)*SHAT33(I,J,K) + &
2._EB*(SHAT12(I,J,K)*SHAT12(I,J,K) + &
SHAT13(I,J,K)*SHAT13(I,J,K) + &
SHAT23(I,J,K)*SHAT23(I,J,K)) ) )
ENDDO
ENDDO
ENDDO
! calculate the grid filtered stress tensor, beta
BETA11 => WORK1
BETA22 => WORK2
BETA33 => WORK3
BETA12 => WORK4
BETA13 => WORK5
BETA23 => WORK6
BETA11 = SS*S11
BETA22 = SS*S22
BETA33 = SS*S33
BETA12 = SS*S12
BETA13 = SS*S13
BETA23 = SS*S23
! ghost values for beta_ij should be filled already
! test filter the grid filtered stress tensor
BETAHAT11 => WORK1
BETAHAT22 => WORK2
BETAHAT33 => WORK3
BETAHAT12 => WORK4
BETAHAT13 => WORK5
BETAHAT23 => WORK6
CALL TEST_FILTER(BETAHAT11,BETA11,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(BETAHAT22,BETA22,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(BETAHAT33,BETA33,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(BETAHAT12,BETA12,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(BETAHAT13,BETA13,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(BETAHAT23,BETA23,-1.E10_EB,1.E10_EB)
! test filter the density
RHOPHAT => WORK7
CALL TEST_FILTER(RHOPHAT,RHOP,RHOMIN,RHOMAX)
! calculate the Mij tensor
M11 => WORK1
M22 => WORK2
M33 => WORK3
M12 => WORK4
M13 => WORK5
M23 => WORK6
DO K = 1,KBAR
DO J = 1,JBAR
DO I = 1,IBAR
TEMP_TERM = ALPHA*RHOPHAT(I,J,K)*SSHAT(I,J,K)
M11(I,J,K) = 2._EB*(BETAHAT11(I,J,K) - TEMP_TERM*SHAT11(I,J,K))
M22(I,J,K) = 2._EB*(BETAHAT22(I,J,K) - TEMP_TERM*SHAT22(I,J,K))
M33(I,J,K) = 2._EB*(BETAHAT33(I,J,K) - TEMP_TERM*SHAT33(I,J,K))
M12(I,J,K) = 2._EB*(BETAHAT12(I,J,K) - TEMP_TERM*SHAT12(I,J,K))
M13(I,J,K) = 2._EB*(BETAHAT13(I,J,K) - TEMP_TERM*SHAT13(I,J,K))
M23(I,J,K) = 2._EB*(BETAHAT23(I,J,K) - TEMP_TERM*SHAT23(I,J,K))
ENDDO
ENDDO
ENDDO
! calculate the Leonard term, Lij
L11 => TURB_WORK4
L22 => TURB_WORK5
L33 => TURB_WORK6
L12 => TURB_WORK7
L13 => TURB_WORK8
L23 => TURB_WORK9
CALL CALC_VARDEN_LEONARD_TERM
! calculate Mij*Lij & Mij*Mij
MM => TURB_WORK1
MMHAT => TURB_WORK1
ML => TURB_WORK2
MLHAT => TURB_WORK2
DO K = 1,KBAR
DO J = 1,JBAR
DO I = 1,IBAR
ML(I,J,K) = M11(I,J,K)*L11(I,J,K) + M22(I,J,K)*L22(I,J,K) + M33(I,J,K)*L33(I,J,K) + &
2._EB*(M12(I,J,K)*L12(I,J,K) + M13(I,J,K)*L13(I,J,K) + M23(I,J,K)*L23(I,J,K))
MM(I,J,K) = M11(I,J,K)*M11(I,J,K) + M22(I,J,K)*M22(I,J,K) + M33(I,J,K)*M33(I,J,K) + &
2._EB*(M12(I,J,K)*M12(I,J,K) + M13(I,J,K)*M13(I,J,K) + M23(I,J,K)*M23(I,J,K))
ENDDO
ENDDO
ENDDO
! extrapolate to ghost
CALL EX2G3D(ML,0._EB,1.E10_EB)
CALL EX2G3D(MM,0._EB,1.E10_EB)
! do some smoothing
CALL TEST_FILTER(MLHAT,ML,0._EB,1.E10_EB)
CALL TEST_FILTER(MMHAT,MM,0._EB,1.E10_EB)
DO K = 1,KBAR
DO J = 1,JBAR
DO I = 1,IBAR
! calculate the local Smagorinsky coefficient
! perform "clipping" in case MLij is negative
IF (MLHAT(I,J,K) < 0._EB) MLHAT(I,J,K) = 0._EB
! calculate the effective viscosity
! handle the case where we divide by zero, note MMHAT is positive semi-definite
IF (MMHAT(I,J,K)<=ZERO_P) THEN
CSD2(I,J,K) = 0._EB
ELSE
CSD2(I,J,K) = MLHAT(I,J,K)/MMHAT(I,J,K) ! (Cs*Delta)**2
ENDIF
END DO
END DO
END DO
END SUBROUTINE VARDEN_DYNSMAG
SUBROUTINE CALC_VARDEN_LEONARD_TERM
REAL(EB), POINTER, DIMENSION(:,:,:) :: L11,L22,L33,L12,L13,L23
REAL(EB), POINTER, DIMENSION(:,:,:) :: RHOP,RHOPHAT
REAL(EB), POINTER, DIMENSION(:,:,:) :: UP,VP,WP
REAL(EB), POINTER, DIMENSION(:,:,:) :: RUU,RVV,RWW,RUV,RUW,RVW
REAL(EB), POINTER, DIMENSION(:,:,:) :: RU,RV,RW
REAL(EB), POINTER, DIMENSION(:,:,:) :: RUU_HAT,RVV_HAT,RWW_HAT,RUV_HAT,RUW_HAT,RVW_HAT
REAL(EB), POINTER, DIMENSION(:,:,:) :: RU_HAT,RV_HAT,RW_HAT
REAL(EB) :: INV_RHOPHAT
INTEGER :: I,J,K
! *****************************************************************************
! CAUTION WHEN MODIFYING: The order in which the tensor components are computed
! is important because we overwrite pointers several times to conserve memory.
! *****************************************************************************
IF (PREDICTOR) THEN
RHOP=>RHO
ELSE
RHOP=>RHOS
ENDIF
RHOPHAT => WORK7
! Compute rho*UiUj
UP => TURB_WORK1 ! will be overwritten by RU
VP => TURB_WORK2
WP => TURB_WORK3
RUU => TURB_WORK4 ! will be overwritten by RUU_HAT
RVV => TURB_WORK5
RWW => TURB_WORK6
RUV => TURB_WORK7
RUW => TURB_WORK8
RVW => TURB_WORK9
RUU = RHOP*UP*UP
RVV = RHOP*VP*VP
RWW = RHOP*WP*WP
RUV = RHOP*UP*VP
RUW = RHOP*UP*WP
RVW = RHOP*VP*WP
! extrapolate to ghost cells
CALL EX2G3D(RUU,0.E10_EB,1.E10_EB)
CALL EX2G3D(RVV,0.E10_EB,1.E10_EB)
CALL EX2G3D(RWW,0.E10_EB,1.E10_EB)
CALL EX2G3D(RUV,-1.E10_EB,1.E10_EB)
CALL EX2G3D(RUW,-1.E10_EB,1.E10_EB)
CALL EX2G3D(RVW,-1.E10_EB,1.E10_EB)
! Test filter rho*UiUj
RUU_HAT => TURB_WORK4 ! will be overwritten by Lij
RVV_HAT => TURB_WORK5
RWW_HAT => TURB_WORK6
RUV_HAT => TURB_WORK7
RUW_HAT => TURB_WORK8
RVW_HAT => TURB_WORK9
CALL TEST_FILTER(RUU_HAT,RUU,0.E10_EB,1.E10_EB)
CALL TEST_FILTER(RVV_HAT,RVV,0.E10_EB,1.E10_EB)
CALL TEST_FILTER(RWW_HAT,RWW,0.E10_EB,1.E10_EB)
CALL TEST_FILTER(RUV_HAT,RUV,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(RUW_HAT,RUW,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(RVW_HAT,RVW,-1.E10_EB,1.E10_EB)
! Compute rho*Ui
RU => TURB_WORK1 ! will be overwritten by RU_HAT
RV => TURB_WORK2
RW => TURB_WORK3
RU = RHOP*UP
RV = RHOP*VP
RW = RHOP*WP
! extrapolate to ghost cells
CALL EX2G3D(RU,-1.E10_EB,1.E10_EB)
CALL EX2G3D(RV,-1.E10_EB,1.E10_EB)
CALL EX2G3D(RW,-1.E10_EB,1.E10_EB)
! Test filter rho*Ui
RU_HAT => TURB_WORK1
RV_HAT => TURB_WORK2
RW_HAT => TURB_WORK3
CALL TEST_FILTER(RU_HAT,RU,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(RV_HAT,RV,-1.E10_EB,1.E10_EB)
CALL TEST_FILTER(RW_HAT,RW,-1.E10_EB,1.E10_EB)
! Compute variable density Leonard stress
L11 => TURB_WORK4
L22 => TURB_WORK5
L33 => TURB_WORK6
L12 => TURB_WORK7
L13 => TURB_WORK8
L23 => TURB_WORK9
DO K = 1,KBAR
DO J = 1,JBAR
DO I = 1,IBAR
INV_RHOPHAT = 1._EB/RHOPHAT(I,J,K)
L11(I,J,K) = RUU_HAT(I,J,K) - RU_HAT(I,J,K)*RU_HAT(I,J,K)*INV_RHOPHAT
L22(I,J,K) = RVV_HAT(I,J,K) - RV_HAT(I,J,K)*RV_HAT(I,J,K)*INV_RHOPHAT
L33(I,J,K) = RWW_HAT(I,J,K) - RW_HAT(I,J,K)*RW_HAT(I,J,K)*INV_RHOPHAT
L12(I,J,K) = RUV_HAT(I,J,K) - RU_HAT(I,J,K)*RV_HAT(I,J,K)*INV_RHOPHAT
L13(I,J,K) = RUW_HAT(I,J,K) - RU_HAT(I,J,K)*RW_HAT(I,J,K)*INV_RHOPHAT
L23(I,J,K) = RVW_HAT(I,J,K) - RV_HAT(I,J,K)*RW_HAT(I,J,K)*INV_RHOPHAT
ENDDO
ENDDO
ENDDO
END SUBROUTINE CALC_VARDEN_LEONARD_TERM
SUBROUTINE EX2G3D(A,A_MIN,A_MAX)
! Second order extrapolation of 3D array to ghost cells
REAL(EB), INTENT(IN) :: A_MIN,A_MAX
REAL(EB), INTENT(INOUT) :: A(0:IBP1,0:JBP1,0:KBP1)
A(0,:,:) = MIN(A_MAX,MAX(A_MIN,2._EB*A(1,:,:)-A(2,:,:)))
A(:,0,:) = MIN(A_MAX,MAX(A_MIN,2._EB*A(:,1,:)-A(:,2,:)))
A(:,:,0) = MIN(A_MAX,MAX(A_MIN,2._EB*A(:,:,1)-A(:,:,2)))
A(IBP1,:,:) = MIN(A_MAX,MAX(A_MIN,2._EB*A(IBAR,:,:)-A(IBM1,:,:)))
A(:,JBP1,:) = MIN(A_MAX,MAX(A_MIN,2._EB*A(:,JBAR,:)-A(:,JBM1,:)))
A(:,:,KBP1) = MIN(A_MAX,MAX(A_MIN,2._EB*A(:,:,KBAR)-A(:,:,KBM1)))
END SUBROUTINE EX2G3D
SUBROUTINE TEST_FILTER(PHIBAR,PHI,PHI_MIN,PHI_MAX)
REAL(EB), INTENT(IN) :: PHI_MIN,PHI_MAX
REAL(EB), INTENT(IN) :: PHI(0:IBP1,0:JBP1,0:KBP1)
REAL(EB), INTENT(OUT) :: PHIBAR(0:IBP1,0:JBP1,0:KBP1)
REAL(EB), POINTER, DIMENSION(:) :: PHI1,PHI2
INTEGER I,J,K
PHI1 => TURB_WORK11
PHI2 => TURB_WORK12
! filter in x:
DO K = 0,KBP1
DO J = 0,JBP1
PHI1(0:IBP1) = PHI(0:IBP1,J,K)
CALL TOPHAT_FILTER_1D(PHI2(0:IBP1),PHI1(0:IBP1),0,IBP1,PHI_MIN,PHI_MAX)
PHIBAR(0:IBP1,J,K) = PHI2(0:IBP1)
ENDDO
ENDDO
IF (.NOT.TWO_D) THEN
! filter in y:
DO K = 0,KBP1
DO I = 0,IBP1
PHI1(0:JBP1) = PHIBAR(I,0:JBP1,K)
CALL TOPHAT_FILTER_1D(PHI2(0:JBP1),PHI1(0:JBP1),0,JBP1,PHI_MIN,PHI_MAX)
PHIBAR(I,0:JBP1,K) = PHI2(0:JBP1)
ENDDO
ENDDO
ENDIF
! filter in z:
DO J = 0,JBP1
DO I = 0,IBP1
PHI1(0:KBP1) = PHIBAR(I,J,0:KBP1)
CALL TOPHAT_FILTER_1D(PHI2(0:KBP1),PHI1(0:KBP1),0,KBP1,PHI_MIN,PHI_MAX)
PHIBAR(I,J,0:KBP1) = PHI2(0:KBP1)
ENDDO
ENDDO
END SUBROUTINE TEST_FILTER
SUBROUTINE TOPHAT_FILTER_1D(UBAR,U,N_LO,N_HI,U_MIN,U_MAX)
INTEGER, INTENT(IN) :: N_LO,N_HI
REAL(EB), INTENT(IN) :: U(N_LO:N_HI),U_MIN,U_MAX
REAL(EB), INTENT(OUT) :: UBAR(N_LO:N_HI)
INTEGER :: J
!REAL(EB),PARAMETER:: W(-1:1) = (/0.25_EB,0.5_EB,0.25_EB/) ! trapezoid rule
!REAL(EB),PARAMETER:: W(-1:1) = (/ONSI,TWTH,ONSI/) ! Simpson's rule
! Filter the u field to obtain ubar
DO J=N_LO+1,N_HI-1
!UBAR(J) = DOT_PRODUCT(W(-1:1),U(J-1:J+1))
UBAR(J) = 0.5_EB*U(J) + 0.25_EB*(U(J-1)+U(J+1))
ENDDO
! set boundary values (not ideal, but fast and simple)
!UBAR(N_LO) = UBAR(N_LO+1)
!UBAR(N_HI) = UBAR(N_HI-1)
UBAR(N_LO) = MIN(U_MAX,MAX(U_MIN,2._EB*UBAR(N_LO+1)-UBAR(N_LO+2)))
UBAR(N_HI) = MIN(U_MAX,MAX(U_MIN,2._EB*UBAR(N_HI-1)-UBAR(N_HI-2)))
END SUBROUTINE TOPHAT_FILTER_1D
SUBROUTINE WERNER_WENGLE_WALL_MODEL(SF,U_TAU,U1,NU,DZ,ROUGHNESS)
REAL(EB), INTENT(OUT) :: SF
REAL(EB), INTENT(IN) :: U1,NU,DZ,ROUGHNESS
REAL(EB), PARAMETER :: A=8.3_EB,B=1._EB/7._EB
REAL(EB), PARAMETER :: Z_PLUS_TURBULENT = 11.81_EB
REAL(EB), PARAMETER :: ALPHA=7.202125273562269_EB !! ALPHA=(1._EB-B)/2._EB*A**((1._EB+B)/(1._EB-B))
REAL(EB), PARAMETER :: BETA=1._EB+B
REAL(EB), PARAMETER :: ETA=(1._EB+B)/A
REAL(EB), PARAMETER :: GAMMA=2._EB/(1._EB+B)
REAL(EB), PARAMETER :: RKAPPA=2.44_EB ! 1./von Karman constant
REAL(EB), PARAMETER :: BTILDE=8.50_EB ! see Pope p. 297
REAL(EB) :: U_TAU,TAU_W,NUODZ,Z_PLUS,TAU_ROUGH
! References (for smooth walls):
!
! Werner, H., Wengle, H. (1991) Large-eddy simulation of turbulent flow over
! and around a cube in a plate channel. 8th Symposium on Turbulent Shear
! Flows, Munich, Germany.
!
! Pierre Sagaut. Large Eddy Simulation for Incompressible Flows: An Introduction.
! Springer, 2001.
!
! Temmerman, L., Leschziner, M.A., Mellen, C.P., and Frohlich, J. (2003)
! Investigation of wall-function approximations and subgrid-scale models in
! Large Eddy Simulation of separated flow in a channel with streamwise
! periodic constrictions. International Journal of Heat and Fluid Flow,
! Vol. 24, No. 2, pp. 157-180.
!
! Breuer, M., Kniazev, B., and Abel, M. (2007) Development of wall models
! for LES of separated flows using statistical evaluations. Computers and
! Fluids, Vol. 36, pp. 817-837.
!
! McDermott, R. (2009) FDS Wall Flows, Part I: Straight Channels, NIST Technical Note.
!
! References (for rough surfaces):
!
! S. B. Pope (2000) Turbulent Flows, Cambridge.
!
! Moeng, C.-H. (1984) A Large-Eddy Simulation Model for the Study of Planetary
! Boundary-Layer Turbulence. Journal of the Atmospheric Sciences, Vol. 41, No. 13,
! pp. 2052-2062.
!
! Stoll, R., Porte-Agel, F. (2008) Large-Eddy Simulation of the Stable Atmospheric
! Boundary Layer using Dynamic Models with Different Averaging Schemes. Boundary-Layer
! Meteorology, 126:1-28.
!
! Comments:
!
! The slip factor (SF) is based on the following approximation to the wall stress
! (note that u0 is the ghost cell value of the streamwise velocity component and
! z is the wall-normal direction):
! tau_w = mu*(u1-u0)/dz = mu*(u1-SF*u1)/dz = mu*u1/dz*(1-SF)
! note that tau_w/rho = nu*u1/dz*(1-SF)
TAU_ROUGH = 0._EB
IF (ROUGHNESS>0._EB) THEN
! Pope (2000)
TAU_ROUGH = ( U1/(RKAPPA*LOG(0.5_EB*DZ/ROUGHNESS)+BTILDE) )**2 ! actually tau_w/rho
ENDIF
! Werner-Wengle
NUODZ = NU/DZ
TAU_W = (ALPHA*(NUODZ)**BETA + ETA*(NUODZ)**B*ABS(U1))**GAMMA ! actually tau_w/rho
TAU_W = MAX(TAU_W,TAU_ROUGH)
U_TAU = SQRT(TAU_W)
Z_PLUS = DZ/(NU/(U_TAU+1.E-10_EB))
IF (Z_PLUS>Z_PLUS_TURBULENT) THEN
SF = 1._EB-TAU_W/(NUODZ*ABS(U1)) ! log layer
ELSE
SF = -1._EB ! viscous sublayer
ENDIF
!! check values...
!IF (Z_PLUS>Z_PLUS_TURBULENT) THEN
! print *,'A = ',A
! print *,'B = ',B
! print *,'ALPHA = ',ALPHA
! print *,'BETA = ',BETA
! print *,'ETA = ',ETA
! print *,'GAMMA = ',GAMMA
! print *,'U1 = ',U1
! print *,'NU/DZ = ',NU_OVER_DZ
! print *,'TAU_W/RHO = ',TAU_W
! print *,'Z_PLUS = ',Z_PLUS
! print *,'SF = ',SF
! print *
!ENDIF
END SUBROUTINE WERNER_WENGLE_WALL_MODEL
REAL(EB) FUNCTION VELTAN2D(U_VELO,U_SURF,NN,DN,DIVU,GRADU,GRADP,TAU_IJ,DT,RRHO,MU,I_VEL)
REAL(EB), INTENT(IN) :: U_VELO(2),U_SURF(2),NN(2),DN,DIVU,GRADU(2,2),GRADP(2),TAU_IJ(2,2),DT,RRHO,MU
INTEGER, INTENT(IN) :: I_VEL
REAL(EB) :: C(2,2),SS(2),SLIP_COEF,ETA,AA,BB,U_STRM_0,DUMMY, &
U_STRM,U_NORM,U_STRM_WALL,U_NORM_WALL,DPDS,DUSDS,DUSDN,TSN
INTEGER :: SUBIT
! Cartesian grid coordinate system orthonormal basis vectors
REAL(EB), DIMENSION(2), PARAMETER :: XX=(/1._EB, 0._EB/),YY=(/0._EB, 1._EB/)
! streamwise unit vector
SS = (/NN(2),-NN(1)/)
! directional cosines (see Pope, Eq. A.11)
C(1,1) = DOT_PRODUCT(XX,SS)
C(1,2) = DOT_PRODUCT(XX,NN)
C(2,1) = DOT_PRODUCT(YY,SS)
C(2,2) = DOT_PRODUCT(YY,NN)
! transform velocity (see Pope, Eq. A.17)
U_STRM = C(1,1)*U_VELO(1) + C(2,1)*U_VELO(2)
U_NORM = C(1,2)*U_VELO(1) + C(2,2)*U_VELO(2)
! transform wall velocity
U_STRM_WALL = C(1,1)*U_SURF(1) + C(2,1)*U_SURF(2)
U_NORM_WALL = C(1,2)*U_SURF(1) + C(2,2)*U_SURF(2)
! transform pressure gradient
DPDS = C(1,1)*GRADP(1) + C(2,1)*GRADP(2)
! transform velocity gradient tensor (Pope A.23)
DUSDS = C(1,1)*C(1,1)*GRADU(1,1) + C(1,1)*C(2,1)*GRADU(1,2) &
+ C(2,1)*C(1,1)*GRADU(2,1) + C(2,1)*C(2,1)*GRADU(2,2)
DUSDN = C(1,1)*C(1,2)*GRADU(1,1) + C(1,1)*C(2,2)*GRADU(1,2) &
+ C(2,1)*C(1,2)*GRADU(2,1) + C(2,1)*C(2,2)*GRADU(2,2)
! transform stress tensor
TSN = C(1,1)*C(1,2)*TAU_IJ(1,1) + C(1,1)*C(2,2)*TAU_IJ(1,2) &
+ C(2,1)*C(1,2)*TAU_IJ(2,1) + C(2,1)*C(2,2)*TAU_IJ(2,2)
! update boundary layer equations
! update wall-normal velocity
U_NORM = U_NORM_WALL + DN*(DIVU-0.5_EB*DUSDS)
! ODE solution
IF (DNS) THEN
ETA = U_NORM + RRHO*MU/DN
AA = -(0.5_EB*DUSDS + TWTH*ETA/DN)
BB = (TWTH*U_STRM_WALL/DN + ONSI*DUSDN)*ETA - (U_NORM*0.5_EB*DUSDN + RRHO*( DPDS + TSN/(2._EB*DN) ))
!AA = -0.5_EB*(DUSDS + ETA/DN)
!BB = 0.5_EB*US_WALL/DN*ETA - (UN*0.5_EB*DUSDN + RRHO*( DPDS + TSN/(2._EB*DN) ))
U_STRM = ((AA*U_STRM + BB)*EXP(AA*DT) - BB)/AA
ELSE
U_STRM_0 = U_STRM
DO SUBIT=1,1
CALL WERNER_WENGLE_WALL_MODEL(SLIP_COEF,DUMMY,U_STRM-U_STRM_WALL,MU*RRHO,DN,0._EB)
!IF (SLIP_COEF< -1._EB .OR. SLIP_COEF>-1._EB) THEN
! PRINT *,SUBIT,'WARNING: SLIP_COEF=',SLIP_COEF
!ENDIF
ETA = RRHO*(1-SLIP_COEF)*MU/(2._EB*DN**2)
AA = -(0.5_EB*DUSDS + TWTH*U_NORM/DN + ETA)
BB = ETA*U_STRM_WALL - (U_NORM*ONTH*DUSDN + RRHO*( DPDS + TSN/(2._EB*DN) ))
U_STRM = ((AA*U_STRM_0 + BB)*EXP(AA*DT) - BB)/AA
ENDDO
ENDIF
! transform velocity back to Cartesian component I_VEL
VELTAN2D = C(I_VEL,1)*U_STRM + C(I_VEL,2)*U_NORM
END FUNCTION VELTAN2D
REAL(EB) FUNCTION VELTAN3D(U_VELO,U_SURF,NN,DN,DIVU,GRADU,GRADP,TAU_IJ,DT,RRHO,MU,I_VEL,ROUGHNESS,U_INT)
USE MATH_FUNCTIONS, ONLY: CROSS_PRODUCT, NORM2
REAL(EB), INTENT(IN) :: U_VELO(3),U_SURF(3),NN(3),DN,DIVU,GRADU(3,3),GRADP(3),TAU_IJ(3,3),DT,RRHO,MU,ROUGHNESS,U_INT
INTEGER, INTENT(IN) :: I_VEL
REAL(EB) :: C(3,3),SS(3),PP(3),SLIP_COEF,ETA,AA,BB,U_STRM_0,DUMMY,U_RELA(3), &
U_STRM,U_ORTH,U_NORM,DPDS,DUSDS,DUSDN,TSN,DUPDP,DUNDN
INTEGER :: SUBIT,I,J
! Cartesian grid coordinate system orthonormal basis vectors
REAL(EB), DIMENSION(3), PARAMETER :: E1=(/1._EB,0._EB,0._EB/),E2=(/0._EB,1._EB,0._EB/),E3=(/0._EB,0._EB,1._EB/)
! find a vector PP in the tangent plane of the surface and orthogonal to U_VELO-U_SURF
U_RELA = U_VELO-U_SURF
CALL CROSS_PRODUCT(PP,NN,U_RELA) ! PP = NN x U_RELA
IF (ABS(NORM2(PP))<=ZERO_P) THEN
! tangent vector is completely arbitrary, just perpendicular to NN
IF (ABS(NN(1))>=ZERO_P .OR. ABS(NN(2))>=ZERO_P) PP = (/NN(2),-NN(1),0._EB/)
IF (ABS(NN(1))<=ZERO_P .AND. ABS(NN(2))<=ZERO_P) PP = (/NN(3),0._EB,-NN(1)/)
ENDIF
PP = PP/NORM2(PP) ! normalize to unit vector
CALL CROSS_PRODUCT(SS,PP,NN) ! define the streamwise unit vector SS
!! check unit normal vectors
!print *,DOT_PRODUCT(SS,SS) ! should be 1
!print *,DOT_PRODUCT(SS,PP) ! should be 0
!print *,DOT_PRODUCT(SS,NN) ! should be 0
!print *,DOT_PRODUCT(PP,PP) ! should be 1
!print *,DOT_PRODUCT(PP,NN) ! should be 0
!print *,DOT_PRODUCT(NN,NN) ! should be 1
!print * ! blank line
! directional cosines (see Pope, Eq. A.11)
C(1,1) = DOT_PRODUCT(E1,SS)
C(1,2) = DOT_PRODUCT(E1,PP)
C(1,3) = DOT_PRODUCT(E1,NN)
C(2,1) = DOT_PRODUCT(E2,SS)
C(2,2) = DOT_PRODUCT(E2,PP)
C(2,3) = DOT_PRODUCT(E2,NN)
C(3,1) = DOT_PRODUCT(E3,SS)
C(3,2) = DOT_PRODUCT(E3,PP)
C(3,3) = DOT_PRODUCT(E3,NN)
! transform velocity (see Pope, Eq. A.17)
U_STRM = C(1,1)*U_RELA(1) + C(2,1)*U_RELA(2) + C(3,1)*U_RELA(3)
U_ORTH = C(1,2)*U_RELA(1) + C(2,2)*U_RELA(2) + C(3,2)*U_RELA(3)
U_NORM = C(1,3)*U_RELA(1) + C(2,3)*U_RELA(2) + C(3,3)*U_RELA(3)
!! check UP, should be zero
!print *, U_ORTH
! transform pressure gradient
DPDS = C(1,1)*GRADP(1) + C(2,1)*GRADP(2) + C(3,1)*GRADP(3)
! transform tensors (Pope A.23)
DUSDS = 0._EB
DUPDP = 0._EB
DUNDN = 0._EB
DUSDN = 0._EB
TSN = 0._EB
DO I=1,3
DO J=1,3
DUSDS = DUSDS + C(I,1)*C(J,1)*GRADU(I,J)
DUPDP = DUPDP + C(I,2)*C(J,2)*GRADU(I,J)
DUNDN = DUNDN + C(I,3)*C(J,3)*GRADU(I,J)
DUSDN = DUSDN + C(I,1)*C(J,3)*GRADU(I,J)
TSN = TSN + C(I,1)*C(J,3)*TAU_IJ(I,J)
ENDDO
ENDDO
!! check trace of transformed tensor
!print *
!print *,DIVU,DUSDS+DUPDP+DUNDN
!print *,GRADU
! update boundary layer equations
! update wall-normal velocity
U_NORM = DN*(DIVU-0.5_EB*DUSDS)
! ODE solution
IF (DNS) THEN
ETA = U_NORM + RRHO*MU/DN
AA = -(0.5_EB*DUSDS + TWTH*ETA/DN)
BB = ONSI*DUSDN*ETA - (U_NORM*0.5_EB*DUSDN + RRHO*( DPDS + TSN/(2._EB*DN) ))
IF (ABS(AA)>=ZERO_P) THEN
U_STRM = ((AA*U_STRM + BB)*EXP(AA*DT) - BB)/AA