MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_physics_hllc.t
Go to the documentation of this file.
2 
3  implicit none
4  public
5 
6  procedure(sub_diffuse_hllcd), pointer :: phys_diffuse_hllcd => null()
7  procedure(sub_get_lcd), pointer :: phys_get_lcd => null()
8  procedure(sub_get_wcd), pointer :: phys_get_wcd => null()
9  procedure(sub_hllc_init_species), pointer :: phys_hllc_init_species => null()
10 
11  abstract interface
12  subroutine sub_diffuse_hllcd(ixI^L,ixO^L,idims,wLC,wRC,fLC,fRC,patchf)
14  integer, intent(in) :: ixI^L,ixO^L,idims
15  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
16  double precision, dimension(ixI^S,1:nwflux),intent(in) :: fLC, fRC
17  integer, dimension(ixI^S), intent(inout) :: patchf
18  end subroutine sub_diffuse_hllcd
19 
20  subroutine sub_get_lcd(wLC,wRC,fLC,fRC,cmin,cmax,idims,ixI^L,ixO^L, &
21  whll,Fhll,lambdaCD,patchf)
23  integer, intent(in) :: ixI^L,ixO^L,idims
24  double precision, dimension(ixI^S,1:nw), intent(in) :: wLC,wRC
25  double precision, dimension(ixI^S,1:nwflux), intent(in) :: fLC,fRC
26  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
27  integer, dimension(ixI^S), intent(inout) :: patchf
28  double precision, dimension(ixI^S,1:nwflux), intent(out) :: Fhll,whll
29  double precision, dimension(ixI^S), intent(out) :: lambdaCD
30  end subroutine sub_get_lcd
31 
32  subroutine sub_get_wcd(wLC,wRC,whll,fRC,fLC,Fhll,patchf,lambdaCD,cmin,cmax,&
33  ixI^L,ixO^L,idims,f)
35  integer, intent(in) :: ixI^L,ixO^L,idims
36  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
37  double precision, dimension(ixI^S,1:nwflux), intent(in) :: whll, Fhll
38  double precision, dimension(ixI^S), intent(in) :: lambdaCD
39  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
40  double precision, dimension(ixI^S,1:nwflux), intent(in) :: fRC,fLC
41  integer, dimension(ixI^S), intent(in) :: patchf
42  double precision, dimension(ixI^S,1:nwflux),intent(out) :: f
43  end subroutine sub_get_wcd
44 
45  subroutine sub_hllc_init_species(ii, rho_, mom, e_)
47  integer, intent(in) :: ii
48  integer, intent(out) :: rho_, mom(1:ndir), e_
49  end subroutine sub_hllc_init_species
50  end interface
51 
52 contains
53 
54  subroutine phys_hllc_check
55  if (.not. associated(phys_diffuse_hllcd)) &
57 
58  if (.not. associated(phys_get_lcd)) &
60 
61  if (.not. associated(phys_get_wcd)) &
63  end subroutine phys_hllc_check
64 
65  ! When method is hllcd or hllcd1 then: this subroutine is to impose enforce
66  ! regions where we AVOID HLLC and use TVDLF instead: this is achieved by setting
67  ! patchf to 4 in certain regions. An additional input parameter is nxdiffusehllc
68  ! which sets the size of the fallback region. This nul version enforces TVDLF
69  ! everywhere!
70  subroutine dummy_diffuse_hllcd(ixI^L,ixO^L,idims,wLC,wRC,fLC,fRC,patchf)
72  integer, intent(in) :: ixI^L,ixO^L,idims
73  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
74  double precision, dimension(ixI^S,1:nwflux),intent(in) :: fLC, fRC
75  integer, dimension(ixI^S), intent(inout) :: patchf
76 
77  patchf(ixo^s) = 4 ! enforce TVDLF everywhere
78  end subroutine dummy_diffuse_hllcd
79 
80  ! Calculate lambda at CD and set the patchf to know the orientation
81  ! of the riemann fan and decide on the flux choice
82  ! We also compute here the HLL flux and w value, for fallback strategy
83  ! In this nul version, we simply compute nothing and ensure TVDLF fallback
84  subroutine dummy_get_lcd(wLC,wRC,fLC,fRC,cmin,cmax,idims,ixI^L,ixO^L, &
85  whll,Fhll,lambdaCD,patchf)
87  integer, intent(in) :: ixI^L,ixO^L,idims
88  double precision, dimension(ixI^S,1:nw), intent(in) :: wLC,wRC
89  double precision, dimension(ixI^S,1:nwflux), intent(in) :: fLC,fRC
90  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
91  integer, dimension(ixI^S), intent(inout) :: patchf
92  double precision, dimension(ixI^S,1:nwflux), intent(out) :: Fhll,whll
93  double precision, dimension(ixI^S), intent(out) :: lambdaCD
94 
95  ! Next must normally be computed
96  fhll(ixo^s,1:nwflux) = zero
97  whll(ixo^s,1:nwflux) = zero
98  lambdacd(ixo^s) = zero
99 
100  ! This actually ensures fallback to TVDLF
101  patchf(ixo^s)=4
102  end subroutine dummy_get_lcd
103 
104  ! compute the intermediate state U*. Only needed where patchf=-1/1. This nul
105  ! version simply nullifies all values
106  subroutine dummy_get_wcd(wLC,wRC,whll,fRC,fLC,Fhll,patchf,lambdaCD,cmin,cmax,&
107  ixI^L,ixO^L,idims,f)
109  integer, intent(in) :: ixI^L,ixO^L,idims
110  double precision, dimension(ixI^S,1:nw), intent(in) :: wRC,wLC
111  double precision, dimension(ixI^S,1:nwflux), intent(in) :: whll, Fhll
112  double precision, dimension(ixI^S), intent(in) :: lambdaCD
113  double precision, dimension(ixI^S), intent(in) :: cmax,cmin
114  double precision, dimension(ixI^S,1:nwflux), intent(in) :: fRC,fLC
115  integer, dimension(ixI^S), intent(in) :: patchf
116  double precision, dimension(ixI^S,1:nwflux),intent(out) :: f
117 
118  ! Next must normally be computed
119  f(ixo^s,1:nwflux) = zero
120  end subroutine dummy_get_wcd
121 
122 end module mod_physics_hllc
This module contains definitions of global parameters and variables and some generic functions/subrou...
procedure(sub_get_wcd), pointer phys_get_wcd
procedure(sub_hllc_init_species), pointer phys_hllc_init_species
subroutine dummy_get_wcd(wLC, wRC, whll, fRC, fLC, Fhll, patchf, lambdaCD, cmin, cmax, ixIL, ixOL, idims, f)
procedure(sub_get_lcd), pointer phys_get_lcd
procedure(sub_diffuse_hllcd), pointer phys_diffuse_hllcd
subroutine dummy_get_lcd(wLC, wRC, fLC, fRC, cmin, cmax, idims, ixIL, ixOL, whll, Fhll, lambdaCD, patchf)
subroutine phys_hllc_check
subroutine dummy_diffuse_hllcd(ixIL, ixOL, idims, wLC, wRC, fLC, fRC, patchf)