MPI-AMRVAC  3.0
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_twofl_roe.t
Go to the documentation of this file.
1 !> Subroutines for Roe-type Riemann solver for HD
3 #include "amrvac.h"
4 
7 
8  implicit none
9  private
10 
11  integer, parameter :: fastRW_ = 3,fastlw_=4,slowrw_=5,slowlw_=6 ! Characteristic
12  integer, parameter :: entroW_ = 8,diverw_=7,alfvrw_=1,alfvlw_=2 ! waves
13 
14  public :: twofl_roe_init
15 
16 contains
17 
18  subroutine twofl_roe_init()
20  integer :: il
21 
25 
26  nworkroe=15
27  allocate(entropycoef(nw))
28 
29  do il = 1, nw
30  select case(il)
31  case(fastrw_,fastlw_,slowrw_,slowlw_)
32  entropycoef(il) = 0.2d0
33  case(alfvrw_,alfvlw_)
34  entropycoef(il) = 0.4d0
35  case default
36  entropycoef(il) = -1.0d0
37  end select
38  end do
39 
40  end subroutine twofl_roe_init
41 
42  ! Eight-wave MHD Riemann solver. See Powell, Notes on the eigensystem, Gombosi
43  ! Calculate the wroe average of primitive variables in wL and wR, assignment:
44  ! rho -> sqrho, m -> v, e -> p, B_idim -> B_idim, B_idir -> beta_idir
45  ! Calculate also alpha_f,alpha_s,c_f,c_s,csound2,dp,rhodv
46  !
47  ! wL,wR,wroe are all interface centered quantities
48  subroutine twofl_average(wL,wR,x,ix^L,idim,wroe,workroe)
50 
51  integer, intent(in) :: ix^L, idim
52  double precision, intent(in) :: wL(ixG^T, nw), wR(ixG^T, nw)
53  double precision, intent(inout) :: wroe(ixG^T, nw)
54  double precision, intent(inout) :: workroe(ixG^T, nworkroe)
55  double precision, intent(in) :: x(ixG^T, 1:^ND)
56 
57  call average2(wl,wr,x,ix^l,idim,wroe,workroe(ixg^t,1),workroe(ixg^t,2), &
58  workroe(ixg^t,3),workroe(ixg^t,4),workroe(ixg^t,5),workroe(ixg^t,6), &
59  workroe(ixg^t,7),workroe(ixg^t,8))
60 
61  end subroutine twofl_average
62 
63  ! Eight-wave MHD Riemann solver. See Powell, Notes on the eigensystem, Gombosi
64  ! Calculate the wroe average of primitive variables in wL and wR, assignment:
65  ! rho -> sqrho, m -> v, e -> p, B_idim -> B_idim, B_idir -> beta_idir
66  ! Calculate also alpha_f,alpha_s,c_f,c_s,csound2,dp,rhodv
67  !
68  ! wL,wR,wroe are all interface centered quantities
69  subroutine average2(wL,wR,x,ix^L,idim,wroe,cfast,cslow,afast,aslow,csound2,dp, &
70  rhodv,tmp)
72  integer :: ix^L,idim,idir,jdir,iw
73  double precision, dimension(ixG^T,nw) :: wL,wR,wroe
74  double precision, intent(in) :: x(ixG^T,1:ndim)
75  double precision, dimension(ixG^T) :: cfast,cslow,afast,aslow,csound2,dp, &
76  rhodv,tmp
77 
78  if (ndir==1) call mpistop("twofl with d=11 is the same as HD")
79 
80  !Averaging primitive variables
81  wroe(ix^s,rho_c_)=half*(wl(ix^s,rho_c_)+wr(ix^s,rho_c_))
82  do idir=1,ndir
83  wroe(ix^s,mom_c(idir))=half*(wl(ix^s,mom_c(idir))/wl(ix^s,rho_c_)+wr(ix^s,mom_c(idir))/wr(ix^s,rho_c_))
84  wroe(ix^s,mag(idir))=half*(wl(ix^s,mag(idir))+wr(ix^s,mag(idir)))
85  end do
86  ! Use afast and aslow for pressures pL and pR
87  call twofl_get_pthermal_c(wl,x,ixg^ll,ix^l,afast)
88  call twofl_get_pthermal_c(wr,x,ixg^ll,ix^l,aslow)
89 
90  if(twofl_eq_energy/=0) then
91  wroe(ix^s,e_c_)=half*(afast(ix^s)+aslow(ix^s))
92  ! dp=pR-pL
93  dp(ix^s)=aslow(ix^s)-afast(ix^s)
94  else
95  dp(ix^s)=aslow(ix^s)-afast(ix^s)
96  end if
97 
98  !CONSERVATIVE rho*dv_idim=dm_idim-v_idim*drho
99  rhodv(ix^s)=wr(ix^s,mom_c(idim))-wl(ix^s,mom_c(idim))-&
100  wroe(ix^s,mom_c(idim))*(wr(ix^s,rho_c_)-wl(ix^s,rho_c_))
101 
102  !Calculate csound2,cfast,cslow,alphafast and alphaslow
103 
104  ! get csound**2
105  if(twofl_eq_energy/=0) then
106  csound2(ix^s)=twofl_gamma*wroe(ix^s,e_c_)/wroe(ix^s,rho_c_)
107  else
108  csound2(ix^s)=twofl_gamma*twofl_adiab*wroe(ix^s,rho_c_)**(twofl_gamma-one)
109  end if
110 
111  ! aa=B**2/rho+a**2
112  cfast(ix^s)=sum(wroe(ix^s,mag(:))**2,dim=ndim+1)/wroe(ix^s,rho_c_)+csound2(ix^s)
113 
114  ! cs**2=0.5*(aa+dsqrt(aa**2-4*a**2*(b_i**2/rho)))
115  cslow(ix^s)=half*(cfast(ix^s)-dsqrt(cfast(ix^s)**2-&
116  4d0*csound2(ix^s)*wroe(ix^s,mag(idim))**2/wroe(ix^s,rho_c_)))
117 
118  ! cf**2=aa-cs**2
119  cfast(ix^s)=cfast(ix^s)-cslow(ix^s)
120 
121  ! alpha_f**2=(a**2-cs**2)/(cf**2-cs**2)
122  afast(ix^s)=(csound2(ix^s)-cslow(ix^s))/(cfast(ix^s)-cslow(ix^s))
123  afast(ix^s)=min(one,max(afast(ix^s),zero))
124 
125  ! alpha_s=dsqrt(1-alpha_f**2)
126  aslow(ix^s)=dsqrt(one-afast(ix^s))
127 
128  ! alpha_f=dsqrt(alpha_f**2)
129  afast(ix^s)=dsqrt(afast(ix^s))
130 
131  ! cf=dsqrt(cf**2)
132  cfast(ix^s)=dsqrt(cfast(ix^s))
133 
134  ! cs=dsqrt(cs**2)
135  cslow(ix^s)=dsqrt(cslow(ix^s))
136 
137  !Replace the primitive variables with more useful quantities:
138  ! rho -> dsqrt(rho)
139  wroe(ix^s,rho_c_)=dsqrt(wroe(ix^s,rho_c_))
140 
141  ! Avoid sgn(b_idim)==0
142  where(dabs(wroe(ix^s,mag(idim)))<smalldouble)&
143  wroe(ix^s,mag(idim))=smalldouble
144  ! B_idir,jdir -> beta_idir,jdir
145  idir=idim+1-ndir*(idim/ndir)
146  if(ndir==2)then
147  where(wroe(ix^s,mag(idir))>=zero)
148  wroe(ix^s,mag(idir))=one
149  elsewhere
150  wroe(ix^s,mag(idir))=-one
151  end where
152  else
153  !beta_j=B_j/dsqrt(B_i**2+B_j**2); beta_i=B_i/dsqrt(B_i**2+B_j**2)
154  jdir=idir+1-ndir*(idir/ndir)
155  tmp(ix^s)=dsqrt(wroe(ix^s,mag(idir))**2+wroe(ix^s,mag(jdir))**2)
156  where(tmp(ix^s)>smalldouble)
157  wroe(ix^s,mag(idir))=wroe(ix^s,mag(idir))/tmp(ix^s)
158  wroe(ix^s,mag(jdir))=wroe(ix^s,mag(jdir))/tmp(ix^s)
159  elsewhere
160  wroe(ix^s,mag(idir))=dsqrt(half)
161  wroe(ix^s,mag(jdir))=dsqrt(half)
162  end where
163  endif
164 
165  end subroutine average2
166 
167  ! Calculate the il-th characteristic speed and the jump in the il-th
168  ! characteristic variable in the idim direction within ixL.
169  ! The eigenvalues and the l=r**(-1) matrix is calculated from wroe.
170  ! jump(il)=Sum_il l(il,iw)*(wR(iw)-wL(iw)), where w are the conservative
171  ! variables. However part of the summation is done in advance and saved into
172  ! bdv,bdb,dp and dv variables. "smalla" contains a lower limit for "a" to be
173  ! used in the entropy fix.
174  !
175  ! All the variables are centered on the cell interface, thus the
176  ! "*C" notation is omitted for sake of brevity.
177  subroutine twofl_get_eigenjump(wL,wR,wroe,x,ix^L,il,idim,smalla,a,jump,workroe)
179 
180  integer, intent(in) :: ix^L,il,idim
181  double precision, dimension(ixG^T,nw) :: wL,wR,wroe
182  double precision, intent(in) :: x(ixG^T,1:ndim)
183  double precision, dimension(ixG^T) :: smalla,a,jump
184  double precision, dimension(ixG^T,nworkroe) :: workroe
185 
186  call geteigenjump2(wl,wr,wroe,x,ix^l,il,idim,smalla,a,jump, &
187  workroe(ixg^t,1),workroe(ixg^t,2), &
188  workroe(ixg^t,3),workroe(ixg^t,4),workroe(ixg^t,5),workroe(ixg^t,6), &
189  workroe(ixg^t,7),workroe(ixg^t,8),workroe(ixg^t,9),workroe(ixg^t,10), &
190  workroe(ixg^t,11),workroe(ixg^t,12),workroe(ixg^t,13))
191 
192  end subroutine twofl_get_eigenjump
193 
194  ! Calculate the il-th characteristic speed and the jump in the il-th
195  ! characteristic variable in the idim direction within ixL.
196  ! The eigenvalues and the l=r**(-1) matrix is calculated from wroe.
197  ! jump(il)=Sum_il l(il,iw)*(wR(iw)-wL(iw)), where w are the conservative
198  ! variables. However part of the summation is done in advance and saved into
199  ! bdv,bdb,dp and dv variables. "smalla" contains a lower limit for "a" to be
200  ! used in the entropy fix.
201  !
202  ! All the variables are centered on the cell interface, thus the
203  ! "*C" notation is omitted for sake of brevity.
204  subroutine geteigenjump2(wL,wR,wroe,x,ix^L,il,idim,smalla,a,jump, &
205  cfast,cslow,afast,aslow,csound2,dp,rhodv,bdv,bdb,cs2L,cs2R,cs2ca2L,cs2ca2R)
207  use mod_tvd
208 
209  integer :: ix^L,il,idim,idir,jdir
210  double precision, dimension(ixG^T,nw) :: wL,wR,wroe
211  double precision, intent(in) :: x(ixG^T,1:ndim)
212  double precision, dimension(ixG^T) :: smalla,a,jump
213  double precision, dimension(ixG^T) :: cfast,cslow,afast,aslow,csound2,dp,rhodv
214  double precision, dimension(ixG^T) :: bdv,bdb
215  double precision, dimension(ixG^T) :: aL,aR,cs2L,cs2R,cs2ca2L,cs2ca2R
216 
217  idir=idim+1-ndir*(idim/ndir)
218  jdir=idir+1-ndir*(idir/ndir)
219 
220  if(il==fastrw_)then
221  !Fast and slow waves use bdv=sqrho**2*sign(bx)*(betay*dvy+betaz*dvz)
222  ! bdb=sqrho*a* (betay*dBy+betaz*dBz)
223  bdv(ix^s)=wroe(ix^s,mag(idir))* &
224  (wr(ix^s,mom_c(idir))/wr(ix^s,rho_c_)-wl(ix^s,mom_c(idir))/wl(ix^s,rho_c_))
225  if(ndir==3)bdv(ix^s)=bdv(ix^s)+wroe(ix^s,mag(jdir))* &
226  (wr(ix^s,mom_c(jdir))/wr(ix^s,rho_c_)-wl(ix^s,mom_c(jdir))/wl(ix^s,rho_c_))
227  bdv(ix^s)=bdv(ix^s)*sign(wroe(ix^s,rho_c_)**2,wroe(ix^s,mag(idim)))
228 
229  bdb(ix^s)=wroe(ix^s,mag(idir))*(wr(ix^s,mag(idir))-wl(ix^s,mag(idir)))
230  if(ndir==3)bdb(ix^s)=bdb(ix^s)+&
231  wroe(ix^s,mag(jdir))*(wr(ix^s,mag(jdir))-wl(ix^s,mag(jdir)))
232  bdb(ix^s)=bdb(ix^s)*dsqrt(csound2(ix^s))*wroe(ix^s,rho_c_)
233  endif
234 
235  if(il==alfvrw_)then
236  !Alfven waves use bdv=0.5*sqrho**2* (betaz*dvy-betay*dvz)
237  ! bdb=0.5*sqrho*sign(bx)*(betaz*dBy-betay*dBz)
238  bdv(ix^s)=wroe(ix^s,mag(jdir))* &
239  (wr(ix^s,mom_c(idir))/wr(ix^s,rho_c_)-wl(ix^s,mom_c(idir))/wl(ix^s,rho_c_)) &
240  -wroe(ix^s,mag(idir))* &
241  (wr(ix^s,mom_c(jdir))/wr(ix^s,rho_c_)-wl(ix^s,mom_c(jdir))/wl(ix^s,rho_c_))
242  bdb(ix^s)=wroe(ix^s,mag(jdir))*(wr(ix^s,mag(idir))-wl(ix^s,mag(idir))) &
243  -wroe(ix^s,mag(idir))*(wr(ix^s,mag(jdir))-wl(ix^s,mag(jdir)))
244  bdv(ix^s)=bdv(ix^s)*half*wroe(ix^s,rho_c_)**2
245  bdb(ix^s)=bdb(ix^s)*half*sign(wroe(ix^s,rho_c_),wroe(ix^s,mag(idim)))
246  endif
247 
248  select case(il)
249  case(fastrw_)
250  a(ix^s)=wroe(ix^s,mom_c(idim))+cfast(ix^s)
251  jump(ix^s)=half/csound2(ix^s)*(&
252  afast(ix^s)*(+cfast(ix^s)*rhodv(ix^s)+dp(ix^s))&
253  +aslow(ix^s)*(-cslow(ix^s)*bdv(ix^s)+bdb(ix^s)))
254  case(fastlw_)
255  a(ix^s)=wroe(ix^s,mom_c(idim))-cfast(ix^s)
256  jump(ix^s)=half/csound2(ix^s)*(&
257  afast(ix^s)*(-cfast(ix^s)*rhodv(ix^s)+dp(ix^s))&
258  +aslow(ix^s)*(+cslow(ix^s)*bdv(ix^s)+bdb(ix^s)))
259  case(slowrw_)
260  a(ix^s)=wroe(ix^s,mom_c(idim))+cslow(ix^s)
261  jump(ix^s)=half/csound2(ix^s)*(&
262  aslow(ix^s)*(+cslow(ix^s)*rhodv(ix^s)+dp(ix^s))&
263  +afast(ix^s)*(+cfast(ix^s)*bdv(ix^s)-bdb(ix^s)))
264  case(slowlw_)
265  a(ix^s)=wroe(ix^s,mom_c(idim))-cslow(ix^s)
266  jump(ix^s)=half/csound2(ix^s)*(&
267  aslow(ix^s)*(-cslow(ix^s)*rhodv(ix^s)+dp(ix^s))&
268  +afast(ix^s)*(-cfast(ix^s)*bdv(ix^s)-bdb(ix^s)))
269  case(entrow_)
270  a(ix^s)=wroe(ix^s,mom_c(idim))
271  jump(ix^s)=wr(ix^s,rho_c_)-wl(ix^s,rho_c_)-dp(ix^s)/csound2(ix^s)
272  case(diverw_)
273  if(divbwave)then
274  a(ix^s)=wroe(ix^s,mom_c(idim))
275  jump(ix^s)=wr(ix^s,mag(idim))-wl(ix^s,mag(idim))
276  else
277  a(ix^s)=zero
278  jump(ix^s)=zero
279  endif
280  case(alfvrw_)
281  a(ix^s)=wroe(ix^s,mom_c(idim))+dabs(wroe(ix^s,mag(idim)))/wroe(ix^s,rho_c_)
282  jump(ix^s)=+bdv(ix^s)-bdb(ix^s)
283  case(alfvlw_)
284  a(ix^s)=wroe(ix^s,mom_c(idim))-dabs(wroe(ix^s,mag(idim)))/wroe(ix^s,rho_c_)
285  jump(ix^s)=-bdv(ix^s)-bdb(ix^s)
286  end select
287 
288  ! Calculate "smalla" or modify "a" based on the "typeentropy" switch
289 
290  select case(typeentropy(il))
291  case('yee')
292  ! Based on Yee JCP 68,151 eq 3.23
293  smalla(ix^s)=entropycoef(il)
294  case('harten','powell', 'ratio')
295  ! Based on Harten & Hyman JCP 50, 235 and Zeeuw & Powell JCP 104,56
296  ! Initialize left and right eigenvalues by velocities
297  al(ix^s)= wl(ix^s,mom_c(idim))/wl(ix^s,rho_c_)
298  ar(ix^s)= wr(ix^s,mom_c(idim))/wr(ix^s,rho_c_)
299  ! Calculate the final "aL" and "aR"
300  select case(il)
301  case(fastrw_)
302  ! These quantities will be used for all the fast and slow waves
303  ! Calculate soundspeed**2 and cs**2+ca**2.
304  call twofl_get_csound2_c_from_conserved(wl,x,ixg^ll,ix^l,cs2l)
305  call twofl_get_csound2_c_from_conserved(wr,x,ixg^ll,ix^l,cs2r)
306  cs2ca2l(ix^s)=cs2l(ix^s)+sum(wl(ix^s,mag(:))**2,dim=ndim+1)/wl(ix^s,rho_c_)
307  cs2ca2r(ix^s)=cs2r(ix^s)+sum(wr(ix^s,mag(:))**2,dim=ndim+1)/wr(ix^s,rho_c_)
308  ! Save the discriminants into cs2L and cs2R
309  cs2l(ix^s)=&
310  dsqrt(cs2ca2l(ix^s)**2-4d0*cs2l(ix^s)*wl(ix^s,mag(idim))**2/wl(ix^s,rho_c_))
311  cs2r(ix^s)=&
312  dsqrt(cs2ca2r(ix^s)**2-4d0*cs2r(ix^s)*wr(ix^s,mag(idim))**2/wr(ix^s,rho_c_))
313 
314  ! The left and right eigenvalues for the fast wave going to right
315  al(ix^s)=al(ix^s) + dsqrt(half*(cs2ca2l(ix^s) + cs2l(ix^s)))
316  ar(ix^s)=ar(ix^s) + dsqrt(half*(cs2ca2r(ix^s) + cs2r(ix^s)))
317  case(fastlw_)
318  al(ix^s)=al(ix^s) - dsqrt(half*(cs2ca2l(ix^s) + cs2l(ix^s)))
319  ar(ix^s)=ar(ix^s) - dsqrt(half*(cs2ca2r(ix^s) + cs2r(ix^s)))
320  case(slowrw_)
321  al(ix^s)=al(ix^s) + dsqrt(half*(cs2ca2l(ix^s) - cs2l(ix^s)))
322  ar(ix^s)=ar(ix^s) + dsqrt(half*(cs2ca2r(ix^s) - cs2r(ix^s)))
323  case(slowlw_)
324  al(ix^s)=al(ix^s) - dsqrt(half*(cs2ca2l(ix^s) - cs2l(ix^s)))
325  ar(ix^s)=ar(ix^s) - dsqrt(half*(cs2ca2r(ix^s) - cs2r(ix^s)))
326  case(entrow_,diverw_)
327  ! These propagate by the velocity
328  case(alfvrw_)
329  ! Store the Alfven speeds into cs2ca2L and cs2ca2R
330  cs2ca2l(ix^s)=dabs(wl(ix^s,mag(idim)))/dsqrt(wl(ix^s,rho_c_))
331  cs2ca2r(ix^s)=dabs(wr(ix^s,mag(idim)))/dsqrt(wr(ix^s,rho_c_))
332 
333  al(ix^s)=al(ix^s) + cs2ca2l(ix^s)
334  ar(ix^s)=ar(ix^s) + cs2ca2r(ix^s)
335  case(alfvlw_)
336  al(ix^s)=al(ix^s) - cs2ca2l(ix^s)
337  ar(ix^s)=ar(ix^s) - cs2ca2r(ix^s)
338  end select
339  end select
340 
341  call entropyfix(ix^l,il,al,ar,a,smalla)
342 
343  end subroutine geteigenjump2
344 
345  ! Multiply q by R(il,iw), where R is the right eigenvalue matrix at wroe
346  subroutine twofl_rtimes(q,w,ix^L,iw,il,idim,rq,workroe)
348 
349  integer, intent(in) :: ix^L, iw, il, idim
350  double precision, intent(in) :: w(ixG^T, nw), q(ixG^T)
351  double precision, intent(inout) :: rq(ixG^T)
352  double precision, intent(inout) :: workroe(ixG^T, nworkroe)
353 
354  call rtimes2(q,w,ix^l,iw,il,idim,rq,&
355  workroe(ixg^t,1),workroe(ixg^t,2), &
356  workroe(ixg^t,3),workroe(ixg^t,4),workroe(ixg^t,5),workroe(ixg^t,6), &
357  workroe(ixg^t,7),workroe(ixg^t,14),workroe(ixg^t,15))
358 
359  end subroutine twofl_rtimes
360 
361  ! Multiply q by R(il,iw), where R is the right eigenvalue matrix at wroe
362  subroutine rtimes2(q,wroe,ix^L,iw,il,idim,rq, &
363  cfast,cslow,afast,aslow,csound2,dp,rhodv,bv,v2a2)
365 
366  integer :: ix^L,iw,il,idim,idir,jdir
367  double precision :: wroe(ixG^T,nw)
368  double precision, dimension(ixG^T) :: q,rq
369  double precision, dimension(ixG^T) :: cfast,cslow,afast,aslow,csound2,dp,rhodv
370  double precision, dimension(ixG^T) :: bv,v2a2
371 
372  idir=idim+1-ndir*(idim/ndir)
373  jdir=idir+1-ndir*(idir/ndir)
374 
375  if(iw == rho_c_) then
376  select case(il)
377  case(fastrw_,fastlw_)
378  rq(ix^s)=q(ix^s)*afast(ix^s)
379  case(slowrw_,slowlw_)
380  rq(ix^s)=q(ix^s)*aslow(ix^s)
381  case(entrow_)
382  rq(ix^s)=q(ix^s)
383  case(diverw_,alfvrw_,alfvlw_)
384  rq(ix^s)=zero
385  end select
386  else if(iw == e_c_) then
387  if(il==fastrw_)then
388  ! Store 0.5*v**2+(2-gamma)/(gamma-1)*a**2
389  v2a2(ix^s)=half*sum(wroe(ix^s,mom_c(:))**2,dim=ndim+1)+ &
390  (two-twofl_gamma)/(twofl_gamma-one)*csound2(ix^s)
391  ! Store sgn(bx)*(betay*vy+betaz*vz) in bv
392  bv(ix^s)=wroe(ix^s,mag(idir))*wroe(ix^s,mom_c(idir))
393  if(ndir==3)bv(ix^s)=bv(ix^s)+wroe(ix^s,mag(jdir))*wroe(ix^s,mom_c(jdir))
394  bv(ix^s)=bv(ix^s)*sign(one,wroe(ix^s,mag(idim)))
395  else if(il==alfvrw_)then
396  !Store betaz*vy-betay*vz in bv
397  bv(ix^s)=(wroe(ix^s,mag(jdir))*wroe(ix^s,mom_c(idir))-&
398  wroe(ix^s,mag(idir))*wroe(ix^s,mom_c(jdir)))
399  endif
400 
401  select case(il)
402  case(fastrw_)
403  rq(ix^s)=q(ix^s)*(-aslow(ix^s)*cslow(ix^s)*bv(ix^s)+afast(ix^s)*&
404  (v2a2(ix^s)+cfast(ix^s)*(cfast(ix^s)+wroe(ix^s,mom_c(idim)))))
405  case(fastlw_)
406  rq(ix^s)=q(ix^s)*(+aslow(ix^s)*cslow(ix^s)*bv(ix^s)+afast(ix^s)*&
407  (v2a2(ix^s)+cfast(ix^s)*(cfast(ix^s)-wroe(ix^s,mom_c(idim)))))
408  case(slowrw_)
409  rq(ix^s)=q(ix^s)*(+afast(ix^s)*cfast(ix^s)*bv(ix^s)+aslow(ix^s)*&
410  (v2a2(ix^s)+cslow(ix^s)*(cslow(ix^s)+wroe(ix^s,mom_c(idim)))))
411  case(slowlw_)
412  rq(ix^s)=q(ix^s)*(-afast(ix^s)*cfast(ix^s)*bv(ix^s)+aslow(ix^s)*&
413  (v2a2(ix^s)+cslow(ix^s)*(cslow(ix^s)-wroe(ix^s,mom_c(idim)))))
414  case(entrow_)
415  rq(ix^s)= q(ix^s)*half*sum(wroe(ix^s,mom_c(:))**2,dim=ndim+1)
416  case(diverw_)
417  if(divbwave)then
418  rq(ix^s)= q(ix^s)*wroe(ix^s,mag(idim))
419  else
420  rq(ix^s)= zero
421  endif
422  case(alfvrw_)
423  rq(ix^s)=+q(ix^s)*bv(ix^s)
424  case(alfvlw_)
425  rq(ix^s)=-q(ix^s)*bv(ix^s)
426  end select
427  else if(any(mom_c(:)==iw)) then
428  if(iw==mom_c(idim))then
429  select case(il)
430  case(fastrw_)
431  rq(ix^s)=q(ix^s)*afast(ix^s)*(wroe(ix^s,iw)+cfast(ix^s))
432  case(fastlw_)
433  rq(ix^s)=q(ix^s)*afast(ix^s)*(wroe(ix^s,iw)-cfast(ix^s))
434  case(slowrw_)
435  rq(ix^s)=q(ix^s)*aslow(ix^s)*(wroe(ix^s,iw)+cslow(ix^s))
436  case(slowlw_)
437  rq(ix^s)=q(ix^s)*aslow(ix^s)*(wroe(ix^s,iw)-cslow(ix^s))
438  case(entrow_)
439  rq(ix^s)=q(ix^s)*wroe(ix^s,iw)
440  case(diverw_,alfvlw_,alfvrw_)
441  rq(ix^s)=zero
442  end select
443  else
444  select case(il)
445  case(fastrw_)
446  rq(ix^s)=q(ix^s)*(afast(ix^s)*wroe(ix^s,iw)-aslow(ix^s)*&
447  cslow(ix^s)*wroe(ix^s,mag(1)-mom_c(1)+iw)*sign(one,wroe(ix^s,mag(idim))))
448  case(fastlw_)
449  rq(ix^s)=q(ix^s)*(afast(ix^s)*wroe(ix^s,iw)+aslow(ix^s)*&
450  cslow(ix^s)*wroe(ix^s,mag(1)-mom_c(1)+iw)*sign(one,wroe(ix^s,mag(idim))))
451  case(slowrw_)
452  rq(ix^s)=q(ix^s)*(aslow(ix^s)*wroe(ix^s,iw)+afast(ix^s)*&
453  cfast(ix^s)*wroe(ix^s,mag(1)-mom_c(1)+iw)*sign(one,wroe(ix^s,mag(idim))))
454  case(slowlw_)
455  rq(ix^s)=q(ix^s)*(aslow(ix^s)*wroe(ix^s,iw)-afast(ix^s)*&
456  cfast(ix^s)*wroe(ix^s,mag(1)-mom_c(1)+iw)*sign(one,wroe(ix^s,mag(idim))))
457  case(entrow_)
458  rq(ix^s)=q(ix^s)*wroe(ix^s,iw)
459  case(diverw_)
460  rq(ix^s)=zero
461  case(alfvrw_)
462  if(iw==mom_c(idir))then
463  rq(ix^s)=+q(ix^s)*wroe(ix^s,mag(jdir))
464  else
465  rq(ix^s)=-q(ix^s)*wroe(ix^s,mag(idir))
466  endif
467  case(alfvlw_)
468  if(iw==mom_c(idir))then
469  rq(ix^s)=-q(ix^s)*wroe(ix^s,mag(jdir))
470  else
471  rq(ix^s)=+q(ix^s)*wroe(ix^s,mag(idir))
472  endif
473  end select
474  end if ! iw=m_idir,m_jdir
475  else if(any(mag(:)==iw)) then
476  if(iw==mag(idim))then
477  if(il==diverw_ .and. divbwave)then
478  rq(ix^s)=q(ix^s)
479  else
480  rq(ix^s)=zero
481  endif
482  else
483  select case(il)
484  case(fastrw_,fastlw_)
485  rq(ix^s)=+q(ix^s)*aslow(ix^s)*dsqrt(csound2(ix^s))*wroe(ix^s,iw)&
486  /wroe(ix^s,rho_c_)
487  case(slowrw_,slowlw_)
488  rq(ix^s)=-q(ix^s)*afast(ix^s)*dsqrt(csound2(ix^s))*wroe(ix^s,iw)&
489  /wroe(ix^s,rho_c_)
490  case(entrow_,diverw_)
491  rq(ix^s)=zero
492  case(alfvrw_,alfvlw_)
493  if(iw==mag(idir))then
494  rq(ix^s)=-q(ix^s)*wroe(ix^s,mag(jdir))&
495  /sign(wroe(ix^s,rho_c_),wroe(ix^s,mag(idim)))
496  else
497  rq(ix^s)=+q(ix^s)*wroe(ix^s,mag(idir))&
498  /sign(wroe(ix^s,rho_c_),wroe(ix^s,mag(idim)))
499  end if
500  end select
501  end if ! iw=b_idir,b_jdir
502  end if
503 
504  end subroutine rtimes2
505 
506 end module mod_twofl_roe
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: comm_lib.t:194
This module contains definitions of global parameters and variables and some generic functions/subrou...
character(len=std_len), dimension(:), allocatable typeentropy
Which type of entropy fix to use with Riemann-type solvers.
integer, parameter ndim
Number of spatial dimensions for grid variables.
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision, dimension(:), allocatable entropycoef
procedure(sub_rtimes), pointer phys_rtimes
procedure(sub_get_eigenjump), pointer phys_get_eigenjump
procedure(sub_average), pointer phys_average
Subroutines for TVD-MUSCL schemes.
Definition: mod_tvd.t:2
subroutine, public entropyfix(ixL, il, aL, aR, a, smalla)
Definition: mod_tvd.t:264
Magneto-hydrodynamics module.
Definition: mod_twofl_phys.t:2
double precision, public twofl_adiab
The adiabatic constant.
subroutine, public twofl_get_csound2_c_from_conserved(w, x, ixIL, ixOL, csound2)
subroutine, public twofl_get_pthermal_c(w, x, ixIL, ixOL, pth)
integer, public e_c_
Index of the energy density (-1 if not present)
integer, dimension(:), allocatable, public mom_c
Indices of the momentum density.
logical, public divbwave
Add divB wave in Roe solver.
integer, dimension(:), allocatable, public mag
Indices of the magnetic field.
integer, public rho_c_
Index of the density (in the w array)
integer, public, protected twofl_eq_energy
double precision, public twofl_gamma
The adiabatic index.
Subroutines for Roe-type Riemann solver for HD.
Definition: mod_twofl_roe.t:2
subroutine average2(wL, wR, x, ixL, idim, wroe, cfast, cslow, afast, aslow, csound2, dp, rhodv, tmp)
Definition: mod_twofl_roe.t:71
subroutine twofl_get_eigenjump(wL, wR, wroe, x, ixL, il, idim, smalla, a, jump, workroe)
subroutine twofl_rtimes(q, w, ixL, iw, il, idim, rq, workroe)
subroutine twofl_average(wL, wR, x, ixL, idim, wroe, workroe)
Definition: mod_twofl_roe.t:49
subroutine geteigenjump2(wL, wR, wroe, x, ixL, il, idim, smalla, a, jump, cfast, cslow, afast, aslow, csound2, dp, rhodv, bdv, bdb, cs2L, cs2R, cs2ca2L, cs2ca2R)
subroutine rtimes2(q, wroe, ixL, iw, il, idim, rq, cfast, cslow, afast, aslow, csound2, dp, rhodv, bv, v2a2)
subroutine, public twofl_roe_init()
Definition: mod_twofl_roe.t:19