MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_mhd_hllc.t
Go to the documentation of this file.
4
5 implicit none
6 private
7
8 public :: mhd_hllc_init
9
10contains
11
12 subroutine mhd_hllc_init()
14
15 phys_diffuse_hllcd => mhd_diffuse_hllcd
16 phys_get_lcd => mhd_get_lcd
17 phys_get_wcd => mhd_get_wcd
18
19 end subroutine mhd_hllc_init
20
21 subroutine mhd_diffuse_hllcd(ixI^L,ixO^L,idim,wLC,wRC,fLC,fRC,patchf)
22 ! when method is hllcd or hllcd1 then:
23 ! this subroutine is to enforce regions where we AVOID HLLC
24 ! and use TVDLF instead: this is achieved by setting patchf to 4 in
25 ! certain regions. An additional input parameter is nxdiffusehllc
26 ! which sets the size of the fallback region.
28
29 integer, intent(in) :: ixi^l,ixo^l,idim
30 double precision, dimension(ixI^S,1:nw), intent(in) :: wrc,wlc
31 double precision, dimension(ixI^S,1:nwflux),intent(in) :: flc, frc
32 integer, dimension(ixI^S), intent(inout) :: patchf
33
34 integer :: ixoo^d,txoo^l
35
36
37 ! In a user-controlled region around any point with flux sign change between
38 ! left and right, ensure fallback to TVDLF
39 {do ixoo^d= ixo^lim^d\}
40 {
41 txoomin^d= max(ixoo^d - nxdiffusehllc*kr(idim,^d), ixomin^d);
42 txoomax^d= min(ixoo^d + nxdiffusehllc*kr(idim,^d), ixomax^d);
43 \}
44 if(abs(patchf(ixoo^d)) == 1 .or. abs(patchf(ixoo^d)) == 4)Then
45 if(any(frc(ixoo^d,1:nwflux)*flc(ixoo^d,1:nwflux)<-smalldouble))Then
46 where(abs(patchf(txoo^s))==1)
47 patchf(txoo^s) = 4
48 endwhere
49 endif
50 endif
51 {enddo^d&\}
52
53 end subroutine mhd_diffuse_hllcd
54
55 subroutine mhd_get_lcd(wLC,wRC,fLC,fRC,cmin,cmax,idim,ixI^L,ixO^L, &
56 whll,Fhll,lambdaCD,patchf)
57 ! Calculate lambda at CD and set the patchf to know the orientation
58 ! of the riemann fan and decide on the flux choice
59 ! We also compute here the HLL flux and w value, for fallback strategy
61
62 integer, intent(in) :: ixi^l,ixo^l,idim
63 double precision, dimension(ixI^S,1:nw), intent(in) :: wlc,wrc
64 double precision, dimension(ixI^S,1:nwflux), intent(in) :: flc,frc
65 double precision, dimension(ixI^S), intent(in) :: cmax,cmin
66 integer , dimension(ixI^S), intent(inout) :: patchf
67 double precision, dimension(ixI^S,1:nwflux), intent(out) :: fhll,whll
68 double precision, dimension(ixI^S), intent(out) :: lambdacd
69
70 logical , dimension(ixO^S) :: cond_patchf
71 double precision :: epsilon
72 integer :: iw,ix^d
73
74 ! on entry, patch is preset to contain values from -2,1,2,4
75 ! -2: take left flux, no computation here
76 ! +2: take right flux, no computation here
77 ! +4: take TVDLF flux, no computation here
78 ! 1: compute the characteristic speed for the CD
79
80 cond_patchf(ixo^s)=(abs(patchf(ixo^s))==1)
81 lambdacd=0.d0
82
83 do iw=1,nwflux
84 where(cond_patchf(ixo^s))
85 !============= compute HLL flux ==============!
86 fhll(ixo^s,iw)= (cmax(ixo^s)*flc(ixo^s,iw)-cmin(ixo^s)*frc(ixo^s,iw) &
87 + cmin(ixo^s)*cmax(ixo^s)*(wrc(ixo^s,iw)-wlc(ixo^s,iw)))&
88 /(cmax(ixo^s)-cmin(ixo^s))
89 !======== compute intermediate HLL state =======!
90 whll(ixo^s,iw) = (cmax(ixo^s)*wrc(ixo^s,iw)-cmin(ixo^s)*wlc(ixo^s,iw)&
91 +flc(ixo^s,iw)-frc(ixo^s,iw))/(cmax(ixo^s)-cmin(ixo^s))
92 end where
93 end do
94
95 ! deduce the characteristic speed at the CD
96 where(cond_patchf(ixo^s))
97 lambdacd(ixo^s)=whll(ixo^s,mom(idim))/whll(ixo^s,rho_)
98 end where
99
100 {do ix^db=ixomin^db,ixomax^db\}
101 if(cond_patchf(ix^d)) then
102 ! double check whether obtained speed is in between min and max speeds given
103 ! and identify in which part of the Riemann fan the time-axis is
104 if(cmin(ix^d)<zero.and.lambdacd(ix^d)>zero&
105 .and.lambdacd(ix^d)<cmax(ix^d)) then
106 patchf(ix^d) = -1
107 else if(cmax(ix^d)>zero.and.lambdacd(ix^d)<zero&
108 .and.lambdacd(ix^d)>cmin(ix^d)) then
109 patchf(ix^d) = 1
110 else if(lambdacd(ix^d)>=cmax(ix^d).or.lambdacd(ix^d) <= cmin(ix^d)) then
111 lambdacd(ix^d) = zero
112 ! we will fall back to HLL flux case in this degeneracy
113 patchf(ix^d) = 3
114 end if
115 end if
116 {end do\}
117
118 where(patchf(ixo^s)== 3)
119 cond_patchf(ixo^s)=.false.
120 end where
121
122 ! handle the specific case where the time axis is exactly on the CD
123 if(any(lambdacd(ixo^s)==zero.and.cond_patchf(ixo^s)))then
124 ! determine which sector (forward or backward) of the Riemann fan is smallest
125 ! and select left or right flux accordingly
126 {do ix^db=ixomin^db,ixomax^db\}
127 if(lambdacd(ix^d)==zero.and.cond_patchf(ix^d)) then
128 if(-cmin(ix^d)>=cmax(ix^d)) then
129 patchf(ix^d) = 1
130 else
131 patchf(ix^d) = -1
132 end if
133 end if
134 {end do\}
135 end if
136
137 ! eigenvalue lambda for contact is near zero: decrease noise by this trick
138 if(flathllc)then
139 epsilon=1.0d-6
140 where(cond_patchf(ixo^s).and. &
141 dabs(lambdacd(ixo^s))/max(cmax(ixo^s),epsilon)< epsilon .and. &
142 dabs(lambdacd(ixo^s))/max(dabs(cmin(ixo^s)),epsilon)< epsilon)
143 lambdacd(ixo^s) = zero
144 end where
145 end if
146
147 end subroutine mhd_get_lcd
148
149 subroutine mhd_get_wcd(wLC,wRC,whll,fRC,fLC,Fhll,patchf,lambdaCD,cmin,cmax,&
150 ixI^L,ixO^L,idim,f)
151 ! compute the intermediate state U*
152 ! only needed where patchf=-1/1
153
154 ! reference Li S., JCP, 203, 2005, 344-357
155 ! reference T. Miyoski, Kusano JCP, 2008, 2005.
157
158 integer, intent(in) :: ixi^l,ixo^l,idim
159 double precision, dimension(ixI^S,1:nw), intent(in) :: wrc,wlc
160 double precision, dimension(ixI^S,1:nwflux), intent(in):: whll, fhll
161 double precision, dimension(ixI^S), intent(in) :: lambdacd
162 double precision, dimension(ixI^S), intent(in) :: cmax,cmin
163 double precision, dimension(ixI^S,1:nwflux), intent(in):: frc,flc
164 double precision, dimension(ixI^S,1:nwflux),intent(out):: f
165 double precision, dimension(ixI^S,1:nw) :: wcd,wsub
166 double precision, dimension(ixI^S,1:nwflux) :: fsub
167 double precision, dimension(ixI^S) :: vsub,cspeed,pcd,vdotbcd
168 integer , dimension(ixI^S), intent(in) :: patchf
169
170 integer :: n, iw, idir,ix^d
171
172 !-------------- auxiliary Speed and array-------------!
173 {do ix^db=ixomin^db,ixomax^db\}
174 if(patchf(ix^d)==1) then
175 cspeed(ix^d)=cmax(ix^d)
176 vsub(ix^d)=wrc(ix^d,mom(idim))/wrc(ix^d,rho_)
177 wsub(ix^d,:)=wrc(ix^d,:)
178 fsub(ix^d,:)=frc(ix^d,:)
179 else if(patchf(ix^d)==-1) then
180 cspeed(ix^d)=cmin(ix^d)
181 vsub(ix^d)=wlc(ix^d,mom(idim))/wlc(ix^d,rho_)
182 wsub(ix^d,:)=wlc(ix^d,:)
183 fsub(ix^d,:)=flc(ix^d,:)
184 end if
185 {end do\}
186
187 {do ix^db=ixomin^db,ixomax^db\}
188 if(abs(patchf(ix^d))==1) then
189 wcd(ix^d,rho_) = wsub(ix^d,rho_)&
190 *(cspeed(ix^d)-vsub(ix^d))/(cspeed(ix^d)-lambdacd(ix^d))
191 do n=1,mhd_n_tracer
192 iw = tracer(n)
193 wcd(ix^d,iw) = wsub(ix^d,iw)*(cspeed(ix^d)-vsub(ix^d))&
194 /(cspeed(ix^d)-lambdacd(ix^d))
195 end do
196 !==== Magnetic field ====!
197 do idir=1,ndir
198 ! case from eq 31
199 wcd(ix^d,mag(idir)) = whll(ix^d,mag(idir))
200 end do
201 !------- Momentum ------!
202 do iw=1, ndir
203 if(iw /= idim)then
204 ! eq. 21 22
205 wcd(ix^d,mom(iw))=(cspeed(ix^d)*wsub(ix^d,mom(iw))-fsub(ix^d,mom(iw))&
206 -wcd(ix^d,mag(idim))*wcd(ix^d,mag(iw)))/&
207 (cspeed(ix^d)-lambdacd(ix^d))
208 else
209 ! eq. 20
210 wcd(ix^d,mom(iw)) = wcd(ix^d,rho_) * lambdacd(ix^d)
211 endif
212 enddo
213 if(mhd_energy) then
214 vdotbcd(ix^d) = sum(whll(ix^d,mom(:))*whll(ix^d,mag(:)))/whll(ix^d,rho_)
215 ! Eq 17
216 pcd(ix^d) = wsub(ix^d,rho_)*(cspeed(ix^d)-vsub(ix^d))&
217 *(lambdacd(ix^d)-vsub(ix^d))&
218 +fsub(ix^d,mom(idim))-wsub(ix^d,mom(idim))*vsub(ix^d)&
219 + wcd(ix^d,mag(idim))**2
220 ! Eq 31
221 wcd(ix^d,e_) = (cspeed(ix^d)*wsub(ix^d,e_) &
222 -fsub(ix^d,e_)+lambdacd(ix^d)*pcd(ix^d)&
223 -vdotbcd(ix^d)*wcd(ix^d,mag(idim)))&
224 /(cspeed(ix^d)-lambdacd(ix^d))
225 end if
226 end if
227 {end do\}
228
229 do iw=1,nwflux
230 if(iw == mag(idim)) then
231 f(ixo^s,iw)=zero
232 else if(mhd_glm .and. iw == psi_) then
233 f(ixo^s,iw)=zero
234 else
235 where(abs(patchf(ixo^s))==1)
236 ! f_i=fsub+lambda (wCD-wSub)
237 f(ixo^s,iw)=fsub(ixo^s,iw)+cspeed(ixo^s)*(wcd(ixo^s,iw)-wsub(ixo^s,iw))
238 endwhere
239 end if
240 end do
241
242 end subroutine mhd_get_wcd
243
244end module mod_mhd_hllc
integer, dimension(:), allocatable, public mag
Indices of the magnetic field.
This module contains definitions of global parameters and variables and some generic functions/subrou...
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision, dimension(:), allocatable, parameter d
subroutine, public mhd_hllc_init()
Magneto-hydrodynamics module.
Definition mod_mhd_phys.t:2
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
integer, public, protected rho_
Index of the density (in the w array)
procedure(sub_get_wcd), pointer phys_get_wcd
procedure(sub_get_lcd), pointer phys_get_lcd
procedure(sub_diffuse_hllcd), pointer phys_diffuse_hllcd