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")
57 hxo^
l=ixo^
l-
kr(idims,^
d);
59 wrp(hxo^s,1:nwflux)=wprim(ixo^s,1:nwflux)
60 wlp(ixo^s,1:nwflux)=wprim(ixo^s,1:nwflux)
63 call reconstruct_lr(ixi^
l,ixo^
l,hxo^
l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
73 wnew(ixo^s,iw)=wnew(ixo^s,iw)-
block%dt(ixo^s)*dtfactor/dxs(idims)* &
74 (flc(ixo^s, iw)-frc(hxo^s, iw))
78 wnew(ixo^s,iw)=wnew(ixo^s,iw)+dxinv(idims)* &
79 (flc(ixo^s, iw)-frc(hxo^s, iw))
85 wnew(ixo^s,iw)=wnew(ixo^s,iw) -
block%dt(ixo^s)*dtfactor * inv_volume &
86 *(
block%surfaceC(ixo^s,idims)*flc(ixo^s, iw) &
87 -
block%surfaceC(hxo^s,idims)*frc(hxo^s, iw))
91 wnew(ixo^s,iw)=wnew(ixo^s,iw) - qdt * inv_volume &
92 *(
block%surfaceC(ixo^s,idims)*flc(ixo^s, iw) &
93 -
block%surfaceC(hxo^s,idims)*frc(hxo^s, iw))
103 dtfactor*dble(idimsmax-idimsmin+1)/dble(
ndim),&
104 ixi^
l,ixo^
l,1,nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
115 qtC,sCT,qt,snew,fC,fE,dxs,x)
124 integer,
intent(in) :: method
125 double precision,
intent(in) :: qdt, dtfactor, qtc, qt, dxs(
ndim)
126 integer,
intent(in) :: ixi^
l, ixo^
l, idims^lim
127 double precision,
dimension(ixI^S,1:ndim),
intent(in) :: x
128 type(state) :: sct, snew
129 double precision,
dimension(ixI^S,1:nwflux,1:ndim) :: fc
130 double precision,
dimension(ixI^S,sdim:3) :: fe
133 double precision,
dimension(ixI^S,1:nw) :: wprim
135 double precision,
dimension(ixI^S,1:nw) :: wlc, wrc
137 double precision,
dimension(ixI^S,1:nw) :: wlp, wrp
138 double precision,
dimension(ixI^S,1:nwflux) :: flc, frc
139 double precision,
dimension(ixI^S,1:number_species) :: cmaxc
140 double precision,
dimension(ixI^S,1:number_species) :: cminc
141 double precision,
dimension(ixI^S) :: hspeed
142 double precision,
dimension(ixO^S) :: inv_volume
143 double precision,
dimension(1:ndim) :: dxinv
144 integer,
dimension(ixI^S) :: patchf
145 integer :: idims, iw, ix^
l, hxo^
l, ixc^
l, ixcr^
l, kxc^
l, kxr^
l, ii
147 type(ct_velocity) :: vcts
149 associate(wct=>sct%w, wnew=>snew%w)
159 ix^
l=ix^
l^ladd2*
kr(idims,^
d);
161 if (ixi^
l^ltix^
l|.or.|.or.) &
162 call mpistop(
"Error in fv : Nonconforming input limits")
171 hxo^
l=ixo^
l-
kr(idims,^
d);
173 kxcmin^
d=iximin^
d; kxcmax^
d=iximax^
d-
kr(idims,^
d);
174 kxr^
l=kxc^
l+
kr(idims,^
d);
182 ixcmax^
d=ixomax^
d; ixcmin^
d=hxomin^
d;
189 wrp(kxc^s,1:
nw)=wprim(kxr^s,1:
nw)
190 wlp(kxc^s,1:
nw)=wprim(kxc^s,1:
nw)
197 call reconstruct_lr(ixi^
l,ixcr^
l,ixcr^
l,idims,wprim,wlc,wrc,wlp,wrp,x,dxs(idims))
215 call phys_get_cbounds(wlc,wrc,wlp,wrp,x,ixi^
l,ixc^
l,idims,hspeed,cmaxc,cminc)
244 call mpistop(
'unkown Riemann flux in finite volume')
253 hxo^
l=ixo^
l-
kr(idims,^
d);
258 fc(ixi^s,iw,idims)=-
block%dt(ixi^s)*dtfactor/dxs(idims)*fc(ixi^s,iw,idims)
262 fc(ixi^s,1:
nwflux,idims)=dxinv(idims)*fc(ixi^s,1:
nwflux,idims)
270 call tvdlimit2(method,qdt,ixi^
l,ixc^
l,ixo^
l,idims,wlc,wrc,wnew,x,fc,dxs)
274 inv_volume = 1.d0/
block%dvolume(ixo^s)
276 hxo^
l=ixo^
l-
kr(idims,^
d);
280 fc(ixi^s,iw,idims)=-
block%dt(ixi^s)*dtfactor*fc(ixi^s,iw,idims)*
block%surfaceC(ixi^s,idims)
281 wnew(ixo^s,iw)=wnew(ixo^s,iw) + (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims)) * &
286 fc(ixi^s,iw,idims)=-qdt*fc(ixi^s,iw,idims)*
block%surfaceC(ixi^s,idims)
287 wnew(ixo^s,iw)=wnew(ixo^s,iw) + (fc(ixo^s,iw,idims)-fc(hxo^s,iw,idims)) * &
293 call tvdlimit2(method,qdt,ixi^
l,ixc^
l,ixo^
l,idims,wlc,wrc,wnew,x,fc,dxs)
298 if (.not.
slab.and.idimsmin==1) &
309 dtfactor*dble(idimsmax-idimsmin+1)/dble(
ndim),&
310 ixi^
l,ixo^
l,1,
nw,qtc,wct,wprim,qt,wnew,x,.false.,active)
318 flc(ixc^s, iw)=half*(flc(ixc^s, iw)+frc(ixc^s, iw))
319 fc(ixc^s,iw,idims)=flc(ixc^s, iw)
323 subroutine get_riemann_flux_tvdlf(iws,iwe)
324 integer,
intent(in) :: iws,iwe
325 double precision :: fac(ixC^S)
327 fac(ixc^s) = -0.5d0*tvdlfeps*cmaxc(ixc^s,ii)
331 flc(ixc^s, iw)=0.5d0*(flc(ixc^s, iw)+frc(ixc^s, iw))
333 if (flux_type(idims, iw) /= flux_no_dissipation)
then
334 flc(ixc^s, iw)=flc(ixc^s, iw) + fac(ixc^s)*(wrc(ixc^s,iw)-wlc(ixc^s,iw))
336 fc(ixc^s,iw,idims)=flc(ixc^s, iw)
339 end subroutine get_riemann_flux_tvdlf
341 subroutine get_riemann_flux_hll(iws,iwe)
342 integer,
intent(in) :: iws,iwe
346 if(flux_type(idims, iw) == flux_tvdlf)
then
348 if(stagger_grid) cycle
349 fc(ixc^s,iw,idims) = -tvdlfeps*half*max(cmaxc(ixc^s,ii),dabs(cminc(ixc^s,ii))) * &
350 (wrc(ixc^s,iw)-wlc(ixc^s,iw))
352 {
do ix^db=ixcmin^db,ixcmax^db\}
353 if(cminc(ix^d,ii) >= zero)
then
354 fc(ix^d,iw,idims)=flc(ix^d,iw)
355 else if(cmaxc(ix^d,ii) <= zero)
then
356 fc(ix^d,iw,idims)=frc(ix^d,iw)
359 fc(ix^d,iw,idims)=(cmaxc(ix^d,ii)*flc(ix^d, iw)-cminc(ix^d,ii)*frc(ix^d,iw)&
360 +cminc(ix^d,ii)*cmaxc(ix^d,ii)*(wrc(ix^d,iw)-wlc(ix^d,iw)))&
361 /(cmaxc(ix^d,ii)-cminc(ix^d,ii))
367 end subroutine get_riemann_flux_hll
369 subroutine get_riemann_flux_hllc(iws,iwe)
370 integer,
intent(in) :: iws, iwe
371 double precision,
dimension(ixI^S,1:nwflux) :: whll, Fhll, fCD
372 double precision,
dimension(ixI^S) :: lambdaCD
374 integer :: rho_, p_, e_, mom(1:ndir)
377 if (
allocated(iw_mom)) mom(:) = iw_mom(:)
380 if(
associated(phys_hllc_init_species))
then
381 call phys_hllc_init_species(ii, rho_, mom(:), e_)
387 where(cminc(ixc^s,1) >= zero)
389 elsewhere(cmaxc(ixc^s,1) <= zero)
393 if(method==fs_hllcd) &
394 call phys_diffuse_hllcd(ixi^l,ixc^l,idims,wlc,wrc,flc,frc,patchf)
397 if(any(patchf(ixc^s)==1)) &
398 call phys_get_lcd(wlc,wrc,flc,frc,cminc(ixi^s,ii),cmaxc(ixi^s,ii),idims,ixi^l,ixc^l, &
399 whll,fhll,lambdacd,patchf)
402 if(any(abs(patchf(ixc^s))== 1))
then
404 call phys_get_wcd(wlc,wrc,whll,frc,flc,fhll,patchf,lambdacd,&
405 cminc(ixi^s,ii),cmaxc(ixi^s,ii),ixi^l,ixc^l,idims,fcd)
409 if (flux_type(idims, iw) == flux_tvdlf)
then
410 flc(ixc^s,iw)=-tvdlfeps*half*max(cmaxc(ixc^s,ii),abs(cminc(ixc^s,ii))) * &
411 (wrc(ixc^s,iw) - wlc(ixc^s,iw))
413 where(patchf(ixc^s)==-2)
414 flc(ixc^s,iw)=flc(ixc^s,iw)
415 elsewhere(abs(patchf(ixc^s))==1)
416 flc(ixc^s,iw)=fcd(ixc^s,iw)
417 elsewhere(patchf(ixc^s)==2)
418 flc(ixc^s,iw)=frc(ixc^s,iw)
419 elsewhere(patchf(ixc^s)==3)
421 flc(ixc^s,iw)=fhll(ixc^s,iw)
422 elsewhere(patchf(ixc^s)==4)
424 flc(ixc^s,iw) = half*((flc(ixc^s,iw)+frc(ixc^s,iw)) &
425 -tvdlfeps * max(cmaxc(ixc^s,ii), dabs(cminc(ixc^s,ii))) * &
426 (wrc(ixc^s,iw)-wlc(ixc^s,iw)))
430 fc(ixc^s,iw,idims)=flc(ixc^s,iw)
433 end subroutine get_riemann_flux_hllc
436 subroutine get_riemann_flux_hlld(iws,iwe)
437 integer,
intent(in) :: iws, iwe
438 double precision,
dimension(ixI^S,1:nwflux) :: w1R,w1L,f1R,f1L,f2R,f2L
439 double precision,
dimension(ixI^S,1:nwflux) :: w2R,w2L
440 double precision,
dimension(ixI^S) :: sm,s1R,s1L,suR,suL,Bx
441 double precision,
dimension(ixI^S) :: pts,ptR,ptL,signBx,r1L,r1R,tmp
443 double precision,
dimension(ixI^S,ndir) :: vRC, vLC
445 double precision,
dimension(ixI^S,ndir) :: BR, BL
446 integer :: ip1,ip2,ip3,idir,ix^D
447 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
449 associate(sr=>cmaxc,sl=>cminc)
467 vrc(ixc^s,:)=wrp(ixc^s,mom(:))
468 vlc(ixc^s,:)=wlp(ixc^s,mom(:))
470 br(ixc^s,:)=wrc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
471 bl(ixc^s,:)=wlc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
473 br(ixc^s,:)=wrc(ixc^s,mag(:))
474 bl(ixc^s,:)=wlc(ixc^s,mag(:))
476 if(stagger_grid)
then
477 bx(ixc^s)=block%ws(ixc^s,ip1)
481 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))
483 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(br(ixc^s,:)**2,dim=ndim+1)
484 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(bl(ixc^s,:)**2,dim=ndim+1)
485 if(iw_equi_rho>0)
then
486 sur(ixc^s) = wrc(ixc^s,rho_)+ block%equi_vars(ixc^s,iw_equi_rho,ip1)
488 sur(ixc^s) = wrc(ixc^s,rho_)
490 sur(ixc^s)=(sr(ixc^s,ii)-vrc(ixc^s,ip1))*sur(ixc^s)
491 if(iw_equi_rho>0)
then
492 sul(ixc^s) = wlc(ixc^s,rho_)+ block%equi_vars(ixc^s,iw_equi_rho,ip1)
494 sul(ixc^s) = wlc(ixc^s,rho_)
496 sul(ixc^s)=(sl(ixc^s,ii)-vlc(ixc^s,ip1))*sul(ixc^s)
498 sm(ixc^s)=(sur(ixc^s)*vrc(ixc^s,ip1)-sul(ixc^s)*vlc(ixc^s,ip1)-&
499 ptr(ixc^s)+ptl(ixc^s))/(sur(ixc^s)-sul(ixc^s))
501 w1r(ixc^s,mom(ip1))=sm(ixc^s)
502 w1l(ixc^s,mom(ip1))=sm(ixc^s)
503 w2r(ixc^s,mom(ip1))=sm(ixc^s)
504 w2l(ixc^s,mom(ip1))=sm(ixc^s)
506 w1r(ixc^s,mag(ip1))=bx(ixc^s)
507 w1l(ixc^s,mag(ip1))=bx(ixc^s)
509 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(wrc(ixc^s,mag(:))**2,dim=ndim+1)
510 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(wlc(ixc^s,mag(:))**2,dim=ndim+1)
514 w1r(ixc^s,rho_)=sur(ixc^s)/(sr(ixc^s,ii)-sm(ixc^s))
515 w1l(ixc^s,rho_)=sul(ixc^s)/(sl(ixc^s,ii)-sm(ixc^s))
519 r1r(ixc^s)=sur(ixc^s)*(sr(ixc^s,ii)-sm(ixc^s))-bx(ixc^s)**2
520 where(abs(r1r(ixc^s))>smalldouble)
521 r1r(ixc^s)=1.d0/r1r(ixc^s)
525 r1l(ixc^s)=sul(ixc^s)*(sl(ixc^s,ii)-sm(ixc^s))-bx(ixc^s)**2
526 where(abs(r1l(ixc^s))>smalldouble)
527 r1l(ixc^s)=1.d0/r1l(ixc^s)
532 w1r(ixc^s,mom(ip2))=vrc(ixc^s,ip2)-bx(ixc^s)*br(ixc^s,ip2)*&
533 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
534 w1l(ixc^s,mom(ip2))=vlc(ixc^s,ip2)-bx(ixc^s)*bl(ixc^s,ip2)*&
535 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
537 w1r(ixc^s,mag(ip2))=(sur(ixc^s)*(sr(ixc^s,ii)-vrc(ixc^s,ip1))-bx(ixc^s)**2)*r1r(ixc^s)
538 w1l(ixc^s,mag(ip2))=(sul(ixc^s)*(sl(ixc^s,ii)-vlc(ixc^s,ip1))-bx(ixc^s)**2)*r1l(ixc^s)
543 w1r(ixc^s,mom(ip3))=vrc(ixc^s,ip3)-bx(ixc^s)*br(ixc^s,ip3)*&
544 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
545 w1l(ixc^s,mom(ip3))=vlc(ixc^s,ip3)-bx(ixc^s)*bl(ixc^s,ip3)*&
546 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
548 w1r(ixc^s,mag(ip3))=br(ixc^s,ip3)*w1r(ixc^s,mag(ip2))
549 w1l(ixc^s,mag(ip3))=bl(ixc^s,ip3)*w1l(ixc^s,mag(ip2))
552 w1r(ixc^s,mag(ip2))=br(ixc^s,ip2)*w1r(ixc^s,mag(ip2))
553 w1l(ixc^s,mag(ip2))=bl(ixc^s,ip2)*w1l(ixc^s,mag(ip2))
556 w1r(ixc^s,mag(:))=w1r(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
557 w1l(ixc^s,mag(:))=w1l(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
562 w1r(ixc^s,p_)=sur(ixc^s)*(sm(ixc^s)-vrc(ixc^s,ip1))+ptr(ixc^s)
564 w1l(ixc^s,p_)=w1r(ixc^s,p_)
567 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)
568 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)
571 w1r(ixc^s,e_)=((sr(ixc^s,ii)-vrc(ixc^s,ip1))*wrc(ixc^s,e_)-ptr(ixc^s)*vrc(ixc^s,ip1)+&
572 w1r(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vrc(ixc^s,:)*wrc(ixc^s,mag(:)),dim=ndim+1)-&
573 sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)))/(sr(ixc^s,ii)-sm(ixc^s))
574 w1l(ixc^s,e_)=((sl(ixc^s,ii)-vlc(ixc^s,ip1))*wlc(ixc^s,e_)-ptl(ixc^s)*vlc(ixc^s,ip1)+&
575 w1l(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vlc(ixc^s,:)*wlc(ixc^s,mag(:)),dim=ndim+1)-&
576 sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)))/(sl(ixc^s,ii)-sm(ixc^s))
579 w1r(ixc^s,e_)=w1r(ixc^s,e_)+(sum(w1r(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
580 sum(wrc(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*vrc(ixc^s,ip1))/(sr(ixc^s,ii)-sm(ixc^s))
581 w1l(ixc^s,e_)=w1l(ixc^s,e_)+(sum(w1l(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
582 sum(wlc(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*vlc(ixc^s,ip1))/(sl(ixc^s,ii)-sm(ixc^s))
585 #if !defined(E_RM_W0) || E_RM_W0 == 1
586 w1r(ixc^s,e_)= w1r(ixc^s,e_) + 1d0/(phys_gamma - 1) * block%equi_vars(ixc^s,iw_equi_p,ip1) * &
587 (sm(ixc^s)-vrc(ixc^s,ip1))/(sr(ixc^s,ii)-sm(ixc^s))
588 w1l(ixc^s,e_)= w1l(ixc^s,e_) + 1d0/(phys_gamma - 1) * block%equi_vars(ixc^s,iw_equi_p,ip1) * &
589 (sm(ixc^s)-vlc(ixc^s,ip1))/(sl(ixc^s,ii)-sm(ixc^s))
591 w1r(ixc^s,e_)= w1r(ixc^s,e_) + phys_gamma /(phys_gamma - 1) * block%equi_vars(ixc^s,iw_equi_p,ip1) * &
592 (sm(ixc^s)-vrc(ixc^s,ip1))/(sr(ixc^s,ii)-sm(ixc^s))
593 w1l(ixc^s,e_)= w1l(ixc^s,e_) + phys_gamma /(phys_gamma - 1) * block%equi_vars(ixc^s,iw_equi_p,ip1) * &
594 (sm(ixc^s)-vlc(ixc^s,ip1))/(sl(ixc^s,ii)-sm(ixc^s))
600 w2r(ixc^s,rho_)=w1r(ixc^s,rho_)
601 w2l(ixc^s,rho_)=w1l(ixc^s,rho_)
602 w2r(ixc^s,mag(ip1))=w1r(ixc^s,mag(ip1))
603 w2l(ixc^s,mag(ip1))=w1l(ixc^s,mag(ip1))
605 r1r(ixc^s)=sqrt(w1r(ixc^s,rho_))
606 r1l(ixc^s)=sqrt(w1l(ixc^s,rho_))
607 tmp(ixc^s)=1.d0/(r1r(ixc^s)+r1l(ixc^s))
608 signbx(ixc^s)=sign(1.d0,bx(ixc^s))
610 s1r(ixc^s)=sm(ixc^s)+abs(bx(ixc^s))/r1r(ixc^s)
611 s1l(ixc^s)=sm(ixc^s)-abs(bx(ixc^s))/r1l(ixc^s)
613 w2r(ixc^s,mom(ip2))=(r1l(ixc^s)*w1l(ixc^s,mom(ip2))+r1r(ixc^s)*w1r(ixc^s,mom(ip2))+&
614 (w1r(ixc^s,mag(ip2))-w1l(ixc^s,mag(ip2)))*signbx(ixc^s))*tmp(ixc^s)
615 w2l(ixc^s,mom(ip2))=w2r(ixc^s,mom(ip2))
617 w2r(ixc^s,mag(ip2))=(r1l(ixc^s)*w1r(ixc^s,mag(ip2))+r1r(ixc^s)*w1l(ixc^s,mag(ip2))+&
618 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip2))-w1l(ixc^s,mom(ip2)))*signbx(ixc^s))*tmp(ixc^s)
619 w2l(ixc^s,mag(ip2))=w2r(ixc^s,mag(ip2))
622 w2r(ixc^s,mom(ip3))=(r1l(ixc^s)*w1l(ixc^s,mom(ip3))+r1r(ixc^s)*w1r(ixc^s,mom(ip3))+&
623 (w1r(ixc^s,mag(ip3))-w1l(ixc^s,mag(ip3)))*signbx(ixc^s))*tmp(ixc^s)
624 w2l(ixc^s,mom(ip3))=w2r(ixc^s,mom(ip3))
626 w2r(ixc^s,mag(ip3))=(r1l(ixc^s)*w1r(ixc^s,mag(ip3))+r1r(ixc^s)*w1l(ixc^s,mag(ip3))+&
627 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip3))-w1l(ixc^s,mom(ip3)))*signbx(ixc^s))*tmp(ixc^s)
628 w2l(ixc^s,mag(ip3))=w2r(ixc^s,mag(ip3))
632 w2r(ixc^s,e_)=w1r(ixc^s,e_)+r1r(ixc^s)*(sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)-&
633 sum(w2r(ixc^s,mom(:))*w2r(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
634 w2l(ixc^s,e_)=w1l(ixc^s,e_)-r1l(ixc^s)*(sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)-&
635 sum(w2l(ixc^s,mom(:))*w2l(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
640 w1r(ixc^s,mom(idir))=w1r(ixc^s,mom(idir))*w1r(ixc^s,rho_)
641 w1l(ixc^s,mom(idir))=w1l(ixc^s,mom(idir))*w1l(ixc^s,rho_)
642 w2r(ixc^s,mom(idir))=w2r(ixc^s,mom(idir))*w2r(ixc^s,rho_)
643 w2l(ixc^s,mom(idir))=w2l(ixc^s,mom(idir))*w2l(ixc^s,rho_)
645 if(iw_equi_rho>0)
then
646 w1r(ixc^s,rho_) = w1r(ixc^s,rho_) - block%equi_vars(ixc^s,iw_equi_rho,ip1)
647 w1l(ixc^s,rho_) = w1l(ixc^s,rho_) - block%equi_vars(ixc^s,iw_equi_rho,ip1)
648 w2r(ixc^s,rho_) = w2r(ixc^s,rho_) - block%equi_vars(ixc^s,iw_equi_rho,ip1)
649 w2l(ixc^s,rho_) = w2l(ixc^s,rho_) - block%equi_vars(ixc^s,iw_equi_rho,ip1)
653 if(flux_type(idims, iw) == flux_special)
then
655 f1l(ixc^s,iw)=flc(ixc^s,iw)
656 f1r(ixc^s,iw)=f1l(ixc^s,iw)
657 f2l(ixc^s,iw)=f1l(ixc^s,iw)
658 f2r(ixc^s,iw)=f1l(ixc^s,iw)
659 else if(flux_type(idims, iw) == flux_hll)
then
661 f1l(ixc^s,iw)=(sr(ixc^s,ii)*flc(ixc^s, iw)-sl(ixc^s,ii)*frc(ixc^s, iw) &
662 +sr(ixc^s,ii)*sl(ixc^s,ii)*(wrc(ixc^s,iw)-wlc(ixc^s,iw)))/(sr(ixc^s,ii)-sl(ixc^s,ii))
663 f1r(ixc^s,iw)=f1l(ixc^s,iw)
664 f2l(ixc^s,iw)=f1l(ixc^s,iw)
665 f2r(ixc^s,iw)=f1l(ixc^s,iw)
668 f1l(ixc^s,iw)=flc(ixc^s,iw)+sl(ixc^s,ii)*(w1l(ixc^s,iw)-wlc(ixc^s,iw))
669 f1r(ixc^s,iw)=frc(ixc^s,iw)+sr(ixc^s,ii)*(w1r(ixc^s,iw)-wrc(ixc^s,iw))
670 f2l(ixc^s,iw)=f1l(ixc^s,iw)+s1l(ixc^s)*(w2l(ixc^s,iw)-w1l(ixc^s,iw))
671 f2r(ixc^s,iw)=f1r(ixc^s,iw)+s1r(ixc^s)*(w2r(ixc^s,iw)-w1r(ixc^s,iw))
676 {
do ix^db=ixcmin^db,ixcmax^db\}
677 if(sl(ix^d,ii)>0.d0)
then
678 fc(ix^d,iws:iwe,ip1)=flc(ix^d,iws:iwe)
679 else if(s1l(ix^d)>=0.d0)
then
680 fc(ix^d,iws:iwe,ip1)=f1l(ix^d,iws:iwe)
681 else if(sm(ix^d)>=0.d0)
then
682 fc(ix^d,iws:iwe,ip1)=f2l(ix^d,iws:iwe)
683 else if(s1r(ix^d)>=0.d0)
then
684 fc(ix^d,iws:iwe,ip1)=f2r(ix^d,iws:iwe)
685 else if(sr(ix^d,ii)>=0.d0)
then
686 fc(ix^d,iws:iwe,ip1)=f1r(ix^d,iws:iwe)
687 else if(sr(ix^d,ii)<0.d0)
then
688 fc(ix^d,iws:iwe,ip1)=frc(ix^d,iws:iwe)
693 end subroutine get_riemann_flux_hlld
697 subroutine get_riemann_flux_hlld_mag2()
702 double precision,
dimension(ixI^S,1:nwflux) :: w1R,w1L,f1R,f1L,f2R,f2L
703 double precision,
dimension(ixI^S,1:nwflux) :: w2R,w2L
704 double precision,
dimension(ixI^S) :: sm,s1R,s1L,suR,suL,Bx
705 double precision,
dimension(ixI^S) :: pts,ptR,ptL,signBx,r1L,r1R,tmp
707 double precision,
dimension(ixI^S,ndir) :: vRC, vLC
709 double precision,
dimension(ixI^S,ndir) :: BR, BL
710 integer :: ip1,ip2,ip3,idir,ix^D
711 double precision :: phiPres, thetaSM, du, dv, dw
712 integer :: ixV^L, ixVb^L, ixVc^L, ixVd^L, ixVe^L, ixVf^L
713 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
714 double precision,
parameter :: aParam = 4d0
722 associate(sr=>cmaxc,sl=>cminc)
729 vrc(ixc^s,:)=wrp(ixc^s,mom(:))
730 vlc(ixc^s,:)=wlp(ixc^s,mom(:))
736 phipres = min(1d0, maxval(max(s1l(ixo^s),s1r(ixo^s)))/maxval(cmaxc(ixo^s,1)))
737 phipres = phipres*(2d0 - phipres)
744 du = minval(wprim(ixv^s,mom(1))-wprim(ixo^s,mom(1)))
780 dv = minval(min(wprim(ixo^s,mom(2))-wprim(ixv^s,mom(2)),&
781 wprim(ixvb^s,mom(2))-wprim(ixo^s,mom(2)),&
782 wprim(ixvc^s,mom(2))-wprim(ixvd^s,mom(2)),&
783 wprim(ixve^s,mom(2))-wprim(ixvc^s,mom(2))&
817 dw = minval(min(wprim(ixo^s,mom(3))-wprim(ixv^s,mom(3)),&
818 wprim(ixvb^s,mom(3))-wprim(ixo^s,mom(3)),&
819 wprim(ixvc^s,mom(3))-wprim(ixvd^s,mom(3)),&
820 wprim(ixve^s,mom(3))-wprim(ixvc^s,mom(3))&
823 thetasm = maxval(cmaxc(ixo^s,1))
825 thetasm = (min(1d0, (thetasm-du)/(thetasm-min(dv,dw))))**aparam
829 br(ixc^s,:)=wrc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
830 bl(ixc^s,:)=wlc(ixc^s,mag(:))+block%B0(ixc^s,:,ip1)
832 br(ixc^s,:)=wrc(ixc^s,mag(:))
833 bl(ixc^s,:)=wlc(ixc^s,mag(:))
838 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(br(ixc^s,:)**2,dim=ndim+1)
839 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(bl(ixc^s,:)**2,dim=ndim+1)
840 sur(ixc^s)=(sr(ixc^s,
index_v_mag)-vrc(ixc^s,ip1))*wrc(ixc^s,rho_)
841 sul(ixc^s)=(sl(ixc^s,
index_v_mag)-vlc(ixc^s,ip1))*wlc(ixc^s,rho_)
843 sm(ixc^s)=(sur(ixc^s)*vrc(ixc^s,ip1)-sul(ixc^s)*vlc(ixc^s,ip1)-&
844 thetasm*(ptr(ixc^s)-ptl(ixc^s)) )/(sur(ixc^s)-sul(ixc^s))
846 w1r(ixc^s,mom(ip1))=sm(ixc^s)
847 w1l(ixc^s,mom(ip1))=sm(ixc^s)
848 w2r(ixc^s,mom(ip1))=sm(ixc^s)
849 w2l(ixc^s,mom(ip1))=sm(ixc^s)
851 w1r(ixc^s,mag(ip1))=bx(ixc^s)
852 w1l(ixc^s,mag(ip1))=bx(ixc^s)
854 ptr(ixc^s)=wrp(ixc^s,p_)+0.5d0*sum(wrc(ixc^s,mag(:))**2,dim=ndim+1)
855 ptl(ixc^s)=wlp(ixc^s,p_)+0.5d0*sum(wlc(ixc^s,mag(:))**2,dim=ndim+1)
859 w1r(ixc^s,rho_)=sur(ixc^s)/(sr(ixc^s,
index_v_mag)-sm(ixc^s))
860 w1l(ixc^s,rho_)=sul(ixc^s)/(sl(ixc^s,
index_v_mag)-sm(ixc^s))
864 r1r(ixc^s)=sur(ixc^s)*(sr(ixc^s,
index_v_mag)-sm(ixc^s))-bx(ixc^s)**2
865 where(r1r(ixc^s)/=0.d0)
866 r1r(ixc^s)=1.d0/r1r(ixc^s)
868 r1l(ixc^s)=sul(ixc^s)*(sl(ixc^s,
index_v_mag)-sm(ixc^s))-bx(ixc^s)**2
869 where(r1l(ixc^s)/=0.d0)
870 r1l(ixc^s)=1.d0/r1l(ixc^s)
873 w1r(ixc^s,mom(ip2))=vrc(ixc^s,ip2)-bx(ixc^s)*br(ixc^s,ip2)*&
874 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
875 w1l(ixc^s,mom(ip2))=vlc(ixc^s,ip2)-bx(ixc^s)*bl(ixc^s,ip2)*&
876 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
878 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)
879 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)
884 w1r(ixc^s,mom(ip3))=vrc(ixc^s,ip3)-bx(ixc^s)*br(ixc^s,ip3)*&
885 (sm(ixc^s)-vrc(ixc^s,ip1))*r1r(ixc^s)
886 w1l(ixc^s,mom(ip3))=vlc(ixc^s,ip3)-bx(ixc^s)*bl(ixc^s,ip3)*&
887 (sm(ixc^s)-vlc(ixc^s,ip1))*r1l(ixc^s)
889 w1r(ixc^s,mag(ip3))=br(ixc^s,ip3)*w1r(ixc^s,mag(ip2))
890 w1l(ixc^s,mag(ip3))=bl(ixc^s,ip3)*w1l(ixc^s,mag(ip2))
893 w1r(ixc^s,mag(ip2))=br(ixc^s,ip2)*w1r(ixc^s,mag(ip2))
894 w1l(ixc^s,mag(ip2))=bl(ixc^s,ip2)*w1l(ixc^s,mag(ip2))
897 w1r(ixc^s,mag(:))=w1r(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
898 w1l(ixc^s,mag(:))=w1l(ixc^s,mag(:))-block%B0(ixc^s,:,ip1)
903 w1r(ixc^s,p_)=(sur(ixc^s)*ptl(ixc^s) - sul(ixc^s)*ptr(ixc^s) +&
904 phipres * sur(ixc^s)*sul(ixc^s)*(vrc(ixc^s,ip1)-vlc(ixc^s,ip1)))/&
905 (sur(ixc^s)-sul(ixc^s))
906 w1l(ixc^s,p_)=w1r(ixc^s,p_)
909 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)
910 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)
913 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)+&
914 w1r(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vrc(ixc^s,:)*wrc(ixc^s,mag(:)),dim=ndim+1)-&
915 sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)))/(sr(ixc^s,
index_v_mag)-sm(ixc^s))
916 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)+&
917 w1l(ixc^s,p_)*sm(ixc^s)+bx(ixc^s)*(sum(vlc(ixc^s,:)*wlc(ixc^s,mag(:)),dim=ndim+1)-&
918 sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)))/(sl(ixc^s,
index_v_mag)-sm(ixc^s))
921 w1r(ixc^s,e_)=w1r(ixc^s,e_)+(sum(w1r(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
922 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))
923 w1l(ixc^s,e_)=w1l(ixc^s,e_)+(sum(w1l(ixc^s,mag(:))*block%B0(ixc^s,:,ip1),dim=ndim+1)*sm(ixc^s)-&
924 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))
929 w2r(ixc^s,rho_)=w1r(ixc^s,rho_)
930 w2l(ixc^s,rho_)=w1l(ixc^s,rho_)
931 w2r(ixc^s,mag(ip1))=w1r(ixc^s,mag(ip1))
932 w2l(ixc^s,mag(ip1))=w1l(ixc^s,mag(ip1))
934 r1r(ixc^s)=sqrt(w1r(ixc^s,rho_))
935 r1l(ixc^s)=sqrt(w1l(ixc^s,rho_))
936 tmp(ixc^s)=1.d0/(r1r(ixc^s)+r1l(ixc^s))
937 signbx(ixc^s)=sign(1.d0,bx(ixc^s))
939 s1r(ixc^s)=sm(ixc^s)+abs(bx(ixc^s))/r1r(ixc^s)
940 s1l(ixc^s)=sm(ixc^s)-abs(bx(ixc^s))/r1l(ixc^s)
942 w2r(ixc^s,mom(ip2))=(r1l(ixc^s)*w1l(ixc^s,mom(ip2))+r1r(ixc^s)*w1r(ixc^s,mom(ip2))+&
943 (w1r(ixc^s,mag(ip2))-w1l(ixc^s,mag(ip2)))*signbx(ixc^s))*tmp(ixc^s)
944 w2l(ixc^s,mom(ip2))=w2r(ixc^s,mom(ip2))
946 w2r(ixc^s,mag(ip2))=(r1l(ixc^s)*w1r(ixc^s,mag(ip2))+r1r(ixc^s)*w1l(ixc^s,mag(ip2))+&
947 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip2))-w1l(ixc^s,mom(ip2)))*signbx(ixc^s))*tmp(ixc^s)
948 w2l(ixc^s,mag(ip2))=w2r(ixc^s,mag(ip2))
951 w2r(ixc^s,mom(ip3))=(r1l(ixc^s)*w1l(ixc^s,mom(ip3))+r1r(ixc^s)*w1r(ixc^s,mom(ip3))+&
952 (w1r(ixc^s,mag(ip3))-w1l(ixc^s,mag(ip3)))*signbx(ixc^s))*tmp(ixc^s)
953 w2l(ixc^s,mom(ip3))=w2r(ixc^s,mom(ip3))
955 w2r(ixc^s,mag(ip3))=(r1l(ixc^s)*w1r(ixc^s,mag(ip3))+r1r(ixc^s)*w1l(ixc^s,mag(ip3))+&
956 r1l(ixc^s)*r1r(ixc^s)*(w1r(ixc^s,mom(ip3))-w1l(ixc^s,mom(ip3)))*signbx(ixc^s))*tmp(ixc^s)
957 w2l(ixc^s,mag(ip3))=w2r(ixc^s,mag(ip3))
961 w2r(ixc^s,e_)=w1r(ixc^s,e_)+r1r(ixc^s)*(sum(w1r(ixc^s,mom(:))*w1r(ixc^s,mag(:)),dim=ndim+1)-&
962 sum(w2r(ixc^s,mom(:))*w2r(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
963 w2l(ixc^s,e_)=w1l(ixc^s,e_)-r1l(ixc^s)*(sum(w1l(ixc^s,mom(:))*w1l(ixc^s,mag(:)),dim=ndim+1)-&
964 sum(w2l(ixc^s,mom(:))*w2l(ixc^s,mag(:)),dim=ndim+1))*signbx(ixc^s)
969 w1r(ixc^s,mom(idir))=w1r(ixc^s,mom(idir))*w1r(ixc^s,rho_)
970 w1l(ixc^s,mom(idir))=w1l(ixc^s,mom(idir))*w1l(ixc^s,rho_)
971 w2r(ixc^s,mom(idir))=w2r(ixc^s,mom(idir))*w2r(ixc^s,rho_)
972 w2l(ixc^s,mom(idir))=w2l(ixc^s,mom(idir))*w2l(ixc^s,rho_)
981 f1l(ixc^s,iw)=flc(ixc^s,iw)
982 f1r(ixc^s,iw)=f1l(ixc^s,iw)
983 f2l(ixc^s,iw)=f1l(ixc^s,iw)
984 f2r(ixc^s,iw)=f1l(ixc^s,iw)
989 f1r(ixc^s,iw)=f1l(ixc^s,iw)
990 f2l(ixc^s,iw)=f1l(ixc^s,iw)
991 f2r(ixc^s,iw)=f1l(ixc^s,iw)
993 f1l(ixc^s,iw)=flc(ixc^s,iw)+sl(ixc^s,
index_v_mag)*(w1l(ixc^s,iw)-wlc(ixc^s,iw))
994 f1r(ixc^s,iw)=frc(ixc^s,iw)+sr(ixc^s,
index_v_mag)*(w1r(ixc^s,iw)-wrc(ixc^s,iw))
995 f2l(ixc^s,iw)=f1l(ixc^s,iw)+s1l(ixc^s)*(w2l(ixc^s,iw)-w1l(ixc^s,iw))
996 f2r(ixc^s,iw)=f1r(ixc^s,iw)+s1r(ixc^s)*(w2r(ixc^s,iw)-w1r(ixc^s,iw))
1001 {
do ix^db=ixcmin^db,ixcmax^db\}
1004 else if(s1l(ix^d)>=0.d0)
then
1006 else if(sm(ix^d)>=0.d0)
then
1008 else if(s1r(ix^d)>=0.d0)
then
1018 end subroutine get_riemann_flux_hlld_mag2
1024 integer,
intent(in) :: ixI^L, ixO^L
1025 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
1026 double precision,
intent(out):: csound(ixI^S)
1027 double precision :: cfast2(ixI^S), AvMinCs2(ixI^S), b2(ixI^S), kmax
1028 double precision :: inv_rho(ixO^S), gamma_A2(ixO^S)
1029 integer :: rho_, p_, e_, mom(1:ndir), mag(1:ndir)
1037 inv_rho=1.d0/w(ixo^s,rho_)
1042 b2(ixo^s) = sum((w(ixo^s, mag(:))+
block%B0(ixo^s,:,
b0i))**2, dim=ndim+1)
1044 b2(ixo^s) = sum(w(ixo^s, mag(:))**2, dim=ndim+1)
1049 avmincs2= w(ixo^s, mag(idims))+
block%B0(ixo^s,idims,
b0i)
1051 avmincs2= w(ixo^s, mag(idims))
1055 csound(ixo^s) = sum(w(ixo^s, mom(:))**2, dim=ndim+1)
1057 cfast2(ixo^s) = b2(ixo^s) * inv_rho+csound(ixo^s)
1058 avmincs2(ixo^s) = cfast2(ixo^s)**2-4.0d0*csound(ixo^s) &
1059 * avmincs2(ixo^s)**2 &
1062 where(avmincs2(ixo^s)<zero)
1063 avmincs2(ixo^s)=zero
1066 avmincs2(ixo^s)=sqrt(avmincs2(ixo^s))
1068 csound(ixo^s) = sqrt(half*(cfast2(ixo^s)+avmincs2(ixo^s)))
1076 subroutine reconstruct_lr(ixI^L,ixL^L,ixR^L,idims,w,wLC,wRC,wLp,wRp,x,dxdim)
1082 integer,
intent(in) :: ixi^
l, ixl^
l, ixr^
l, idims
1083 double precision,
intent(in) :: dxdim
1085 double precision,
dimension(ixI^S,1:nw) :: w
1087 double precision,
dimension(ixI^S,1:nw) :: wlc, wrc
1089 double precision,
dimension(ixI^S,1:nw) :: wlp, wrp
1090 double precision,
dimension(ixI^S,1:ndim) :: x
1092 integer :: jxr^
l, ixc^
l, jxc^
l, iw
1093 double precision :: ldw(ixi^s), rdw(ixi^s), dwc(ixi^s)
1094 double precision :: a2max
1098 call mp5limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp)
1100 call weno3limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,1)
1102 call weno3limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,2)
1104 call weno5limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,1)
1106 call weno5nmlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,1)
1108 call weno5limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,2)
1110 call weno5nmlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,2)
1112 call weno5limiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,3)
1114 call weno5nmlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp,3)
1116 call weno5cu6limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp)
1118 call teno5adlimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp)
1120 call weno7limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp,1)
1122 call weno7limiter(ixi^
l,ixl^
l,idims,w,wlp,wrp,2)
1124 call venklimiter(ixi^
l,ixl^
l,idims,dxdim,w,wlp,wrp)
1130 call ppmlimiter(ixi^
l,
ixm^
ll,idims,w,w,wlp,wrp)
1136 jxr^
l=ixr^
l+
kr(idims,^
d);
1137 ixcmax^
d=jxrmax^
d; ixcmin^
d=ixlmin^
d-
kr(idims,^
d);
1138 jxc^
l=ixc^
l+
kr(idims,^
d);
1141 w(ixcmin^
d:jxcmax^
d,iw)=dlog10(w(ixcmin^
d:jxcmax^
d,iw))
1142 wlp(ixl^s,iw)=dlog10(wlp(ixl^s,iw))
1143 wrp(ixr^s,iw)=dlog10(wrp(ixr^s,iw))
1146 dwc(ixc^s)=w(jxc^s,iw)-w(ixc^s,iw)
1162 call mpistop(
"idims is wrong in mod_limiter")
1168 wlp(ixl^s,iw)=wlp(ixl^s,iw)+half*ldw(ixl^s)
1169 wrp(ixr^s,iw)=wrp(ixr^s,iw)-half*rdw(jxr^s)
1172 w(ixcmin^
d:jxcmax^
d,iw)=10.0d0**w(ixcmin^
d:jxcmax^
d,iw)
1173 wlp(ixl^s,iw)=10.0d0**wlp(ixl^s,iw)
1174 wrp(ixr^s,iw)=10.0d0**wrp(ixr^s,iw)
1183 wlc(ixl^s,1:nwflux) = wlp(ixl^s,1:nwflux)
1184 wrc(ixr^s,1:nwflux) = wrp(ixr^s,1:nwflux)
1188 wlp(ixl^s,nwflux+1:nwflux+nwaux) = wlc(ixl^s,nwflux+1:nwflux+nwaux)
1189 wrp(ixr^s,nwflux+1:nwflux+nwaux) = wrc(ixr^s,nwflux+1:nwflux+nwaux)
subroutine get_hlld2_modif_c(w, x, ixIL, ixOL, csound)
Calculate fast magnetosonic wave speed.
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
integer, 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 slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
integer, parameter fs_hllc
double precision, dimension(ndim) a2max_global
global largest a2 for schmid scheme
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
integer, parameter flux_tvdlf
Indicates the flux should be treated with tvdlf.
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
integer, parameter flux_hll
Indicates the flux should be treated with hll.
procedure(sub_add_source_geom), pointer phys_add_source_geom
integer, parameter flux_special
Indicates the flux should be specially treated.
integer, dimension(:, :), allocatable flux_type
Array per direction per variable, which can be used to specify that certain fluxes have to be treated...
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_face_to_center), pointer phys_face_to_center
procedure(sub_modify_wlr), pointer phys_modify_wlr
procedure(sub_get_h_speed), pointer phys_get_h_speed
logical phys_energy
Solve energy equation or not.
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 iw_mom
Indices of the momentum density.
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, dimension(:), allocatable, protected iw_mag
Indices of the magnetic field components.
integer iw_rho
Index of the (gas) density.
integer nwflux
Number of flux variables.
integer index_v_mag
index of the var whose velocity appears in the induction eq.
integer iw_e
Index of the energy density.