-
Notifications
You must be signed in to change notification settings - Fork 38
/
flux.f90
1424 lines (1424 loc) · 47.8 KB
/
flux.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 flux_mod
!
!.. Use Statements ..
use module_kind_types
use eqn_idx, only : nec,nmx,nmy,nmz,nmb,nme,nee,ntk,ntl,nq
!
implicit none
!
private
!
public :: interp_cv_to_flxpts
public :: interp_dusp_to_flxpts
public :: flux_gradient_at_solpts_visc
public :: correct_interface_flux_visc
public :: discontinuous_gradient_at_solpts
public :: correct_gradients_at_solpts
public :: correct_gradients_br2
!
contains
!
!###############################################################################
!
subroutine interp_cv_to_flxpts()
!
use order_mod, only : maxSP,geom_solpts
use geovar, only : ncell,cell
use interpolation_mod, only : interp
use flowvar, only : usp,faceusp
!
use parallel_mod, only : exchange_faceusp
!
!.. Local Scalars ..
integer :: n,nc,n1,n2,np,k,m
integer :: nf,kf,lf,sf
integer :: this_geom,this_order
integer :: face_geom,face_order
!
!.. Local Arrays ..
real(wp), dimension(1:maxSP,1:nq) :: tusp
!
!.. Local Parameters ..
character(len=*), parameter :: pname = "interp_cv_to_flxpts"
!
continue
!
call debug_timer(entering_procedure,pname)
!
! Loop through each interior cell and interpolate the conserved
! variables at the solution points to the flux points of each cell
!
do nc = 1,ncell
!
n1 = cell(nc)%beg_sp ! beginning index of sol-pts for this cell
n2 = cell(nc)%end_sp ! ending index of sol-pts for this cell
np = n2-n1+1 ! number of sol-pts in this cell
!
this_geom = cell(nc)%geom ! geometry type of this cell
this_order = cell(nc)%order ! solution order of this cell
!
! Create transpose of usp array to reduce cache misses
!
do k = 1,np
do m = 1,nq
tusp(k,m) = usp(m,k+n1-1)
end do
end do
!
! Now interpolate the interior solution gradients to the flux points
!
! Loop over the faces of this cell
!
do n = 1,size(cell(nc)%face)
!
nf = cell(nc)%face(n)%idx ! index of the current face
sf = cell(nc)%face(n)%side ! side of face that cell nc is on
!
face_geom = cell(nc)%face(n)%geom ! geom of current face
face_order = cell(nc)%face(n)%order ! order of current face
!
do kf = 1,geom_solpts(face_geom,face_order)
!
! kf: index of the flx-pt local to this face
! lf: index of flx-pt within interp%mat
lf = cell(nc)%face(n)%get_fp_idx(kf)
!
faceusp(nf)%v(1:nq,kf,sf) = interp(this_geom,this_order)% &
toFace(face_order)% &
dot_mat( lf , tusp(1:np,1:nq) )
!do m = 1,nq
! faceusp(nf)%v(m,kf,sf) = interp(this_geom,this_order)% &
! toFace(face_order)% &
! dot( lf , tusp(1:np,m) )
!end do
!
end do
!
end do
!
end do
!
! Exchange the left and right states for faces on communication boundaries
!
if (ncpu > 1) call exchange_faceusp(faceusp)
!
call debug_timer(leaving_procedure,pname)
!
end subroutine interp_cv_to_flxpts
!
!###############################################################################
!
subroutine interp_dusp_to_flxpts()
!
use order_mod, only : maxSP,geom_solpts
use geovar, only : nr,ncell,cell
use interpolation_mod, only : interp
use flowvar, only : dusp,facedusp
!
use parallel_mod, only : exchange_facedusp
!
!.. Local Scalars ..
integer :: n,nc,n1,n2,np,k,l,m
integer :: nf,kf,lf,sf
integer :: this_geom,this_order
integer :: face_geom,face_order
!
!.. Local Arrays ..
real(wp), dimension(1:maxSP,1:nr,1:nq) :: tdusp
!
!.. Local Parameters ..
character(len=*), parameter :: pname = "interp_dusp_to_flxpts"
!
continue
!
call debug_timer(entering_procedure,pname)
!
! Loop through each interior cell and interpolate the conserved
! variables at the solution points to the flux points of each cell
!
do nc = 1,ncell
!
n1 = cell(nc)%beg_sp ! beginning index of sol-pts for this cell
n2 = cell(nc)%end_sp ! ending index of sol-pts for this cell
np = n2-n1+1 ! number of sol-pts in this cell
!
this_geom = cell(nc)%geom ! geometry type of this cell
this_order = cell(nc)%order ! solution order of this cell
!
! Create transpose of dusp array to reduce cache misses
!
do k = 1,np
do m = 1,nq
do l = 1,nr
tdusp(k,l,m) = dusp(l,m,k+n1-1)
end do
end do
end do
!
! Now interpolate the interior solution gradients to the flux points
!
! Loop over the faces of this cell
!
do n = 1,size(cell(nc)%face)
!
nf = cell(nc)%face(n)%idx ! index of the current face
sf = cell(nc)%face(n)%side ! side of face that cell nc is on
!
face_geom = cell(nc)%face(n)%geom ! geom of current face
face_order = cell(nc)%face(n)%order ! order of current face
!
do kf = 1,geom_solpts(face_geom,face_order)
!
! kf: index of the flx-pt local to this face
! lf: index of flx-pt within interp%mat
lf = cell(nc)%face(n)%get_fp_idx(kf)
!
facedusp(nf)%v(1:nr,1:nq,kf,sf) = interp(this_geom,this_order)% &
toFace(face_order)% &
dot_mat2( lf , tdusp(1:np,1:nr,1:nq) )
!do m = 1,nq
! do l = 1,nr
! facedusp(nf)%v(l,m,kf,sf) = interp(this_geom,this_order)% &
! toFace(face_order)% &
! dot( lf , tdusp(1:np,l,m) )
! end do
!end do
!
end do
!
end do
!
end do
!
! Exchange the left and right states for faces on communication boundaries
!
if (ncpu > 1) call exchange_facedusp(facedusp)
!
call debug_timer(leaving_procedure,pname)
!
end subroutine interp_dusp_to_flxpts
!
!###############################################################################
!###############################################################################
!###############################################################################
!
! Procedures for computing the divergence of inviscid and viscous fluxes
!
!###############################################################################
!###############################################################################
!###############################################################################
!
subroutine flux_gradient_at_solpts_visc(time_step,rk_step)
!
!... this subroutine get derivative of fluxes (no interaction)
!... for each cell, get f, g, ftd, gtd at the solution points
!
!.. Use Statements ..
use order_mod, only : maxSP,geom_solpts
use geovar, only : nr,ncell,using_quadrature
use geovar, only : cell
use interpolation_mod, only : interp
use derivatives_mod, only : deriv
use metrics_mod, only : metrics_dt,geofa
use flowvar, only : usp,dusp,residual,faceflx,faceusp
use ovar, only : governing_equations
use ovar, only : output_dir,dump_fluxes
!
!.. Formal Arguments ..
integer, intent(in) :: time_step
integer, intent(in) :: rk_step
!
!.. Local Scalars ..
integer :: k,l,m,n,nc,n1,n2,np
integer :: nf,sf,kf,lf
integer :: this_geom,this_order
integer :: face_geom,face_order
integer :: ierr,iosp
!
character(len=100) :: flux_fmt
character(len=300) :: spf_fname
!
!.. Local Arrays ..
real(wp), dimension(1:nr) :: unrm
real(wp), dimension(1:nq,1:nr) :: vflx
real(wp), dimension(1:nq,1:maxSP) :: lcl_res
real(wp), dimension(1:maxSP,1:nq,1:nr) :: rsflx
real(wp), dimension(1:maxSP,1:nq,1:nr) :: xyflx
!
!.. Local Parameters ..
character(len=*), parameter :: pname = "flux_gradient_at_solpts_visc"
!
continue
!
call debug_timer(entering_procedure,pname)
!
if (dump_fluxes) then
write (flux_fmt,4) nr,nr
write (spf_fname,1) trim(output_dir),time_step,rk_step
spf_fname = trim(adjustl(spf_fname))
open (newunit=iosp,file=trim(spf_fname),status="replace",action="write", &
position="rewind",form="formatted", &
iostat=ierr,iomsg=error_message)
call io_error(pname,trim(spf_fname),1,__LINE__,__FILE__,ierr,error_message)
end if
!
! Initialize faceflx to zero before we start
!
#ifdef SPECIAL_FOR_PGI
call faceflx%zero_face_var
#else
call faceflx%zero
#endif
!
! Loop over the cells and get the derivative of the
! fluxes at each solution point within each cell
!
do nc = 1,ncell
!
n1 = cell(nc)%beg_sp ! beginning index of sol-pts for this cell
n2 = cell(nc)%end_sp ! ending index of sol-pts for this cell
np = n2-n1+1 ! number of sol-pts in this cell
!
this_geom = cell(nc)%geom ! geometry type of this cell
this_order = cell(nc)%order ! solution order of this cell
!
! #########
! #########
! ### 1 ###
! #########
! #########
!
! Get the flux at each quadrature point in the current cell
!
do k = 1,np
!
n = k + n1-1
!
! Inviscid fluxes
!
do l = 1,nr
xyflx(k,:,l) = invs_flx_cv( usp(:,n) , l ) ! x-directional flux
end do
!
if (governing_equations == NavierStokes_Eqns) then
!
! Viscous fluxes
!
vflx(1:nq,1:nr) = visc_flx_cv( usp(1:nq,n) , dusp(1:nr,1:nq,n) )
!
! Subtract the viscous flux from the inviscid flux
!
xyflx(k,1:nq,1:nr) = xyflx(k,1:nq,1:nr) - vflx(1:nq,1:nr)
!
end if
!
! Transform the fluxes in the x-y physical
! domain to the r-s computational domain
!
do m = 1,nq
do l = 1,nr
rsflx(k,m,l) = dot(xyflx(k,m,1:nr),metrics_dt(nc)%met(k,l,1:nr)) * &
metrics_dt(nc)%jac(k)
end do
end do
!
end do
!
if (dump_fluxes) then
write (iosp,2) nc
do k = 1,np
write (iosp,3) k,k+n1-1
do m = 1,nq
write (iosp,flux_fmt) (xyflx(k,m,l),l=1,nr), (rsflx(k,m,l),l=1,nr)
end do
end do
end if
!
! #########
! #########
! ### 2 ###
! #########
! #########
!
! Compute the divergence of the discontinuous fluxes
!
! Zero the residual array for the current cell
!
lcl_res = zero
!
! Loop through each coordinate direction and add the divergence of
! the flux to the local residual array for the current cell
!
do m = 1,nq
lcl_res(m,1:np) = deriv(this_geom,this_order)%div( rsflx(1:np,m,1:nr) )
end do
!
! Copy the local residual for the current cell into the residual array
!
residual(1:nq,n1:n2) = lcl_res(1:nq,1:np)
!
! #########
! #########
! ### 3 ###
! #########
! #########
!
! Interpolate the fluxes to the interface flux points and dot this
! interpolated vector with the face normal at each flux point
! NOTE: We dont need to exchange this information for communication and
! periodic boundary faces since this flux is only used within the
! current cell.
!
! Loop over the faces of this cell
!
do n = 1,size(cell(nc)%face)
!
nf = cell(nc)%face(n)%idx ! index of the current face
sf = cell(nc)%face(n)%side ! side of face that cell nc is on
!
face_geom = cell(nc)%face(n)%geom ! geom of current face
face_order = cell(nc)%face(n)%order ! order of current face
!
do kf = 1,geom_solpts(face_geom,face_order)
!
! kf: index of the flx-pt local to this face
! lf: index of flx-pt within interp%mat
lf = cell(nc)%face(n)%get_fp_idx(kf)
!
! Get the unit normal for this face flux/solution point
!
unrm(:) = unit_vector( geofa(nf)%v(:,kf) )
!
vflx(1:nq,1:nr) = interp(this_geom,this_order)% &
toFace(face_order)% &
dot_mat2( lf , xyflx(1:np,1:nq,1:nr) )
!do l = 1,nr
! do m = 1,nq
! vflx(m,l) = interp(this_geom,this_order)% &
! toFace(face_order)% &
! dot( lf , xyflx(1:np,m,l) )
! end do
!end do
!
do m = 1,nq
faceflx(nf)%v(m,kf,sf) = dot( vflx(m,1:nr) , unrm(1:nr) )
end do
!
end do
!
end do
!
end do
!
if (dump_fluxes) then
close (iosp,iostat=ierr,iomsg=error_message)
call io_error(pname,trim(spf_fname),2,__LINE__,__FILE__,ierr,error_message)
end if
!
call debug_timer(leaving_procedure,pname)
!
! Format Statements
!
1 format (a,"flux_files/sp_",i0,".",i0,"_fluxes.dat")
2 format (100("#"),/,15x,"CELL ",i0,/,100("#"))
3 format (100("-"),/,5x,"SP ",i0," (",i0,")",/,100("-"))
4 format ("(",i0,"es15.6,5x,",i0,"es15.6)")
!
end subroutine flux_gradient_at_solpts_visc
!
!###############################################################################
!
subroutine correct_interface_flux_visc(time_step,rk_step)
!
use generic_types_mod, only : matrix
use order_mod, only : maxFP,geom_solpts
use geovar, only : nr,nfbnd,nface
use geovar, only : bface,face,cell
use ovar, only : invs_flux_method,governing_equations
use ovar, only : output_dir,dump_fluxes
use metrics_mod, only : geofa
use correction_mod, only : correct
use flowvar, only : usp,residual,faceflx,faceusp,facedusp
!
!.. Formal Arguments ..
integer, intent(in) :: time_step
integer, intent(in) :: rk_step
!
!.. Local Scalars ..
integer :: i,k,l,m,n,nf,r1,r2,l1,l2
integer :: kl,kr,iu,nfp
integer :: ierr,iofp
integer :: face_geom,face_order
integer :: host_cell,host_geom,host_order
integer :: left_cell,left_geom,left_order
integer :: right_cell,right_geom,right_order
real(wp) :: ds,dsl,dsr,cmat_val
!
character(len=300) :: fpf_fname
!
!.. Local Arrays ..
real(wp), dimension(1:nr) :: unrm
real(wp), dimension(1:nq) :: iflx,vflx
real(wp), dimension(1:nr,1:maxFP) :: rnrm
real(wp), dimension(1:nq,1:maxFP) :: fjl,fjr
real(wp), dimension(1:nq,1:maxFP) :: lflx,rflx
real(wp), dimension(1:nq,1:maxFP) :: luqp,ruqp,fncq
real(wp), dimension(1:nr,1:nq,1:maxFP) :: lduqp,rduqp
!
!.. Local Parameters ..
character(len=*), parameter :: pname = "correct_interface_flux_visc"
!
continue
!
call debug_timer(entering_procedure,pname)
!
if (dump_fluxes) then
write (fpf_fname,2) trim(output_dir),time_step,rk_step
fpf_fname = trim(adjustl(fpf_fname))
open (newunit=iofp,file=trim(fpf_fname),status="replace",action="write", &
position="rewind",form="formatted", &
iostat=ierr,iomsg=error_message)
call io_error(pname,trim(fpf_fname),1,__LINE__,__FILE__,ierr,error_message)
end if
!
vflx = zero
!
! ##############################################################
! ##############################################################
! ##############################################################
! ########## BOUNDARY FACES ##########
! ##############################################################
! ##############################################################
! ##############################################################
!
! Loop through the boundary faces to apply the correction
! due to the boundary fluxes.
!
boundary_fluxes: do nf = 1,nfbnd
!
face_geom = face(nf)%geom ! geometry of face
face_order = face(nf)%order ! order of face
!
nfp = geom_solpts(face_geom,face_order) ! number of flux points on this face
!
host_cell = face(nf)%left%cell ! host cell on boundary face
host_geom = cell(host_cell)%geom ! geometry of host cell
host_order = cell(host_cell)%order ! order of host cell
!
l1 = cell(host_cell)%beg_sp ! beginning index for sol-pts in left cell
l2 = cell(host_cell)%end_sp ! ending index for sol-pts in left cell
!
! value of interior flux polynomial at FP
lflx(1:nq,1:nfp) = faceflx(nf)%v(1:nq,1:nfp,1)
!
! FP solution from left cell on face
luqp(1:nq,1:nfp) = faceusp(nf)%v(1:nq,1:nfp,1)
! FP solution from right cell on face
ruqp(1:nq,1:nfp) = faceusp(nf)%v(1:nq,1:nfp,2)
!
if (governing_equations == NavierStokes_Eqns) then
! Solution gradient at FP from left cell on face
lduqp(1:nr,1:nq,1:nfp) = facedusp(nf)%v(1:nr,1:nq,1:nfp,1)
! Solution gradient at FP from right cell on face
rduqp(1:nr,1:nq,1:nfp) = facedusp(nf)%v(1:nr,1:nq,1:nfp,2)
end if
!
rnrm(1:nr,1:nfp) = geofa(nf)%v(1:nr,1:nfp) ! Area face normals at FP
!
! Get the common flux at each quadrature flux point
!
do k = 1,nfp
!
! Get the unit normal at this flux point
!
unrm(:) = unit_vector( rnrm(:,k) )
!
! Get the inviscid upwind flux at the current flux point
!
if (any(bface(1,nf) == bc_walls)) then
!
! This is a wall boundary so ensure that the inviscid
! boundary flux is only due to the pressure term.
!
iflx(:) = invs_wall_flx_cv( unrm(:) , luqp(:,k) )
!
! Compute the viscous wall boundary flux
!
if (governing_equations == NavierStokes_Eqns) then
vflx(:) = visc_wall_flx_cv( unrm(:) , ruqp(:,k) , rduqp(:,:,k) )
end if
!
else
!
! For all other boundaries, compute the boundary
! fluxes in the same way as the interior faces.
!
if (invs_flux_method == invs_flux_ROE) then
iflx(:) = upwind_flux_cv ( unrm(:) , luqp(:,k) , ruqp(:,k) )
else if (invs_flux_method == invs_flux_HLLC) then
iflx(:) = hllc_upwind_flux_cv ( unrm(:) , luqp(:,k) , ruqp(:,k) )
else if (invs_flux_method == invs_flux_LDFSS) then
iflx(:) = ldfss_upwind_flux_cv ( unrm(:) , luqp(:,k) , ruqp(:,k) )
else if (invs_flux_method == invs_flux_Rotated_RHLL) then
if (nr == 2) then
iflx(:) = rotated_rhll_upwind_flux_cv_2d ( unrm(:) , luqp(:,k) , &
ruqp(:,k) )
else
iflx(:) = rotated_rhll_upwind_flux_cv_3d ( unrm(:) , luqp(:,k) , &
ruqp(:,k) )
end if
end if
!
! Get the common viscous flux at the current flux point
!
if (governing_equations == NavierStokes_Eqns) then
vflx(:) = common_visc_flx_cv( unrm(:), luqp(:,k), lduqp(:,:,k), &
ruqp(:,k), rduqp(:,:,k) )
end if
!
end if
!
! Get the total common boundary flux at this flux point
!
fncq(:,k) = iflx(:) - vflx(:)
!
end do
!
if (dump_fluxes) then
write (iofp,3) nf,host_cell
do k = 1,nfp
write (iofp,4) k,face(nf)%left%get_fp_idx(face_geom,face_order,k)
do m = 1,nq
write (iofp,7) fncq(m,k), lflx(m,k)
end do
end do
end if
!
! Compute the flux jump at the interface
!
fjl(1:nq,1:nfp) = fncq(1:nq,1:nfp) - lflx(1:nq,1:nfp)
!
! Apply the flux correction due to the flux jump
! NOTE: We only need to apply the correction to the left cell because
! this is a boundary face and the right cell is a ghost cell.
!
do k = 1,nfp
!
! Get the index for the correction function
! corresponding to the current interface flux point
!
kl = face(nf)%left%get_fp_idx(face_geom,face_order,k)
!
! Get the contribution to the face length/area from the current flux
! point for the left cell
!
ds = norm2( rnrm(:,k) )
dsl = ds * cl(host_geom)
!
! Add the flux correction to the residuals of the sol-pts in cell
! host_cell
!
do i = correct(host_geom,host_order)%row_ptr(kl)+1, &
correct(host_geom,host_order)%row_ptr(kl+1)
n = l1 - 1 + correct(host_geom,host_order)%row_idx(i)
cmat_val = dsl * correct(host_geom,host_order)%row_val(i)
do m = 1,nq
residual(m,n) = residual(m,n) + fjl(m,k) * cmat_val
end do
end do
!
end do
!
end do boundary_fluxes
!
! ##############################################################
! ##############################################################
! ##############################################################
! ########## INTERIOR FACES ##########
! ##############################################################
! ##############################################################
! ##############################################################
!
! Loop over the interior faces and apply the correction from
! the fluxes at each face point to all the cell solution points
!
interior_fluxes: do nf = nfbnd+1,nface
!
face_geom = face(nf)%geom ! geometry of face
face_order = face(nf)%order ! order of face
!
nfp = geom_solpts(face_geom,face_order) ! number of flux points on this face
!
left_cell = face(nf)%left%cell ! left cell on face
right_cell = face(nf)%right%cell ! right cell on face
!
left_geom = cell(left_cell)%geom
left_order = cell(left_cell)%order
!
right_geom = cell(right_cell)%geom
right_order = cell(right_cell)%order
!
l1 = cell(left_cell)%beg_sp ! beginning index for sol-pts in left cell
l2 = cell(left_cell)%end_sp ! ending index for sol-pts in left cell
!
r1 = cell(right_cell)%beg_sp ! beginning index for sol-pts in right cell
r2 = cell(right_cell)%end_sp ! ending index for sol-pts in right cell
!
lflx(1:nq,1:nfp) = faceflx(nf)%v(1:nq,1:nfp,1)
rflx(1:nq,1:nfp) = faceflx(nf)%v(1:nq,1:nfp,2)
!
luqp(1:nq,1:nfp) = faceusp(nf)%v(1:nq,1:nfp,1)
ruqp(1:nq,1:nfp) = faceusp(nf)%v(1:nq,1:nfp,2)
!
if (governing_equations == NavierStokes_Eqns) then
lduqp(1:nr,1:nq,1:nfp) = facedusp(nf)%v(1:nr,1:nq,1:nfp,1)
rduqp(1:nr,1:nq,1:nfp) = facedusp(nf)%v(1:nr,1:nq,1:nfp,2)
end if
!
rnrm(1:nr,1:nfp) = geofa(nf)%v(1:nr,1:nfp)
!
! Get the common flux at each quadrature flux point
!
do k = 1,nfp
!
! Get the unit normal at this flux point
!
unrm(:) = unit_vector( rnrm(:,k) )
!
! Get the inviscid upwind flux at the current flux point
!
if (invs_flux_method == invs_flux_ROE) then
iflx(:) = upwind_flux_cv ( unrm(:) , luqp(:,k) , ruqp(:,k) )
else if (invs_flux_method == invs_flux_HLLC) then
iflx(:) = hllc_upwind_flux_cv ( unrm(:) , luqp(:,k) , ruqp(:,k) )
else if (invs_flux_method == invs_flux_LDFSS) then
iflx(:) = ldfss_upwind_flux_cv ( unrm(:) , luqp(:,k) , ruqp(:,k) )
else if (invs_flux_method == invs_flux_Rotated_RHLL) then
if (nr == 2) then
iflx(:) = rotated_rhll_upwind_flux_cv_2d ( unrm(:) , luqp(:,k) , &
ruqp(:,k) )
else
iflx(:) = rotated_rhll_upwind_flux_cv_3d ( unrm(:) , luqp(:,k) , &
ruqp(:,k) )
end if
end if
!
! Get the common viscous flux at the current flux point
!
if (governing_equations == NavierStokes_Eqns) then
vflx(:) = common_visc_flx_cv( unrm(:), luqp(:,k), lduqp(:,:,k), &
ruqp(:,k), rduqp(:,:,k) )
end if
!
! Get the total common boundary flux at this flux point
!
fncq(:,k) = iflx(:) - vflx(:)
!
end do
!
if (dump_fluxes) then
write (iofp,5) nf,left_cell,right_cell
do k = 1,nfp
write (iofp,6) k,face(nf)%left%get_fp_idx(face_geom,face_order,k), &
face(nf)%right%get_fp_idx(face_geom,face_order,k)
do m = 1,nq
write (iofp,7) fncq(m,k), lflx(m,k), rflx(m,k)
end do
end do
end if
!
! Compute the flux jump at the interface
!
fjl(1:nq,1:nfp) = fncq(1:nq,1:nfp) - lflx(1:nq,1:nfp)
fjr(1:nq,1:nfp) = fncq(1:nq,1:nfp) - rflx(1:nq,1:nfp)
!
! Apply the flux correction due to the flux jump
!
do k = 1,nfp
!
! Get the index for the correction function
! corresponding to the current interface flux point
!
kl = face(nf)%left%get_fp_idx(face_geom,face_order,k)
kr = face(nf)%right%get_fp_idx(face_geom,face_order,k)
!
! Get the face length for the left and right cell
!
ds = norm2( rnrm(:,k) )
dsl = ds * cl(left_geom)
dsr = ds * cl(right_geom)
!
! Add the flux correction to the residuals of the sol-pts in cell
! left_cell
!
do i = correct(left_geom,left_order)%row_ptr(kl)+1, &
correct(left_geom,left_order)%row_ptr(kl+1)
n = l1 - 1 + correct(left_geom,left_order)%row_idx(i)
cmat_val = dsl * correct(left_geom,left_order)%row_val(i)
do m = 1,nq
residual(m,n) = residual(m,n) + fjl(m,k) * cmat_val
end do
end do
!
! Add the flux correction to the residuals of the sol-pts in cell
! right_cell
! NOTE: only do this for interior cells, i.e. skip this for ghost cells
!
do i = correct(right_geom,right_order)%row_ptr(kr)+1, &
correct(right_geom,right_order)%row_ptr(kr+1)
n = r1 - 1 + correct(right_geom,right_order)%row_idx(i)
cmat_val = dsr * correct(right_geom,right_order)%row_val(i)
do m = 1,nq
residual(m,n) = residual(m,n) - fjr(m,k) * cmat_val
end do
end do
!
end do
!
end do interior_fluxes
!
if (dump_fluxes) then
close (iofp,iostat=ierr,iomsg=error_message)
call io_error(pname,trim(fpf_fname),2,__LINE__,__FILE__,ierr,error_message)
end if
!
call debug_timer(leaving_procedure,pname)
!
! Format Statements
!
1 format (a,i0)
2 format (a,"flux_files/fp_",i0,".",i0,"_fluxes.dat")
3 format (100("#"),/, &
15x,"BOUNDARY FACE ",i0," : HC=",i0,/, &
100("#"))
4 format (100("-"),/, &
5x,"FP ",i0," : HP=",i0,/, &
100("-"))
5 format (100("#"),/, &
15x,"INTERIOR FACE ",i0," : LC=",i0,", RC=",i0,/, &
100("#"))
6 format (100("-"),/, &
5x,"FP ",i0," : LP=",i0,", RP=",i0,/, &
100("-"))
7 format (es16.6,5x,es15.6,:,es15.6)
!
end subroutine correct_interface_flux_visc
!
!###############################################################################
!###############################################################################
!###############################################################################
!
! Procedures for getting the conservative variable gradients
!
!###############################################################################
!###############################################################################
!###############################################################################
!
subroutine discontinuous_gradient_at_solpts()
!
! This routine computes the discontinuous gradient at each solution point
!
!.. Use Statements ..
use order_mod, only : maxpts
use geovar, only : nr,ncell,cell
use metrics_mod, only : metrics_dt
use derivatives_mod, only : deriv
use flowvar, only : usp,dusp
!
!.. Local Scalars ..
integer :: k,l,m,n,nc,n1,n2,np
integer :: this_geom,this_order
!
!.. Local Arrays ..
!real(wp), dimension(1:nr) :: dudr
real(wp), dimension(1:nr,1:nq) :: dudr
real(wp), dimension(1:maxpts,1:nq) :: tusp
!
!.. Local Parameters ..
character(len=*), parameter :: pname = "discontinuous_gradient_at_solpts"
!
continue
!
call debug_timer(entering_procedure,pname)
!
! Loop over the cells and get the derivative of
! the primitive variables at each solution point
!
do nc = 1,ncell
!
n1 = cell(nc)%beg_sp ! beginning index for solpts in cell nc
n2 = cell(nc)%end_sp ! ending index for solpts in cell nc
np = n2-n1+1 ! number of solution points in this cell
!
this_geom = cell(nc)%geom
this_order = cell(nc)%order
!
! Store the transpose of usp for this cell to reduce cache misses
!
do k = 1,np
do m = 1,nq
tusp(k,m) = usp(m,k+n1-1)
end do
end do
!
! Compute the gradient of the conservative variables at each solution point
!
do k = 1,np
!
n = k+n1-1 ! total index of solution point on current partition
!
dudr(1:nr,1:nq) = deriv(this_geom,this_order)%grad( k , tusp(1:np,1:nq) )
!
do m = 1,nq
do l = 1,nr
dusp(l,m,n) = dot( dudr(1:nr,m) , metrics_dt(nc)%met(k,1:nr,l) )
end do
end do
!
!do m = 1,nq
! !
! ! Compute the derivatives in r-s space
! !
! do l = 1,nr
! dudr(l) = deriv(this_geom,this_order)% &
! d(l)% &
! dot( k , tusp(1:np,m) )
! end do
! !
! ! Transform the derivatives back to x-y space
! !
! do l = 1,nr
! dusp(l,m,n) = dot( dudr(1:nr) , metrics_dt(nc)%met(k,1:nr,l) )
! end do
! !
!end do
!
end do
!
end do
!
call debug_timer(leaving_procedure,pname)
!
end subroutine discontinuous_gradient_at_solpts
!
!###############################################################################
!
subroutine correct_gradients_at_solpts()
!
!.. Use Statements ..
use generic_types_mod, only : matrix
use order_mod, only : maxpts,geom_solpts
use geovar, only : nfbnd,nface,nr,n_solpts,ncell
use geovar, only : face,cell,bface
use metrics_mod, only : geofa,metrics_dt
use correction_mod, only : correct
use flowvar, only : dusp,faceusp
use ovar, only : walls_are_exact
use ovar, only : LDG_beta
!
!.. Local Scalars ..
integer :: i,k,l,m,n,nk,nf,r1,r2,l1,l2,n1,n2,nc
integer :: kl,kr,nfp
integer :: face_geom,face_order
integer :: host_cell,host_geom,host_order
integer :: left_cell,left_geom,left_order
integer :: right_cell,right_geom,right_order
real(wp) :: cbias,hfpb,hfmb,cmat_val
real(wp) :: correct_n,correct_m,ds,dsl,dsr
!
!.. Local Arrays ..
real(wp), dimension(1:nr) :: unrm
real(wp), dimension(1:nq) :: ucom
real(wp), dimension(1:nq) :: ljfp,rjfp
real(wp), dimension(1:nq,1:maxpts) :: lufp,rufp
!
!.. Local Parameters ..
character(len=*), parameter :: pname = "correct_gradients_at_solpts"
!
continue
!
call debug_timer(entering_procedure,pname)
!
! #######################
! #######################
! ### BOUNDARY FACES ###
! #######################
! #######################
!
! Loop through the boundary faces and correct the gradient of the
! continuous solution polynomial in the cells adjacent to the boundary.
!
boundary_faces: do nf = 1,nfbnd
!
face_geom = face(nf)%geom ! geometry of face
face_order = face(nf)%order ! order of face
!
nfp = geom_solpts(face_geom,face_order) ! number of flux points on this face
!
host_cell = face(nf)%left%cell ! left cell on face
!
host_geom = cell(host_cell)%geom
host_order = cell(host_cell)%order
!
l1 = cell(host_cell)%beg_sp ! beginning index for sol-pts in left cell
l2 = cell(host_cell)%end_sp ! ending index for sol-pts in left cell
!
! For all the flux points on the current face, extract the conservative
! conservative variables for the left and right cells from faceusp
!
lufp(1:nq,1:nfp) = faceusp(nf)%v(1:nq,1:nfp,1)
rufp(1:nq,1:nfp) = faceusp(nf)%v(1:nq,1:nfp,2)
!
! Correct the gradients of the conservative variables using the
! solution jumps at the interface and the correction function
!
! Use the LDG flux formulation for computing the common interface solutions
! u_com = half*(u_left + u_right) - &
! (u_left - u_right)*dot(LDG_beta,unit_normal)
!
do k = 1,nfp
!
! Get the index for the correction function
! corresponding to the current interface flux point
!
kl = face(nf)%left%get_fp_idx(face_geom,face_order,k)
!
! Get the unit normal at this flux point
!
unrm(1:nr) = unit_vector( geofa(nf)%v(1:nr,k) )
!
! Get the face length for the current face
!
dsl = cl(host_geom) * norm2( geofa(nf)%v(1:nr,k) )
!
! Compute the LDG directional parameters
! NOTE: cbias is used for all boundary faces that allow for flow
! across the interface, e.g, inflow/outflow, cpu/partition
! boundary, or periodic boundary. These boundaries are
! treated the same as an interior face except that we dont
! need to worry about correcting the right cell.
!