19 subroutine hancock(qdt,dtfactor,ixI^L,ixO^L,idims^LIM,qtC,sCT,qt,snew,dxs,x)
25 integer,
intent(in) :: ixi^
l, ixo^
l, idims^lim
26 double precision,
intent(in) :: qdt, dtfactor,qtc, qt, dxs(
ndim), x(ixi^s,1:
ndim)
27 type(state) :: sct, snew
29 double precision,
dimension(ixI^S,1:nw) :: wprim, wlc, wrc
31 double precision,
dimension(ixI^S,1:nw) :: wlp, wrp
32 double precision,
dimension(ixO^S) :: inv_volume
33 double precision :: flc(ixi^s, nwflux), frc(ixi^s, nwflux)
34 double precision :: dxinv(1:
ndim)
35 integer :: idims, iw, ix^
l, hxo^
l
38 associate(wct=>sct%w,wnew=>snew%w)
42 ix^
l=ix^
l^laddkr(idims,^
d);
44 if (ixi^
l^ltix^
l|.or.|.or.) &
45 call mpistop(
"Error in Hancock: Nonconforming input limits")
59 hxo^
l=ixo^
l-
kr(idims,^
d);
61 wrp(hxo^s,1:nw_recon)=wprim(ixo^s,1:nw_recon)
62 wlp(ixo^s,1:nw_recon)=wprim(ixo^s,1:nw_recon)
65 call reconstruct_lr(ixi^
l,ixo^
l,hxo^
l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
75 wnew(ixo^s,iw)=wnew(ixo^s,iw)-
block%dt(ixo^s)*dtfactor/dxs(idims)* &
76 (flc(ixo^s, iw)-frc(hxo^s, iw))
80 wnew(ixo^s,iw)=wnew(ixo^s,iw)+dxinv(idims)* &
81 (flc(ixo^s, iw)-frc(hxo^s, iw))
87 wnew(ixo^s,iw)=wnew(ixo^s,iw) -
block%dt(ixo^s)*dtfactor * inv_volume &
88 *(
block%surfaceC(ixo^s,idims)*flc(ixo^s, iw) &
89 -
block%surfaceC(hxo^s,idims)*frc(hxo^s, iw))
93 wnew(ixo^s,iw)=wnew(ixo^s,iw) - qdt * inv_volume &
94 *(
block%surfaceC(ixo^s,idims)*flc(ixo^s, iw) &
95 -
block%surfaceC(hxo^s,idims)*frc(hxo^s, iw))
105 dtfactor*dble(idimsmax-idimsmin+1)/dble(
ndim),&
106 ixi^
l,ixo^
l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
117 qtC,sCT,qt,snew,fC,fE,dxs,x)
126 integer,
intent(in) :: method
127 double precision,
intent(in) :: qdt, dtfactor, qtc, qt, dxs(
ndim)
128 integer,
intent(in) :: ixi^
l, ixo^
l, idims^lim
129 double precision,
dimension(ixI^S,1:ndim),
intent(in) :: x
130 double precision,
dimension(ixI^S,1:nwflux,1:ndim) :: fc
131 double precision,
dimension(ixI^S,sdim:3) :: fe
132 type(state) :: sct, snew
135 double precision,
dimension(ixI^S,1:nw) :: wprim
137 double precision,
dimension(ixI^S,1:nw) :: wlc, wrc
139 double precision,
dimension(ixI^S,1:nw) :: wlp, wrp
140 double precision,
dimension(ixI^S,1:nwflux) :: flc, frc
141 double precision,
dimension(ixI^S,1:number_species) :: cmaxc
142 double precision,
dimension(ixI^S,1:number_species) :: cminc
143 double precision,
dimension(ixI^S) :: hspeed
144 double precision,
dimension(ixO^S) :: inv_volume
145 double precision,
dimension(1:ndim) :: dxinv
146 integer :: idims, iw, ix^
d, hx^
d, ix^
l, hxo^
l, ixc^
l, ixcr^
l, kxc^
l, kxr^
l, ii
148 type(ct_velocity) :: vcts
150 associate(wct=>sct%w, wnew=>snew%w)
156 ix^
l=ix^
l^ladd2*
kr(idims,^
d);
158 if (ixi^
l^ltix^
l|.or.|.or.) &
159 call mpistop(
"Error in fv : Nonconforming input limits")
168 hxo^
l=ixo^
l-
kr(idims,^
d);
170 kxcmin^
d=iximin^
d; kxcmax^
d=iximax^
d-
kr(idims,^
d);
171 kxr^
l=kxc^
l+
kr(idims,^
d);
179 ixcmax^
d=ixomax^
d; ixcmin^
d=hxomin^
d;
186 wrp(kxc^s,1:
nw)=wprim(kxr^s,1:
nw)
187 wlp(kxc^s,1:
nw)=wprim(kxc^s,1:
nw)
193 call reconstruct_lr(ixi^
l,ixcr^
l,ixcr^
l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
211 call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^
l,ixc^
l,idims,hspeed,cmaxc,cminc)
240 call mpistop(
'unkown Riemann flux in finite volume')
251 hxomin^
d=ixomin^
d-hx^
d\
253 {
do ix^db=hxomin^db,ixomax^db\}
254 fc(ix^
d,iw,idims)=
block%dt(ix^
d)*dxinv(idims)*fc(ix^
d,iw,idims)
256 {
do ix^db=ixomin^db,ixomax^db\}
257 wnew(ix^d,iw)=wnew(ix^d,iw)+fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims)
261 if(method==fs_tvdmu) &
262 call tvdlimit2(method,qdt,ixi^l,ixc^l,ixo^l,idims,wlc,wrc,wnew,x,fc,dxs)
265 dxinv(1:ndim)=-qdt/dxs(1:ndim)
268 hxomin^d=ixomin^d-hx^d\
270 {
do ix^db=hxomin^db,ixomax^db\}
271 fc(ix^d,iw,idims)=dxinv(idims)*fc(ix^d,iw,idims)
273 {
do ix^db=ixomin^db,ixomax^db\}
274 wnew(ix^d,iw)=wnew(ix^d,iw)+fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims)
278 if(method==fs_tvdmu) &
279 call tvdlimit2(method,qdt,ixi^l,ixc^l,ixo^l,idims,wlc,wrc,wnew,x,fc,dxs)
283 inv_volume(ixo^s) = 1.d0/block%dvolume(ixo^s)
284 if(local_timestep)
then
287 hxomin^d=ixomin^d-hx^d\
289 {
do ix^db=hxomin^db,ixomax^db\}
290 fc(ix^d,iw,idims)=-block%dt(ix^d)*dtfactor*fc(ix^d,iw,idims)*block%surfaceC(ix^d,idims)
292 {
do ix^db=ixomin^db,ixomax^db\}
293 wnew(ix^d,iw)=wnew(ix^d,iw)+(fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims))*inv_volume(ix^d)
297 if (method==fs_tvdmu) &
298 call tvdlimit2(method,qdt,ixi^l,ixc^l,ixo^l,idims,wlc,wrc,wnew,x,fc,dxs)
303 hxomin^d=ixomin^d-hx^d\
305 {
do ix^db=hxomin^db,ixomax^db\}
306 fc(ix^d,iw,idims)=-qdt*fc(ix^d,iw,idims)*block%surfaceC(ix^d,idims)
308 {
do ix^db=ixomin^db,ixomax^db\}
309 wnew(ix^d,iw)=wnew(ix^d,iw)+(fc(ix^d,iw,idims)-fc(ix^d-hx^d,iw,idims))*inv_volume(ix^d)
316 if (.not.slab.and.idimsmin==1) &
317 call phys_add_source_geom(qdt,dtfactor,ixi^l,ixo^l,wct,wprim,wnew,x)
319 if(stagger_grid)
call phys_face_to_center(ixo^l,snew)
322 if(fix_small_values)
then
323 call phys_handle_small_values(.false.,wnew,x,ixi^l,ixo^l,
'multi-D finite_volume')
327 if(stagger_grid) snew%ws=sct%ws
331 call addsource2(qdt*dble(idimsmax-idimsmin+1)/dble(ndim),&
332 dtfactor*dble(idimsmax-idimsmin+1)/dble(ndim),&
333 ixi^l,ixo^l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
340 fc(ixc^s,iw,idims)=half*(flc(ixc^s,iw)+frc(ixc^s,iw))
344 subroutine get_riemann_flux_tvdlf(iws,iwe)
345 integer,
intent(in) :: iws,iwe
348 double precision :: fac(ixC^S),phi
350 fac(ixc^s) = -0.5d0*tvdlfeps*cmaxc(ixc^s,ii)
352 fc(ixc^s,iw,idims)=0.5d0*(flc(ixc^s, iw)+frc(ixc^s, iw))
354 if(flux_type(idims, iw) /= flux_no_dissipation)
then
355 if(flux_adaptive_diffusion)
then
356 {
do ix^db=ixcmin^db,ixcmax^db\}
357 jx^d=ix^d+kr(idims,^d)\
360 if(((wrc(ix^d,iw)-wlc(ix^d,iw))*(sct%w(jx^d,iw)-sct%w(ix^d,iw))) .gt. 1.e-18)
then
361 phi=min((wrc(ix^d,iw)-wlc(ix^d,iw))**2/((sct%w(jx^d,iw)-sct%w(ix^d,iw))**2+1.e-18),one)
365 fc(ix^d,iw,idims)=fc(ix^d,iw,idims)+fac(ix^d)*(wrc(ix^d,iw)-wlc(ix^d,iw))*phi
368 {
do ix^db=ixcmin^db,ixcmax^db\}
369 fc(ix^d,iw,idims)=fc(ix^d,iw,idims)+fac(ix^d)*(wrc(ix^d,iw)-wlc(ix^d,iw))
375 end subroutine get_riemann_flux_tvdlf
377 subroutine get_riemann_flux_hll(iws,iwe)
378 integer,
intent(in) :: iws,iwe
380 double precision :: phi
382 if(flux_adaptive_diffusion)
then
384 if(flux_type(idims, iw) == flux_tvdlf)
then
385 if(stagger_grid)
then
387 fc(ixc^s,iw,idims)=0.d0
389 fc(ixc^s,iw,idims)=-tvdlfeps*half*max(cmaxc(ixc^s,ii),dabs(cminc(ixc^s,ii)))*&
390 (wrc(ixc^s,iw)-wlc(ixc^s,iw))
393 {
do ix^db=ixcmin^db,ixcmax^db\}
394 if(cminc(ix^d,ii) >= zero)
then
395 fc(ix^d,iw,idims)=flc(ix^d,iw)
396 else if(cmaxc(ix^d,ii) <= zero)
then
397 fc(ix^d,iw,idims)=frc(ix^d,iw)
400 phi=max(abs(cmaxc(ix^d,ii)),abs(cminc(ix^d,ii)))/(cmaxc(ix^d,ii)-cminc(ix^d,ii))
401 fc(ix^d,iw,idims)=(cmaxc(ix^d,ii)*flc(ix^d, iw)-cminc(ix^d,ii)*frc(ix^d,iw)&
402 +phi*cminc(ix^d,ii)*cmaxc(ix^d,ii)*(wrc(ix^d,iw)-wlc(ix^d,iw)))&
403 /(cmaxc(ix^d,ii)-cminc(ix^d,ii))
410 if(flux_type(idims, iw) == flux_tvdlf)
then
411 if(stagger_grid)
then
413 fc(ixc^s,iw,idims)=0.d0
415 fc(ixc^s,iw,idims)=-tvdlfeps*half*max(cmaxc(ixc^s,ii),dabs(cminc(ixc^s,ii)))*&
416 (wrc(ixc^s,iw)-wlc(ixc^s,iw))
419 {
do ix^db=ixcmin^db,ixcmax^db\}
420 if(cminc(ix^d,ii) >= zero)
then
421 fc(ix^d,iw,idims)=flc(ix^d,iw)
422 else if(cmaxc(ix^d,ii) <= zero)
then
423 fc(ix^d,iw,idims)=frc(ix^d,iw)
425 fc(ix^d,iw,idims)=(cmaxc(ix^d,ii)*flc(ix^d, iw)-cminc(ix^d,ii)*frc(ix^d,iw)&
426 +cminc(ix^d,ii)*cmaxc(ix^d,ii)*(wrc(ix^d,iw)-wlc(ix^d,iw)))&
427 /(cmaxc(ix^d,ii)-cminc(ix^d,ii))
433 end subroutine get_riemann_flux_hll
435 subroutine get_riemann_flux_hllc(iws,iwe)
436 integer,
intent(in) :: iws, iwe
437 double precision,
dimension(ixI^S,1:nwflux) :: whll, Fhll, fCD
438 double precision,
dimension(ixI^S) :: lambdaCD
440 integer,
dimension(ixI^S) :: patchf
441 integer :: rho_, p_, e_, mom(1:ndir)
444 if (
allocated(iw_mom)) mom(:) = iw_mom(:)
447 if(
associated(phys_hllc_init_species))
then
448 call phys_hllc_init_species(ii, rho_, mom(:), e_)
454 where(cminc(ixc^s,1) >= zero)
456 elsewhere(cmaxc(ixc^s,1) <= zero)
460 if(method==fs_hllcd) &
461 call phys_diffuse_hllcd(ixi^l,ixc^l,idims,wlc,wrc,flc,frc,patchf)
464 if(any(patchf(ixc^s)==1)) &
465 call phys_get_lcd(wlc,wrc,flc,frc,cminc(ixi^s,ii),cmaxc(ixi^s,ii),idims,ixi^l,ixc^l, &
466 whll,fhll,lambdacd,patchf)
469 if(any(abs(patchf(ixc^s))== 1))
then
471 call phys_get_wcd(wlc,wrc,whll,frc,flc,fhll,patchf,lambdacd,&
472 cminc(ixi^s,ii),cmaxc(ixi^s,ii),ixi^l,ixc^l,idims,fcd)
476 if (flux_type(idims, iw) == flux_tvdlf)
then
477 flc(ixc^s,iw)=-tvdlfeps*half*max(cmaxc(ixc^s,ii),abs(cminc(ixc^s,ii))) * &
478 (wrc(ixc^s,iw) - wlc(ixc^s,iw))
480 where(patchf(ixc^s)==-2)
481 flc(ixc^s,iw)=flc(ixc^s,iw)
482 elsewhere(abs(patchf(ixc^s))==1)
483 flc(ixc^s,iw)=fcd(ixc^s,iw)
484 elsewhere(patchf(ixc^s)==2)
485 flc(ixc^s,iw)=frc(ixc^s,iw)
486 elsewhere(patchf(ixc^s)==3)
488 flc(ixc^s,iw)=fhll(ixc^s,iw)
489 elsewhere(patchf(ixc^s)==4)
491 flc(ixc^s,iw) = half*((flc(ixc^s,iw)+frc(ixc^s,iw)) &
492 -tvdlfeps * max(cmaxc(ixc^s,ii), dabs(cminc(ixc^s,ii))) * &
493 (wrc(ixc^s,iw)-wlc(ixc^s,iw)))
497 fc(ixc^s,iw,idims)=flc(ixc^s,iw)
500 end subroutine get_riemann_flux_hllc
503 subroutine get_riemann_flux_hlld(iws,iwe)
504 integer,
intent(in) :: iws, iwe
505 double precision,
dimension(ixI^S,1:nwflux) :: w1R,w1L,w2R,w2L
506 double precision,
dimension(ixI^S) :: sm,s1R,s1L,suR,suL,Bx
507 double precision,
dimension(ixI^S) :: pts,ptR,ptL,signBx,r1L,r1R,tmp
509 double precision,
dimension(ixI^S,ndir) :: BR, BL
510 integer :: ip1,ip2,ip3,idir,ix^D,^C&b^C_,^C&m^C_
511 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
513 associate(sr=>cmaxc,sl=>cminc)
516 ^c&mom(^c)=iw_mom(^c)\
518 ^c&mag(^c)=iw_mag(^c)\
526 br(ixc^s,:)=wrc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
527 bl(ixc^s,:)=wlc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
529 br(ixc^s,:)=wrc(ixc^s,mag(:))
530 bl(ixc^s,:)=wlc(ixc^s,mag(:))
532 if(stagger_grid)
then
533 bx(ixc^s)=block%ws(ixc^s,ip1)
537 bx(ixc^s)=(sr(ixc^s,ii)*br(ixc^s,ip1)-sl(ixc^s,ii)*bl(ixc^s,ip1))/(sr(ixc^s,ii)-sl(ixc^s,ii))
540 do ix^db=ixcmin^db,ixcmax^db\}
541 ptr(ix^d)=wrp(ix^d,p_)+0.5d0*(^c&br(ix^d,^c)**2+)
542 ptl(ix^d)=wlp(ix^d,p_)+0.5d0*(^c&bl(ix^d,^c)**2+)
543 if(iw_equi_rho>0)
then
544 sur(ix^d)=(sr(ix^d,ii)-wrp(ix^d,mom(ip1)))*(wrc(ix^d,rho_)+block%equi_vars(ix^d,iw_equi_rho,ip1))
545 sul(ix^d)=(sl(ix^d,ii)-wlp(ix^d,mom(ip1)))*(wlc(ix^d,rho_)+block%equi_vars(ix^d,iw_equi_rho,ip1))
547 sur(ix^d)=(sr(ix^d,ii)-wrp(ix^d,mom(ip1)))*wrc(ix^d,rho_)
548 sul(ix^d)=(sl(ix^d,ii)-wlp(ix^d,mom(ip1)))*wlc(ix^d,rho_)
551 sm(ix^d)=(sur(ix^d)*wrp(ix^d,mom(ip1))-sul(ix^d)*wlp(ix^d,mom(ip1))-&
552 ptr(ix^d)+ptl(ix^d))/(sur(ix^d)-sul(ix^d))
554 w1r(ix^d,mom(ip1))=sm(ix^d)
555 w1l(ix^d,mom(ip1))=sm(ix^d)
556 w2r(ix^d,mom(ip1))=sm(ix^d)
557 w2l(ix^d,mom(ip1))=sm(ix^d)
559 w1r(ix^d,mag(ip1))=bx(ix^d)
560 w1l(ix^d,mag(ip1))=bx(ix^d)
562 ptr(ix^d)=wrp(ix^d,p_)+0.5d0*(^c&wrc(ix^d,b^c_)**2+)
563 ptl(ix^d)=wlp(ix^d,p_)+0.5d0*(^c&wlc(ix^d,b^c_)**2+)
566 w1r(ix^d,rho_)=sur(ix^d)/(sr(ix^d,ii)-sm(ix^d))
567 w1l(ix^d,rho_)=sul(ix^d)/(sl(ix^d,ii)-sm(ix^d))
570 r1r(ix^d)=sur(ix^d)*(sr(ix^d,ii)-sm(ix^d))-bx(ix^d)**2
571 if(r1r(ix^d)/=0.d0) r1r(ix^d)=1.d0/r1r(ix^d)
572 r1l(ix^d)=sul(ix^d)*(sl(ix^d,ii)-sm(ix^d))-bx(ix^d)**2
573 if(r1l(ix^d)/=0.d0) r1l(ix^d)=1.d0/r1l(ix^d)
575 w1r(ix^d,mom(ip2))=wrp(ix^d,mom(ip2))-bx(ix^d)*br(ix^d,ip2)*&
576 (sm(ix^d)-wrp(ix^d,mom(ip1)))*r1r(ix^d)
577 w1l(ix^d,mom(ip2))=wlp(ix^d,mom(ip2))-bx(ix^d)*bl(ix^d,ip2)*&
578 (sm(ix^d)-wlp(ix^d,mom(ip1)))*r1l(ix^d)
580 w1r(ix^d,mag(ip2))=(sur(ix^d)*(sr(ix^d,ii)-wrp(ix^d,mom(ip1)))-bx(ix^d)**2)*r1r(ix^d)
581 w1l(ix^d,mag(ip2))=(sul(ix^d)*(sl(ix^d,ii)-wlp(ix^d,mom(ip1)))-bx(ix^d)**2)*r1l(ix^d)
586 w1r(ix^d,mom(ip3))=wrp(ix^d,mom(ip3))-bx(ix^d)*br(ix^d,ip3)*&
587 (sm(ix^d)-wrp(ix^d,mom(ip1)))*r1r(ix^d)
588 w1l(ix^d,mom(ip3))=wlp(ix^d,mom(ip3))-bx(ix^d)*bl(ix^d,ip3)*&
589 (sm(ix^d)-wlp(ix^d,mom(ip1)))*r1l(ix^d)
591 w1r(ix^d,mag(ip3))=br(ix^d,ip3)*w1r(ix^d,mag(ip2))
592 w1l(ix^d,mag(ip3))=bl(ix^d,ip3)*w1l(ix^d,mag(ip2))
595 w1r(ix^d,mag(ip2))=br(ix^d,ip2)*w1r(ix^d,mag(ip2))
596 w1l(ix^d,mag(ip2))=bl(ix^d,ip2)*w1l(ix^d,mag(ip2))
599 ^c&w1r(ix^d,b^c_)=w1r(ix^d,b^c_)-block%B0(ix^d,^c,ip1)\
600 ^c&w1l(ix^d,b^c_)=w1l(ix^d,b^c_)-block%B0(ix^d,^c,ip1)\
605 w1r(ix^d,p_)=sur(ix^d)*(sm(ix^d)-wrp(ix^d,mom(ip1)))+ptr(ix^d)
606 w1l(ix^d,p_)=w1r(ix^d,p_)
609 w1r(ix^d,p_)=w1r(ix^d,p_)+(^c&block%B0(ix^d,^c,ip1)*(wrc(ix^d,b^c_)-w1r(ix^d,b^c_))+)
610 w1l(ix^d,p_)=w1l(ix^d,p_)+(^c&block%B0(ix^d,^c,ip1)*(wlc(ix^d,b^c_)-w1l(ix^d,b^c_))+)
613 w1r(ix^d,e_)=((sr(ix^d,ii)-wrp(ix^d,mom(ip1)))*wrc(ix^d,e_)-ptr(ix^d)*wrp(ix^d,mom(ip1))+&
614 w1r(ix^d,p_)*sm(ix^d)+bx(ix^d)*((^c&wrp(ix^d,m^c_)*wrc(ix^d,b^c_)+)-&
615 (^c&w1r(ix^d,m^c_)*w1r(ix^d,b^c_)+)))/(sr(ix^d,ii)-sm(ix^d))
616 w1l(ix^d,e_)=((sl(ix^d,ii)-wlp(ix^d,mom(ip1)))*wlc(ix^d,e_)-ptl(ix^d)*wlp(ix^d,mom(ip1))+&
617 w1l(ix^d,p_)*sm(ix^d)+bx(ix^d)*((^c&wlp(ix^d,m^c_)*wlc(ix^d,b^c_)+)-&
618 (^c&w1l(ix^d,m^c_)*w1l(ix^d,b^c_)+)))/(sl(ix^d,ii)-sm(ix^d))
621 w1r(ix^d,e_)=w1r(ix^d,e_)+((^c&w1r(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*sm(ix^d)-&
622 (^c&wrc(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*wrp(ix^d,mom(ip1)))/(sr(ix^d,ii)-sm(ix^d))
623 w1l(ix^d,e_)=w1l(ix^d,e_)+((^c&w1l(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*sm(ix^d)-&
624 (^c&wlc(ix^d,b^c_)*block%B0(ix^d,^c,ip1)+)*wlp(ix^d,mom(ip1)))/(sl(ix^d,ii)-sm(ix^d))
627 w1r(ix^d,e_)=w1r(ix^d,e_)+1d0/(phys_gamma-1)*block%equi_vars(ix^d,iw_equi_p,ip1)*&
628 (sm(ix^d)-wrp(ix^d,mom(ip1)))/(sr(ix^d,ii)-sm(ix^d))
629 w1l(ix^d,e_)=w1l(ix^d,e_)+1d0/(phys_gamma-1)*block%equi_vars(ix^d,iw_equi_p,ip1)*&
630 (sm(ix^d)-wlp(ix^d,mom(ip1)))/(sl(ix^d,ii)-sm(ix^d))
635 w2r(ix^d,rho_)=w1r(ix^d,rho_)
636 w2l(ix^d,rho_)=w1l(ix^d,rho_)
637 w2r(ix^d,mag(ip1))=w1r(ix^d,mag(ip1))
638 w2l(ix^d,mag(ip1))=w1l(ix^d,mag(ip1))
639 r1r(ix^d)=sqrt(w1r(ix^d,rho_))
640 r1l(ix^d)=sqrt(w1l(ix^d,rho_))
641 tmp(ix^d)=1.d0/(r1r(ix^d)+r1l(ix^d))
642 signbx(ix^d)=sign(1.d0,bx(ix^d))
644 s1r(ix^d)=sm(ix^d)+abs(bx(ix^d))/r1r(ix^d)
645 s1l(ix^d)=sm(ix^d)-abs(bx(ix^d))/r1l(ix^d)
647 w2r(ix^d,mom(ip2))=(r1l(ix^d)*w1l(ix^d,mom(ip2))+r1r(ix^d)*w1r(ix^d,mom(ip2))+&
648 (w1r(ix^d,mag(ip2))-w1l(ix^d,mag(ip2)))*signbx(ix^d))*tmp(ix^d)
649 w2l(ix^d,mom(ip2))=w2r(ix^d,mom(ip2))
651 w2r(ix^d,mag(ip2))=(r1l(ix^d)*w1r(ix^d,mag(ip2))+r1r(ix^d)*w1l(ix^d,mag(ip2))+&
652 r1l(ix^d)*r1r(ix^d)*(w1r(ix^d,mom(ip2))-w1l(ix^d,mom(ip2)))*signbx(ix^d))*tmp(ix^d)
653 w2l(ix^d,mag(ip2))=w2r(ix^d,mag(ip2))
656 w2r(ix^d,mom(ip3))=(r1l(ix^d)*w1l(ix^d,mom(ip3))+r1r(ix^d)*w1r(ix^d,mom(ip3))+&
657 (w1r(ix^d,mag(ip3))-w1l(ix^d,mag(ip3)))*signbx(ix^d))*tmp(ix^d)
658 w2l(ix^d,mom(ip3))=w2r(ix^d,mom(ip3))
660 w2r(ix^d,mag(ip3))=(r1l(ix^d)*w1r(ix^d,mag(ip3))+r1r(ix^d)*w1l(ix^d,mag(ip3))+&
661 r1l(ix^d)*r1r(ix^d)*(w1r(ix^d,mom(ip3))-w1l(ix^d,mom(ip3)))*signbx(ix^d))*tmp(ix^d)
662 w2l(ix^d,mag(ip3))=w2r(ix^d,mag(ip3))
666 w2r(ix^d,e_)=w1r(ix^d,e_)+r1r(ix^d)*((^c&w1r(ix^d,m^c_)*w1r(ix^d,b^c_)+)-&
667 (^c&w2r(ix^d,m^c_)*w2r(ix^d,b^c_)+))*signbx(ix^d)
668 w2l(ix^d,e_)=w1l(ix^d,e_)-r1l(ix^d)*((^c&w1l(ix^d,m^c_)*w1l(ix^d,b^c_)+)-&
669 (^c&w2l(ix^d,m^c_)*w2l(ix^d,b^c_)+))*signbx(ix^d)
673 ^c&w1r(ix^d,m^c_)=w1r(ix^d,m^c_)*w1r(ix^d,rho_)\
674 ^c&w1l(ix^d,m^c_)=w1l(ix^d,m^c_)*w1l(ix^d,rho_)\
675 ^c&w2r(ix^d,m^c_)=w2r(ix^d,m^c_)*w2r(ix^d,rho_)\
676 ^c&w2l(ix^d,m^c_)=w2l(ix^d,m^c_)*w2l(ix^d,rho_)\
677 if(iw_equi_rho>0)
then
678 w1r(ix^d,rho_)=w1r(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
679 w1l(ix^d,rho_)=w1l(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
680 w2r(ix^d,rho_)=w2r(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
681 w2l(ix^d,rho_)=w2l(ix^d,rho_)-block%equi_vars(ix^d,iw_equi_rho,ip1)
686 if(flux_type(idims, iw)==flux_special)
then
689 do ix^db=ixcmin^db,ixcmax^db\}
690 fc(ix^d,iw,ip1)=flc(ix^d,iw)
692 else if(flux_type(idims, iw)==flux_hll)
then
695 do ix^db=ixcmin^db,ixcmax^db\}
696 fc(ix^d,iw,ip1)=(sr(ix^d,ii)*flc(ix^d,iw)-sl(ix^d,ii)*frc(ix^d,iw) &
697 +sr(ix^d,ii)*sl(ix^d,ii)*(wrc(ix^d,iw)-wlc(ix^d,iw)))/(sr(ix^d,ii)-sl(ix^d,ii))
704 do ix^db=ixcmin^db,ixcmax^db\}
705 if(sl(ix^d,ii)>0.d0)
then
706 fc(ix^d,iw,ip1)=flc(ix^d,iw)
707 else if(s1l(ix^d)>=0.d0)
then
708 fc(ix^d,iw,ip1)=flc(ix^d,iw)+sl(ix^d,ii)*(w1l(ix^d,iw)-wlc(ix^d,iw))
709 else if(sm(ix^d)>=0.d0)
then
710 fc(ix^d,iw,ip1)=flc(ix^d,iw)+sl(ix^d,ii)*(w1l(ix^d,iw)-wlc(ix^d,iw))+&
711 s1l(ix^d)*(w2l(ix^d,iw)-w1l(ix^d,iw))
712 else if(s1r(ix^d)>=0.d0)
then
713 fc(ix^d,iw,ip1)=frc(ix^d,iw)+sr(ix^d,ii)*(w1r(ix^d,iw)-wrc(ix^d,iw))+&
714 s1r(ix^d)*(w2r(ix^d,iw)-w1r(ix^d,iw))
715 else if(sr(ix^d,ii)>=0.d0)
then
716 fc(ix^d,iw,ip1)=frc(ix^d,iw)+sr(ix^d,ii)*(w1r(ix^d,iw)-wrc(ix^d,iw))
717 else if(sr(ix^d,ii)<0.d0)
then
718 fc(ix^d,iw,ip1)=frc(ix^d,iw)
725 end subroutine get_riemann_flux_hlld
729 subroutine get_riemann_flux_hlld_mag2(iws,iwe)
731 integer,
intent(in) :: iws, iwe
733 double precision,
dimension(ixI^S,1:nwflux) :: w1R,w1L,f1R,f1L,f2R,f2L
734 double precision,
dimension(ixI^S,1:nwflux) :: w2R,w2L
735 double precision,
dimension(ixI^S) :: sm,s1R,s1L,suR,suL,Bx
736 double precision,
dimension(ixI^S) :: pts,ptR,ptL,signBx,r1L,r1R,tmp
738 double precision,
dimension(ixI^S,ndir) :: vRC, vLC
740 double precision,
dimension(ixI^S,ndir) :: BR, BL
741 integer :: ip1,ip2,ip3,idir,ix^D
742 double precision :: phiPres, thetaSM, du, dv, dw
743 integer :: ixV^L, ixVb^L, ixVc^L, ixVd^L, ixVe^L, ixVf^L
744 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
745 double precision,
parameter :: aParam = 4d0
753 associate(sr=>cmaxc,sl=>cminc)
760 vrc(ixc^s,:)=wrp(ixc^s,mom(:))
761 vlc(ixc^s,:)=wlp(ixc^s,mom(:))
764 call get_hlld2_modif_c(wlp,x,ixi^l,ixo^l,s1l)
765 call get_hlld2_modif_c(wrp,x,ixi^l,ixo^l,s1r)
767 phipres = min(1d0, maxval(max(s1l(ixo^s),s1r(ixo^s)))/maxval(cmaxc(ixo^s,1)))
768 phipres = phipres*(2d0 - phipres)
775 du = minval(wprim(ixv^s,mom(1))-wprim(ixo^s,mom(1)))
811 dv = minval(min(wprim(ixo^s,mom(2))-wprim(ixv^s,mom(2)),&
812 wprim(ixvb^s,mom(2))-wprim(ixo^s,mom(2)),&
813 wprim(ixvc^s,mom(2))-wprim(ixvd^s,mom(2)),&
814 wprim(ixve^s,mom(2))-wprim(ixvc^s,mom(2))&
848 dw = minval(min(wprim(ixo^s,mom(3))-wprim(ixv^s,mom(3)),&
849 wprim(ixvb^s,mom(3))-wprim(ixo^s,mom(3)),&
850 wprim(ixvc^s,mom(3))-wprim(ixvd^s,mom(3)),&
851 wprim(ixve^s,mom(3))-wprim(ixvc^s,mom(3))&
854 thetasm = maxval(cmaxc(ixo^s,1))
856 thetasm = (min(1d0, (thetasm-du)/(thetasm-min(dv,dw))))**aparam
860 br(ixc^s,:)=wrc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
861 bl(ixc^s,:)=wlc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
863 br(ixc^s,:)=wrc(ixc^s,mag(:))
864 bl(ixc^s,:)=wlc(ixc^s,mag(:))
867 bx(ixc^s)=(sr(ixc^s,index_v_mag)*br(ixc^s,ip1)-sl(ixc^s,index_v_mag)*bl(ixc^s,ip1)-&
868 flc(ixc^s,mag(ip1))-frc(ixc^s,mag(ip1)))/(sr(ixc^s,index_v_mag)-sl(ixc^s,index_v_mag))
869 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(br(ixc^s,:)**2,dim=ndim+1)
870 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(bl(ixc^s,:)**2,dim=ndim+1)
871 sur(ixc^s)=(sr(ixc^s,index_v_mag)-vrc(ixc^s,ip1))*wrc(ixc^s,rho_)
872 sul(ixc^s)=(sl(ixc^s,index_v_mag)-vlc(ixc^s,ip1))*wlc(ixc^s,rho_)
874 sm(ixc^s)=(sur(ixc^s)*vrc(ixc^s,ip1)-sul(ixc^s)*vlc(ixc^s,ip1)-&
875 thetasm*(ptr(ixc^s)-ptl(ixc^s)) )/(sur(ixc^s)-sul(ixc^s))
877 w1r(ixc^s,mom(ip1))=sm(ixc^s)
878 w1l(ixc^s,mom(ip1))=sm(ixc^s)
879 w2r(ixc^s,mom(ip1))=sm(ixc^s)
880 w2l(ixc^s,mom(ip1))=sm(ixc^s)
882 w1r(ixc^s,mag(ip1))=bx(ixc^s)
883 w1l(ixc^s,mag(ip1))=bx(ixc^s)
885 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(wrc(ixc^s,mag(:))**2,dim=ndim+1)
886 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(wlc(ixc^s,mag(:))**2,dim=ndim+1)
890 w1r(ixc^s,rho_)=sur(ixc^s)/(sr(ixc^s,index_v_mag)-sm(ixc^s))
891 w1l(ixc^s,rho_)=sul(ixc^s)/(sl(ixc^s,index_v_mag)-sm(ixc^s))
895 r1r(ixc^s)=sur(ixc^s)*(sr(ixc^s,index_v_mag)-sm(ixc^s))-bx(ixc^s)**2
896 where(r1r(ixc^s)/=0.d0)
897 r1r(ixc^s)=1.d0/r1r(ixc^s)
899 r1l(ixc^s)=sul(ixc^s)*(sl(ixc^s,index_v_mag)-sm(ixc^s))-bx(ixc^s)**2
900 where(r1l(ixc^s)/=0.d0)
901 r1l(ixc^s)=1.d0/r1l(ixc^s)
904 w1r(ixc^s,mom(ip2))=vrc(ixc^s,ip2)-bx(ixc^s)*br(ixc^s,ip2)*&
905 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
906 w1l(ixc^s,mom(ip2))=vlc(ixc^s,ip2)-bx(ixc^s)*bl(ixc^s,ip2)*&
907 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
909 w1r(ixc^s,mag(ip2))=(sur(ixc^s)*(sr(ixc^s,index_v_mag)-vrc(ixc^s,ip1))-bx(ixc^s)**2)*r1r(ixc^s)
910 w1l(ixc^s,mag(ip2))=(sul(ixc^s)*(sl(ixc^s,index_v_mag)-vlc(ixc^s,ip1))-bx(ixc^s)**2)*r1l(ixc^s)
915 w1r(ixc^s,mom(ip3))=vrc(ixc^s,ip3)-bx(ixc^s)*br(ixc^s,ip3)*&
916 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
917 w1l(ixc^s,mom(ip3))=vlc(ixc^s,ip3)-bx(ixc^s)*bl(ixc^s,ip3)*&
918 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
920 w1r(ixc^s,mag(ip3))=br(ixc^s,ip3)*w1r(ixc^s,mag(ip2))
921 w1l(ixc^s,mag(ip3))=bl(ixc^s,ip3)*w1l(ixc^s,mag(ip2))
924 w1r(ixc^s,mag(ip2))=br(ixc^s,ip2)*w1r(ixc^s,mag(ip2))
925 w1l(ixc^s,mag(ip2))=bl(ixc^s,ip2)*w1l(ixc^s,mag(ip2))
928 w1r(ixc^s,mag(:))=w1r(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
929 w1l(ixc^s,mag(:))=w1l(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
934 w1r(ixc^s,p_)=(sur(ixc^s)*ptl(ixc^s) - sul(ixc^s)*ptr(ixc^s) +&
935 phipres * sur(ixc^s)*sul(ixc^s)*(vrc(ixc^s,ip1)-vlc(ixc^s,ip1)))/&
936 (sur(ixc^s)-sul(ixc^s))
937 w1l(ixc^s,p_)=w1r(ixc^s,p_)
940 w1r(ixc^s,p_)=w1r(ixc^s,p_)+sum(block%B0(ixc^s,:,ip1)*(wrc(ixc^s,mag(:))-w1r(ixc^s,mag(:))),dim=ndim+1)
941 w1l(ixc^s,p_)=w1l(ixc^s,p_)+sum(block%B0(ixc^s,:,ip1)*(wlc(ixc^s,mag(:))-w1l(ixc^s,mag(:))),dim=ndim+1)
944 w1r(ixc^s,e_)=((sr(ixc^s,index_v_mag)-vrc(ixc^s,ip1))*wrc(ixc^s,e_)-ptr(ixc^s)*vrc(ixc^s,ip1)+&
945 w1r(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vrc(ixc^s,:)*wrc(ixc^s,mag(:)),dim=ndim+1)-&
946 sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)))/(sr(ixc^s,index_v_mag)-sm(ixc^s))
947 w1l(ixc^s,e_)=((sl(ixc^s,index_v_mag)-vlc(ixc^s,ip1))*wlc(ixc^s,e_)-ptl(ixc^s)*vlc(ixc^s,ip1)+&
948 w1l(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vlc(ixc^s,:)*wlc(ixc^s,mag(:)),dim=ndim+1)-&
949 sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)))/(sl(ixc^s,index_v_mag)-sm(ixc^s))
952 w1r(ixc^s,e_)=w1r(ixc^s,e_)+(sum(w1r(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
953 sum(wrc(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*vrc(ixc^s,ip1))/(sr(ixc^s,index_v_mag)-sm(ixc^s))
954 w1l(ixc^s,e_)=w1l(ixc^s,e_)+(sum(w1l(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
955 sum(wlc(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*vlc(ixc^s,ip1))/(sl(ixc^s,index_v_mag)-sm(ixc^s))
960 w2r(ixc^s,rho_)=w1r(ixc^s,rho_)
961 w2l(ixc^s,rho_)=w1l(ixc^s,rho_)
962 w2r(ixc^s,mag(ip1))=w1r(ixc^s,mag(ip1))
963 w2l(ixc^s,mag(ip1))=w1l(ixc^s,mag(ip1))
965 r1r(ixc^s)=sqrt(w1r(ixc^s,rho_))
966 r1l(ixc^s)=sqrt(w1l(ixc^s,rho_))
967 tmp(ixc^s)=1.d0/(r1r(ixc^s)+r1l(ixc^s))
968 signbx(ixc^s)=sign(1.d0,bx(ixc^s))
970 s1r(ixc^s)=sm(ixc^s)+abs(bx(ixc^s))/r1r(ixc^s)
971 s1l(ixc^s)=sm(ixc^s)-abs(bx(ixc^s))/r1l(ixc^s)
973 w2r(ixc^s,mom(ip2))=(r1l(ixc^s)*w1l(ixc^s,mom(ip2))+r1r(ixc^s)*w1r(ixc^s,mom(ip2))+&
974 (w1r(ixc^s,mag(ip2))-w1l(ixc^s,mag(ip2)))*signbx(ixc^s))*tmp(ixc^s)
975 w2l(ixc^s,mom(ip2))=w2r(ixc^s,mom(ip2))
977 w2r(ixc^s,mag(ip2))=(r1l(ixc^s)*w1r(ixc^s,mag(ip2))+r1r(ixc^s)*w1l(ixc^s,mag(ip2))+&
978 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip2))-w1l(ixc^s,mom(ip2)))*signbx(ixc^s))*tmp(ixc^s)
979 w2l(ixc^s,mag(ip2))=w2r(ixc^s,mag(ip2))
982 w2r(ixc^s,mom(ip3))=(r1l(ixc^s)*w1l(ixc^s,mom(ip3))+r1r(ixc^s)*w1r(ixc^s,mom(ip3))+&
983 (w1r(ixc^s,mag(ip3))-w1l(ixc^s,mag(ip3)))*signbx(ixc^s))*tmp(ixc^s)
984 w2l(ixc^s,mom(ip3))=w2r(ixc^s,mom(ip3))
986 w2r(ixc^s,mag(ip3))=(r1l(ixc^s)*w1r(ixc^s,mag(ip3))+r1r(ixc^s)*w1l(ixc^s,mag(ip3))+&
987 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip3))-w1l(ixc^s,mom(ip3)))*signbx(ixc^s))*tmp(ixc^s)
988 w2l(ixc^s,mag(ip3))=w2r(ixc^s,mag(ip3))
992 w2r(ixc^s,e_)=w1r(ixc^s,e_)+r1r(ixc^s)*(sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)-&
993 sum(w2r(ixc^s,mom(:))*w2r(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
994 w2l(ixc^s,e_)=w1l(ixc^s,e_)-r1l(ixc^s)*(sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)-&
995 sum(w2l(ixc^s,mom(:))*w2l(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
1000 w1r(ixc^s,mom(idir))=w1r(ixc^s,mom(idir))*w1r(ixc^s,rho_)
1001 w1l(ixc^s,mom(idir))=w1l(ixc^s,mom(idir))*w1l(ixc^s,rho_)
1002 w2r(ixc^s,mom(idir))=w2r(ixc^s,mom(idir))*w2r(ixc^s,rho_)
1003 w2l(ixc^s,mom(idir))=w2l(ixc^s,mom(idir))*w2l(ixc^s,rho_)
1009 if(stagger_grid .and. flux_type(idims, iw) == flux_tvdlf) cycle
1010 if(flux_type(idims, iw) == flux_special)
then
1012 f1l(ixc^s,iw)=flc(ixc^s,iw)
1013 f1r(ixc^s,iw)=f1l(ixc^s,iw)
1014 f2l(ixc^s,iw)=f1l(ixc^s,iw)
1015 f2r(ixc^s,iw)=f1l(ixc^s,iw)
1016 else if(flux_type(idims, iw) == flux_hll)
then
1018 f1l(ixc^s,iw)=(sr(ixc^s,index_v_mag)*flc(ixc^s, iw)-sl(ixc^s,index_v_mag)*frc(ixc^s, iw) &
1019 +sr(ixc^s,index_v_mag)*sl(ixc^s,index_v_mag)*(wrc(ixc^s,iw)-wlc(ixc^s,iw)))/(sr(ixc^s,index_v_mag)-sl(ixc^s,index_v_mag))
1020 f1r(ixc^s,iw)=f1l(ixc^s,iw)
1021 f2l(ixc^s,iw)=f1l(ixc^s,iw)
1022 f2r(ixc^s,iw)=f1l(ixc^s,iw)
1024 f1l(ixc^s,iw)=flc(ixc^s,iw)+sl(ixc^s,index_v_mag)*(w1l(ixc^s,iw)-wlc(ixc^s,iw))
1025 f1r(ixc^s,iw)=frc(ixc^s,iw)+sr(ixc^s,index_v_mag)*(w1r(ixc^s,iw)-wrc(ixc^s,iw))
1026 f2l(ixc^s,iw)=f1l(ixc^s,iw)+s1l(ixc^s)*(w2l(ixc^s,iw)-w1l(ixc^s,iw))
1027 f2r(ixc^s,iw)=f1r(ixc^s,iw)+s1r(ixc^s)*(w2r(ixc^s,iw)-w1r(ixc^s,iw))
1032 {
do ix^db=ixcmin^db,ixcmax^db\}
1033 if(sl(ix^d,index_v_mag)>0.d0)
then
1034 fc(ix^d,iws:iwe,ip1)=flc(ix^d,iws:iwe)
1035 else if(s1l(ix^d)>=0.d0)
then
1036 fc(ix^d,iws:iwe,ip1)=f1l(ix^d,iws:iwe)
1037 else if(sm(ix^d)>=0.d0)
then
1038 fc(ix^d,iws:iwe,ip1)=f2l(ix^d,iws:iwe)
1039 else if(s1r(ix^d)>=0.d0)
then
1040 fc(ix^d,iws:iwe,ip1)=f2r(ix^d,iws:iwe)
1041 else if(sr(ix^d,index_v_mag)>=0.d0)
then
1042 fc(ix^d,iws:iwe,ip1)=f1r(ix^d,iws:iwe)
1043 else if(sr(ix^d,index_v_mag)<0.d0)
then
1044 fc(ix^d,iws:iwe,ip1)=frc(ix^d,iws:iwe)
1049 end subroutine get_riemann_flux_hlld_mag2
1052 subroutine get_hlld2_modif_c(w,x,ixI^L,ixO^L,csound)
1055 integer,
intent(in) :: ixI^L, ixO^L
1056 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
1057 double precision,
intent(out):: csound(ixI^S)
1058 double precision :: cfast2(ixI^S), AvMinCs2(ixI^S), b2(ixI^S), kmax
1059 double precision :: inv_rho(ixO^S), gamma_A2(ixO^S)
1060 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
1068 inv_rho=1.d0/w(ixo^s,rho_)
1073 b2(ixo^s) = sum((w(ixo^s, mag(:))+
block%B0(ixo^s,:,
b0i))**2, dim=ndim+1)
1075 b2(ixo^s) = sum(w(ixo^s, mag(:))**2, dim=ndim+1)
1080 avmincs2= w(ixo^s, mag(idims))+
block%B0(ixo^s,idims,
b0i)
1082 avmincs2= w(ixo^s, mag(idims))
1086 csound(ixo^s) = sum(w(ixo^s, mom(:))**2, dim=ndim+1)
1088 cfast2(ixo^s) = b2(ixo^s) * inv_rho+csound(ixo^s)
1089 avmincs2(ixo^s) = cfast2(ixo^s)**2-4.0d0*csound(ixo^s) &
1090 * avmincs2(ixo^s)**2 &
1093 where(avmincs2(ixo^s)<zero)
1094 avmincs2(ixo^s)=zero
1097 avmincs2(ixo^s)=sqrt(avmincs2(ixo^s))
1099 csound(ixo^s) = sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s)))
1101 end subroutine get_hlld2_modif_c
1107 subroutine reconstruct_lr(ixI^L,ixL^L,ixR^L,idims,w,wLC,wRC,wLp,wRp,x,dxdim)
1113 integer,
intent(in) :: ixi^
l, ixl^
l, ixr^
l, idims
1114 double precision,
intent(in) :: dxdim
1116 double precision,
dimension(ixI^S,1:nw) :: w
1118 double precision,
dimension(ixI^S,1:nw) :: wlc, wrc
1120 double precision,
dimension(ixI^S,1:nw) :: wlp, wrp
1121 double precision,
dimension(ixI^S,1:ndim) :: x
1123 integer :: jxr^
l, ixc^
l, jxc^
l, iw
1124 double precision :: ldw(ixi^s), rdw(ixi^s), dwc(ixi^s)
1125 double precision :: a2max
1129 call mp5limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp)
1131 call weno3limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,1)
1133 call weno3limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,2)
1135 call weno5limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,1)
1137 call weno5nmlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,1)
1139 call weno5limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,2)
1141 call weno5nmlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,2)
1143 call weno5limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,3)
1145 call weno5nmlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,3)
1147 call weno5cu6limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp)
1149 call teno5adlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp)
1151 call weno7limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp,1)
1153 call weno7limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp,2)
1155 call venklimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp)
1161 call ppmlimiter(ixi^
l,
ixm^
ll,idims,w,w,wlp,wrp)
1167 jxr^
l=ixr^
l+
kr(idims,^
d);
1168 ixcmax^
d=jxrmax^
d; ixcmin^
d=ixlmin^
d-
kr(idims,^
d);
1169 jxc^
l=ixc^
l+
kr(idims,^
d);
1172 w(ixcmin^
d:jxcmax^
d,iw)=dlog10(w(ixcmin^
d:jxcmax^
d,iw))
1173 wlp(ixl^s,iw)=dlog10(wlp(ixl^s,iw))
1174 wrp(ixr^s,iw)=dlog10(wrp(ixr^s,iw))
1177 dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
1193 call mpistop(
"idims is wrong in mod_limiter")
1199 wlp(ixl^s,iw)=wlp(ixl^s,iw)+half*ldw(ixl^s)
1200 wrp(ixr^s,iw)=wrp(ixr^s,iw)-half*rdw(jxr^s)
1203 w(ixcmin^
d:jxcmax^
d,iw)=10.0d0**w(ixcmin^
d:jxcmax^
d,iw)
1204 wlp(ixl^s,iw)=10.0d0**wlp(ixl^s,iw)
1205 wrp(ixr^s,iw)=10.0d0**wrp(ixr^s,iw)
1214 wlc(ixl^s,1:nw_recon)=wlp(ixl^s,1:nw_recon)
1215 wrc(ixr^s,1:nw_recon)=wrp(ixr^s,1:nw_recon)
1219 wlp(ixl^s,nw_recon+1:nw_recon+nwaux)=wlc(ixl^s,nw_recon+1:nw_recon+nwaux)
1220 wrp(ixr^s,nw_recon+1:nw_recon+nwaux)=wrc(ixr^s,nw_recon+1:nw_recon+nwaux)
subroutine get_riemann_flux_tvdmu()
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module with finite volume methods for fluxes.
subroutine, public finite_volume(method, qdt, dtfactor, ixIL, ixOL, idimsLIM, qtC, sCT, qt, snew, fC, fE, dxs, x)
finite volume method
subroutine, public reconstruct_lr(ixIL, ixLL, ixRL, idims, w, wLC, wRC, wLp, wRp, x, dxdim)
Determine the upwinded wLC(ixL) and wRC(ixR) from w. the wCT is only used when PPM is exploited.
subroutine, public hancock(qdt, dtfactor, ixIL, ixOL, idimsLIM, qtC, sCT, qt, snew, dxs, x)
The non-conservative Hancock predictor for TVDLF.
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
logical h_correction
If true, do H-correction to fix the carbuncle problem at grid-aligned shocks.
integer, parameter fs_tvdlf
integer ixghi
Upper index of grid block arrays.
integer, parameter fs_hlld
integer, dimension(3, 3) kr
Kronecker delta tensor.
logical, dimension(:), allocatable loglimit
integer, parameter fs_hllcd
integer, parameter ndim
Number of spatial dimensions for grid variables.
logical stagger_grid
True for using stagger grid.
integer, parameter fs_tvdmu
integer b0i
background magnetic field location indicator
double precision, dimension(:), allocatable, parameter d
logical local_timestep
each cell has its own timestep or not
logical need_global_a2max
global value for schmid scheme
integer ixm
the mesh range of a physical block without ghost cells
logical slab
Cartesian geometry or not.
integer, parameter fs_hll
flux schemes
logical b0field
split magnetic field as background B0 field
integer, dimension(:), allocatable type_limiter
Type of slope limiter used for reconstructing variables on cell edges.
integer nghostcells
Number of ghost cells surrounding a grid.
logical fix_small_values
fix small values with average or replace methods
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
double precision, dimension(^nd) a2max_global
global largest a2 for schmid scheme
integer, parameter fs_hllc
integer, parameter ixglo
Lower index of grid block arrays (always 1)
Module with slope/flux limiters.
integer, parameter limiter_weno5cu6
integer, parameter limiter_mpweno7
integer, parameter limiter_teno5ad
integer, parameter limiter_weno3
integer, parameter limiter_ppm
subroutine dwlimiter2(dwC, ixIL, ixCL, idims, typelim, ldw, rdw, a2max)
Limit the centered dwC differences within ixC for iw in direction idim. The limiter is chosen accordi...
integer, parameter limiter_wenozp5
integer, parameter limiter_weno5
integer, parameter limiter_wenoz5
integer, parameter limiter_wenoz5nm
integer, parameter limiter_weno7
integer, parameter limiter_wenozp5nm
integer, parameter limiter_mp5
integer, parameter limiter_venk
integer, parameter limiter_weno5nm
integer, parameter limiter_wenoyc3
This module defines the procedures of a physics module. It contains function pointers for the various...
procedure(sub_get_ct_velocity), pointer phys_get_ct_velocity
procedure(sub_convert), pointer phys_to_primitive
procedure(sub_small_values), pointer phys_handle_small_values
procedure(sub_get_flux), pointer phys_get_flux
procedure(sub_get_cbounds), pointer phys_get_cbounds
procedure(sub_add_source_geom), pointer phys_add_source_geom
integer phys_wider_stencil
To use wider stencils in flux calculations. A value of 1 will extend it by one cell in both direction...
procedure(sub_update_faces), pointer phys_update_faces
procedure(sub_convert), pointer phys_to_conserved
procedure(sub_modify_wlr), pointer phys_modify_wlr
procedure(sub_get_h_speed), pointer phys_get_h_speed
Module for handling split source terms (split from the fluxes)
subroutine, public addsource2(qdt, dtfactor, ixIL, ixOL, iwLIM, qtC, wCT, wCTprim, qt, w, x, qsourcesplit, src_active)
Add source within ixO for iws: w=w+qdt*S[wCT].
Subroutines for TVD-MUSCL schemes.
subroutine, public tvdlimit2(method, qdt, ixIL, ixICL, ixOL, idims, wL, wR, wnew, x, fC, dxs)
Module with all the methods that users can customize in AMRVAC.
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
integer, dimension(:), allocatable start_indices
the indices in 1:nwflux array are assumed consecutive for each species this array should be of size n...
integer, dimension(:), allocatable stop_indices
the indices in 1:nwflux array are assumed consecutive for each species this array should be of size n...
integer nwflux
Number of flux variables.
integer index_v_mag
index of the var whose velocity appears in the induction eq.