9 integer :: soundRW_ = -1
10 integer :: soundLW_ = -1
11 integer :: entropW_ = -1
12 integer :: shearW0_ = -1
49 if (iw == soundrw_ .or. iw == soundlw_)
then
63 integer,
intent(in) :: ix^L, idim
64 double precision,
intent(in) :: wL(ixG^T, nw), wR(ixG^T, nw)
65 double precision,
intent(inout) :: wroe(ixG^T, nw)
66 double precision,
intent(inout) :: workroe(ixG^T, nworkroe)
67 double precision,
intent(in) :: x(ixG^T, 1:^ND)
71 workroe(ix^s, 1) = sqrt(wl(ix^s,
rho_))
72 workroe(ix^s, 2) = sqrt(wr(ix^s,
rho_))
75 wroe(ix^s,
rho_) = workroe(ix^s, 1)*workroe(ix^s, 2)
78 workroe(ix^s, 1) = workroe(ix^s, 1)/workroe(ix^s, 2)
82 wroe(ix^s,
mom(idir)) = (wl(ix^s,
mom(idir))/wl(ix^s,
rho_) * workroe(ix^s, 1)+&
83 wr(ix^s,
mom(idir))/wr(ix^s,
rho_))/(one+workroe(ix^s, 1))
89 wroe(ix^s,
e_) = (workroe(ix^s, 2)+wl(ix^s,
e_))/wl(ix^s,
rho_)
93 workroe(ix^s, 2) = (workroe(ix^s, 2)+wr(ix^s,
e_))/wr(ix^s,
rho_)
94 wroe(ix^s,
e_) = (wroe(ix^s,
e_)*workroe(ix^s, 1) + workroe(ix^s, 2))/(one+workroe(ix^s, 1))
97 subroutine average2(wL,wR,x,ix^L,idim,wroe,tmp,tmp2)
104 integer :: ix^L,idim,idir
105 double precision,
dimension(ixG^T,nw) :: wL,wR,wroe
106 double precision,
dimension(ixG^T,ndim),
intent(in) :: x
107 double precision,
dimension(ixG^T) :: tmp,tmp2
117 subroutine rhd_get_eigenjump(wL,wR,wroe,x,ix^L,il,idim,smalla,a,jump,workroe)
120 integer,
intent(in) :: ix^L,il,idim
121 double precision,
dimension(ixG^T,nw):: wL,wR,wroe
122 double precision,
dimension(ixG^T,ndim),
intent(in) :: x
123 double precision,
dimension(ixG^T) :: smalla,a,jump
124 double precision,
dimension(ixG^T,nworkroe) :: workroe
126 call geteigenjump2(wl,wr,wroe,x,ix^l,il,idim,smalla,a,jump, &
127 workroe(ixg^t,1),workroe(ixg^t,2),workroe(ixg^t,3))
131 subroutine geteigenjump2(wL,wR,wroe,x,ix^L,il,idim,smalla,a,jump,&
132 csound,dpperc2,dvperc)
142 integer :: ix^L,il,idim,idir
143 double precision,
dimension(ixG^T,nw) :: wL,wR,wroe
144 double precision,
dimension(ixG^T,ndim),
intent(in) :: x
145 double precision,
dimension(ixG^T) :: smalla,a,jump,tmp,tmp2
146 double precision,
dimension(ixG^T) :: csound,dpperc2,dvperc
147 double precision :: kin_en(ixG^T)
151 kin_en(ix^s) = 0.5d0 * sum(wroe(ix^s,
mom(:))**2, dim=^nd+1)
152 csound(ix^s)=(
rhd_gamma-one)*(wroe(ix^s,
e_) - kin_en(ix^s))
154 csound(ix^s)=max(
rhd_gamma*smalldouble/wroe(ix^s,
rho_),csound(ix^s))
160 dpperc2(ix^s)=(tmp2(ix^s)-tmp(ix^s))/csound(ix^s)
163 csound(ix^s)=sqrt(csound(ix^s))
166 dvperc(ix^s)=(wr(ix^s,
mom(idim))/wr(ix^s,
rho_)-&
167 wl(ix^s,
mom(idim))/wl(ix^s,
rho_))/csound(ix^s)
171 if (il == soundrw_)
then
172 a(ix^s)=wroe(ix^s,
mom(idim))+csound(ix^s)
173 jump(ix^s)=half*(dpperc2(ix^s)+wroe(ix^s,
rho_)*dvperc(ix^s))
174 else if (il == soundlw_)
then
175 a(ix^s)=wroe(ix^s,
mom(idim))-csound(ix^s)
176 jump(ix^s)=half*(dpperc2(ix^s)-wroe(ix^s,
rho_)*dvperc(ix^s))
177 else if (il == entropw_)
then
178 a(ix^s)=wroe(ix^s,
mom(idim))
179 jump(ix^s)=-dpperc2(ix^s)+wr(ix^s,
rho_)-wl(ix^s,
rho_)
182 idir=il-shearw0_;
if(idir>=idim)idir=idir+1
183 a(ix^s)=wroe(ix^s,
mom(idim))
184 jump(ix^s)=wroe(ix^s,
rho_)*&
185 (wr(ix^s,
mom(idir))/wr(ix^s,
rho_)-wl(ix^s,
mom(idir))/wl(ix^s,
rho_))
196 case(
'harten',
'powell')
198 if (il == soundrw_)
then
200 tmp(ix^s)=wl(ix^s,
mom(idim))/wl(ix^s,
rho_)&
203 tmp2(ix^s)=wr(ix^s,
mom(idim))/wr(ix^s,
rho_)&
205 else if (il == soundlw_)
then
207 tmp(ix^s)=wl(ix^s,
mom(idim))/wl(ix^s,
rho_)&
210 tmp2(ix^s)=wr(ix^s,
mom(idim))/wr(ix^s,
rho_)&
213 tmp(ix^s) =wl(ix^s,
mom(idim))/wl(ix^s,
rho_)
214 tmp2(ix^s)=wr(ix^s,
mom(idim))/wr(ix^s,
rho_)
226 integer,
intent(in) :: ix^L,iw,il,idim
227 double precision,
intent(in) :: wroe(ixG^T,nw)
228 double precision,
intent(in) :: q(ixG^T)
229 double precision,
intent(inout) :: rq(ixG^T)
230 double precision,
intent(inout) :: workroe(ixG^T,nworkroe)
232 call rtimes2(q,wroe,ix^l,iw,il,idim,rq,workroe(ixg^t,1))
235 subroutine rtimes2(q,wroe,ix^L,iw,il,idim,rq,csound)
241 integer :: ix^L,iw,il,idim,idir
242 double precision :: wroe(ixG^T,nw)
243 double precision,
dimension(ixG^T) :: q,rq,csound
246 shearwave=il>shearw0_
250 idir=il-shearw0_;
if(idir>=idim)idir=idir+1
259 else if (iw ==
e_)
then
260 if (il == soundrw_)
then
261 rq(ix^s)=q(ix^s)*(wroe(ix^s,
e_)+wroe(ix^s,
mom(idim))*csound(ix^s))
262 else if (il == soundlw_)
then
263 rq(ix^s)=q(ix^s)*(wroe(ix^s,
e_)-wroe(ix^s,
mom(idim))*csound(ix^s))
264 else if (il == entropw_)
then
265 rq(ix^s)=q(ix^s) * 0.5d0 * sum(wroe(ix^s,
mom(:))**2, dim=^nd+1)
267 rq(ix^s)=q(ix^s)*wroe(ix^s,
mom(idir))
270 if(iw==
mom(idim))
then
271 if (il == soundrw_)
then
272 rq(ix^s)=q(ix^s)*(wroe(ix^s,
mom(idim))+csound(ix^s))
273 else if (il == soundlw_)
then
274 rq(ix^s)=q(ix^s)*(wroe(ix^s,
mom(idim))-csound(ix^s))
275 else if (il == entropw_)
then
276 rq(ix^s)=q(ix^s)*wroe(ix^s,
mom(idim))
282 if(iw==
mom(idir))
then
288 rq(ix^s)=q(ix^s)*wroe(ix^s,iw)
301 integer,
intent(in) :: ix^L, idim
302 double precision,
intent(in) :: wL(ixG^T, nw), wR(ixG^T, nw)
303 double precision,
intent(inout) :: wroe(ixG^T, nw)
304 double precision,
intent(inout) :: workroe(ixG^T, nworkroe)
305 double precision,
intent(in) :: x(ixG^T, 1:^ND)
307 call average2_iso(wl,wr,x,ix^l,idim,wroe,workroe(ixg^t,1))
318 integer:: ix^L,idim,idir
319 double precision,
dimension(ixG^T,nw):: wL,wR,wroe
320 double precision,
intent(in) :: x(ixG^T,1:ndim)
321 double precision,
dimension(ixG^T):: tmp
329 wroe(ix^s,
mom(idir)) = half * (wl(ix^s,
mom(idir))/wl(ix^s,
rho_) + &
330 wr(ix^s,
mom(idir))/wr(ix^s,
rho_))
332 case (
'roe',
'default')
336 tmp(ix^s)=sqrt(wl(ix^s,
rho_)/wr(ix^s,
rho_))
338 wroe(ix^s,
mom(idir))=(wl(ix^s,
mom(idir))/wl(ix^s,
rho_)*tmp(ix^s)+&
339 wr(ix^s,
mom(idir))/wr(ix^s,
rho_))/(one+tmp(ix^s))
345 subroutine rhd_get_eigenjump_iso(wL,wR,wroe,x,ix^L,il,idim,smalla,a,jump,workroe)
354 integer,
intent(in) :: ix^L,il,idim
355 double precision,
dimension(ixG^T,nw):: wL,wR,wroe
356 double precision,
intent(in) :: x(ixG^T,1:ndim)
357 double precision,
dimension(ixG^T) :: smalla,a,jump
358 double precision,
dimension(ixG^T,nworkroe) :: workroe
360 call geteigenjump2_iso(wl,wr,wroe,x,ix^l,il,idim,smalla,a,jump,workroe(ixg^t,1))
364 subroutine geteigenjump2_iso(wL,wR,wroe,x,ix^L,il,idim,smalla,a,jump,csound)
374 integer:: ix^L,il,idim,idir
375 double precision,
dimension(ixG^T,nw):: wL,wR,wroe
376 double precision,
intent(in) :: x(ixG^T,1:ndim)
377 double precision,
dimension(ixG^T) :: smalla,a,jump,tmp,tmp2
378 double precision,
dimension(ixG^T) :: csound
379 DOUBLE PRECISION,
PARAMETER:: qsmall=1.
d-6
387 if (il == soundrw_)
then
388 a(ix^s)=wroe(ix^s,
mom(idim))+csound(ix^s)
389 jump(ix^s)=half*((one-wroe(ix^s,
mom(idim))/csound(ix^s))*&
391 +(wr(ix^s,
mom(idim))-wl(ix^s,
mom(idim)))/csound(ix^s))
392 else if (il == soundlw_)
then
393 a(ix^s)=wroe(ix^s,
mom(idim))-csound(ix^s)
394 jump(ix^s)=half*((one+wroe(ix^s,
mom(idim))/csound(ix^s))*&
396 -(wr(ix^s,
mom(idim))-wl(ix^s,
mom(idim)))/csound(ix^s))
399 idir=il-shearw0_;
if(idir>=idim)idir=idir+1
400 a(ix^s)=wroe(ix^s,
mom(idim))
401 jump(ix^s)=-wroe(ix^s,
mom(idir))*(wr(ix^s,
rho_)-wl(ix^s,
rho_))&
402 +(wr(ix^s,
mom(idir))-wl(ix^s,
mom(idir)))
404 case (
'roe',
'default')
408 where(abs(wl(ix^s,
rho_)-wr(ix^s,
rho_))<=qsmall*(wl(ix^s,
rho_)+wr(ix^s,
rho_)))
411 csound(ix^s) = sqrt(
rhd_gamma*(tmp2(ix^s)-tmp(ix^s))/&
416 if (il == soundrw_)
then
417 a(ix^s)=wroe(ix^s,
mom(idim))+csound(ix^s)
418 jump(ix^s)=half*((wr(ix^s,
rho_)-wl(ix^s,
rho_))+&
419 wroe(ix^s,
rho_)/csound(ix^s)*(wr(ix^s,
mom(idim))/wr(ix^s,
rho_)-&
420 wl(ix^s,
mom(idim))/wl(ix^s,
rho_)))
421 else if (il == soundlw_)
then
422 a(ix^s)=wroe(ix^s,
mom(idim))-csound(ix^s)
423 jump(ix^s)=half*((wr(ix^s,
rho_)-wl(ix^s,
rho_))-&
424 wroe(ix^s,
rho_)/csound(ix^s)*(wr(ix^s,
mom(idim))/wr(ix^s,
rho_)-&
425 wl(ix^s,
mom(idim))/wl(ix^s,
rho_)))
428 idir=il-shearw0_;
if(idir>=idim)idir=idir+1
429 a(ix^s)=wroe(ix^s,
mom(idim))
430 jump(ix^s)=wroe(ix^s,
rho_)*(wr(ix^s,
mom(idir))/wr(ix^s,
rho_)-&
431 wl(ix^s,
mom(idir))/wl(ix^s,
rho_))
441 case(
'harten',
'powell')
443 if (il == soundrw_)
then
445 tmp(ix^s) =wl(ix^s,
mom(idim))/wl(ix^s,
rho_)&
448 tmp2(ix^s)=wr(ix^s,
mom(idim))/wr(ix^s,
rho_)&
450 else if (il == soundlw_)
then
452 tmp(ix^s) =wl(ix^s,
mom(idim))/wl(ix^s,
rho_)&
455 tmp2(ix^s)=wr(ix^s,
mom(idim))/wr(ix^s,
rho_)&
458 tmp(ix^s) =wl(ix^s,
mom(idim))/wl(ix^s,
rho_)
459 tmp2(ix^s)=wr(ix^s,
mom(idim))/wr(ix^s,
rho_)
473 integer,
intent(in) :: ix^L,iw,il,idim
474 double precision,
intent(in) :: wroe(ixG^T,nw)
475 double precision,
intent(in) :: q(ixG^T)
476 double precision,
intent(inout) :: rq(ixG^T)
477 double precision,
intent(inout) :: workroe(ixG^T,nworkroe)
480 call rtimes2_iso(q,wroe,ix^l,iw,il,idim,rq,workroe(ixg^t,1))
490 integer:: ix^L,iw,il,idim,idir
491 double precision:: wroe(ixG^T,nw)
492 double precision,
dimension(ixG^T):: q,rq,csound
495 if (il == soundrw_ .or. il == soundlw_)
then
500 else if(iw==
mom(idim))
then
501 if (il == soundrw_)
then
502 rq(ix^s)=q(ix^s)*(wroe(ix^s,
mom(idim))+csound(ix^s))
503 else if (il == soundlw_)
then
504 rq(ix^s)=q(ix^s)*(wroe(ix^s,
mom(idim))-csound(ix^s))
509 if (il == soundrw_ .or. il == soundlw_)
then
510 rq(ix^s)=q(ix^s)*wroe(ix^s,iw)
513 idir=il-shearw0_;
if(idir>=idim)idir=idir+1
514 if(iw==
mom(idir))
then
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, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision, dimension(:), allocatable entropycoef
character(len=std_len) typeaverage
procedure(sub_rtimes), pointer phys_rtimes
procedure(sub_get_eigenjump), pointer phys_get_eigenjump
procedure(sub_average), pointer phys_average
Radiation-Hydrodynamics physics module Module aims at solving the Hydrodynamic equations toghether wi...
integer, public, protected rho_
Index of the density (in the w array)
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
logical, public, protected rhd_energy
Whether an energy equation is used.
integer, public, protected e_
Index of the energy density (-1 if not present)
subroutine, public rhd_get_pthermal(w, x, ixIL, ixOL, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
double precision, public rhd_gamma
The adiabatic index.
Module with Roe-type Riemann solver for hydrodynamics.
subroutine geteigenjump2(wL, wR, wroe, x, ixL, il, idim, smalla, a, jump, csound, dpperc2, dvperc)
subroutine average2_iso(wL, wR, x, ixL, idim, wroe, tmp)
subroutine rtimes2(q, wroe, ixL, iw, il, idim, rq, csound)
subroutine rhd_average(wL, wR, x, ixL, idim, wroe, workroe)
Calculate the Roe average of w, assignment of variables: rho -> rho, m -> v, e -> h.
subroutine average2(wL, wR, x, ixL, idim, wroe, tmp, tmp2)
subroutine, public rhd_roe_init()
subroutine geteigenjump2_iso(wL, wR, wroe, x, ixL, il, idim, smalla, a, jump, csound)
subroutine rhd_get_eigenjump_iso(wL, wR, wroe, x, ixL, il, idim, smalla, a, jump, workroe)
subroutine rhd_rtimes(q, wroe, ixL, iw, il, idim, rq, workroe)
Multiply q by R(il,iw), where R is the right eigenvalue matrix at wroe.
subroutine rhd_rtimes_iso(q, wroe, ixL, iw, il, idim, rq, workroe)
subroutine rhd_get_eigenjump(wL, wR, wroe, x, ixL, il, idim, smalla, a, jump, workroe)
Calculate the il-th characteristic speed and the jump in the il-th characteristic variable in the idi...
subroutine rhd_average_iso(wL, wR, x, ixL, idim, wroe, workroe)
subroutine rtimes2_iso(q, wroe, ixL, iw, il, idim, rq, csound)
Subroutines for TVD-MUSCL schemes.
subroutine, public entropyfix(ixL, il, aL, aR, a, smalla)