51 integer,
intent(in) :: ixI^L, ixO^L
52 double precision,
intent(in) :: w(ixI^S,nw)
53 double precision,
intent(in) :: x(ixI^S,1:ndim)
54 double precision,
intent(out):: res(ixI^S)
64 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_from_eint => null()
65 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_from_conserved => null()
66 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_equi => null()
72 logical :: has_equi=.false.
76 double precision :: tc_k_para
79 double precision :: tc_k_perp
82 character(len=std_len) :: tc_slope_limiter
84 logical :: tc_constant=.false.
86 logical :: tc_perpendicular=.false.
89 logical :: tc_saturate=.true.
107 double precision,
intent(in) :: phys_gamma
119 subroutine read_mhd_params(fl)
124 end subroutine read_mhd_params
129 fl%tc_slope_limiter=
'MC'
135 call read_mhd_params(fl)
137 if(fl%tc_k_para==0.d0 .and. fl%tc_k_perp==0.d0)
then
149 if(
mype .eq. 0) print*,
"Spitzer MHD: par: ",fl%tc_k_para, &
150 " ,perp: ",fl%tc_k_perp
152 fl%tc_constant=.true.
166 subroutine read_hd_params(fl)
171 end subroutine read_hd_params
177 call read_hd_params(fl)
178 if(fl%tc_k_para==0.d0 )
then
186 if(
mype .eq. 0) print*,
"Spitzer HD par: ",fl%tc_k_para
198 integer,
intent(in) :: ixi^
l, ixo^
l
199 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
200 double precision,
intent(in) :: w(ixi^s,1:nw)
201 double precision :: dtnew
203 double precision :: dxinv(1:
ndim),mf(ixi^s,1:
ndir)
204 double precision :: tmp2(ixi^s),tmp(ixi^s),te(ixi^s),b2(ixi^s)
205 double precision :: dtdiff_tcond,maxtmp2
208 ^
d&dxinv(^
d)=one/
dx^
d;
211 call fl%get_temperature_from_conserved(w,x,ixi^
l,ixo^
l,te)
214 if(fl%tc_constant)
then
215 tmp(ixo^s)=fl%tc_k_para
217 call fl%get_rho(w,x,ixi^
l,ixo^
l,tmp2)
218 tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)/tmp2(ixo^s)
223 mf(ixo^s,:)=w(ixo^s,iw_mag(:))+
block%B0(ixo^s,:,0)
225 mf(ixo^s,:)=w(ixo^s,iw_mag(:))
228 b2(ixo^s)=sum(mf(ixo^s,:)**2,dim=
ndim+1)
230 where(b2(ixo^s)/=0.d0)
231 ^
d&mf(ixo^s,^
d)=mf(ixo^s,^
d)**2/b2(ixo^s);
233 ^
d&mf(ixo^s,^
d)=1.d0;
236 if(fl%tc_saturate) b2(ixo^s)=22.d0*dsqrt(te(ixo^s))
239 tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idim)
240 if(fl%tc_saturate)
then
241 where(tmp2(ixo^s)>b2(ixo^s))
242 tmp2(ixo^s)=b2(ixo^s)
245 maxtmp2=maxval(tmp2(ixo^s))
246 if(maxtmp2==0.d0) maxtmp2=smalldouble
248 dtdiff_tcond=1.d0/
tc_gamma_1/(maxtmp2*dxinv(idim)**2)
250 dtnew=min(dtnew,dtdiff_tcond)
252 dtnew=dtnew/dble(
ndim)
261 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
262 double precision,
intent(in) :: x(ixi^s,1:
ndim)
263 double precision,
intent(in) :: w(ixi^s,1:nw)
264 double precision,
intent(inout) :: wres(ixi^s,1:nw)
265 double precision,
intent(in) :: my_dt
266 logical,
intent(in) :: fix_conserve_at_step
270 double precision :: qd(ixi^s)
271 double precision :: rho(ixi^s),te(ixi^s)
272 double precision :: qvec(ixi^s,1:
ndim)
273 double precision,
allocatable,
dimension(:^D&,:) :: qvec_equi
275 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
276 double precision :: alpha,dxinv(
ndim)
277 integer :: idims,idir,ix^
d,ix^
l,ixc^
l,ixa^
l,ixb^
l
289 call fl%get_temperature_from_eint(w, x, ixi^
l, ixi^
l, te)
290 call fl%get_rho(w, x, ixi^
l, ixi^
l, rho)
293 allocate(qvec_equi(ixi^s,1:
ndim))
294 call fl%get_temperature_equi(w, x, ixi^
l, ixi^
l, te)
295 call fl%get_rho_equi(w, x, ixi^
l, ixi^
l, rho)
298 qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
300 deallocate(qvec_equi)
306 qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
307 ixb^
l=ixo^
l-
kr(idims,^
d);
308 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
312 qvec(ix^s,idims)=qvec(ix^s,idims)*
block%surfaceC(ix^s,idims)
313 ixb^
l=ixo^
l-
kr(idims,^
d);
314 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
316 qd(ixo^s)=qd(ixo^s)/
block%dvolume(ixo^s)
319 if(fix_conserve_at_step)
then
320 allocate(fluxall(ixi^s,1,1:
ndim))
321 fluxall(ixi^s,1,1:
ndim)=my_dt*qvec(ixi^s,1:
ndim)
326 wres(ixo^s,fl%e_)=qd(ixo^s)
333 integer,
intent(in) :: ixI^L, ixO^L
334 double precision,
intent(in) :: x(ixI^S,1:ndim)
335 double precision,
intent(in) :: w(ixI^S,1:nw)
337 double precision,
intent(in) :: rho(ixI^S),Te(ixI^S)
338 double precision,
intent(in) :: alpha
339 double precision,
intent(out) :: qvec(ixI^S,1:ndim)
342 double precision :: qd(ixI^S)
344 double precision,
dimension(ixI^S,1:ndir) :: mf,Bc,Bcf
345 double precision,
dimension(ixI^S,1:ndim) :: gradT
346 double precision,
dimension(ixI^S) :: ka,kaf,ke,kef,qdd,qe,Binv,minq,maxq,Bnorm
347 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
348 integer,
dimension(ndim) :: lowindex
349 integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L
356 mf(ixi^s,:)=w(ixi^s,iw_mag(:))+
block%B0(ixi^s,:,0)
358 mf(ixi^s,:)=w(ixi^s,iw_mag(:))
361 binv(ix^s)=dsqrt(sum(mf(ix^s,:)**2,dim=ndim+1))
362 where(binv(ix^s)/=0.d0)
363 binv(ix^s)=1.d0/binv(ix^s)
369 mf(ix^s,idims)=mf(ix^s,idims)*binv(ix^s)
372 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
376 ixamin^d=ixcmin^d+ix^d;
377 ixamax^d=ixcmax^d+ix^d;
378 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)+mf(ixa^s,1:ndim)
380 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)*0.5d0**ndim
385 ixbmax^d=ixmax^d-kr(idims,^d);
386 call gradientc(te,ixi^l,ixb^l,idims,minq)
387 gradt(ixb^s,idims)=minq(ixb^s)
389 if(fl%tc_constant)
then
390 if(fl%tc_perpendicular)
then
391 ka(ixc^s)=fl%tc_k_para-fl%tc_k_perp
392 ke(ixc^s)=fl%tc_k_perp
394 ka(ixc^s)=fl%tc_k_para
400 where(minq(ix^s) < block%wextra(ix^s,fl%Tcoff_))
401 minq(ix^s)=block%wextra(ix^s,fl%Tcoff_)
403 minq(ix^s)=fl%tc_k_para*sqrt(minq(ix^s)**5)
405 minq(ix^s)=fl%tc_k_para*sqrt(te(ix^s)**5)
409 ixbmin^d=ixcmin^d+ix^d;
410 ixbmax^d=ixcmax^d+ix^d;
411 ka(ixc^s)=ka(ixc^s)+minq(ixb^s)
414 ka(ixc^s)=0.5d0**ndim*ka(ixc^s)
416 if(fl%tc_perpendicular)
then
417 minq(ix^s)=fl%tc_k_perp*rho(ix^s)**2*binv(ix^s)**2/dsqrt(te(ix^s))
420 ixbmin^d=ixcmin^d+ix^d;
421 ixbmax^d=ixcmax^d+ix^d;
422 ke(ixc^s)=ke(ixc^s)+minq(ixb^s)
425 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
426 where(ke(ixc^s)<ka(ixc^s))
427 ka(ixc^s)=ka(ixc^s)-ke(ixc^s)
434 if(fl%tc_slope_limiter==
'no')
then
440 if({ ix^d==0 .and. ^d==idims | .or.})
then
441 ixbmin^d=ixcmin^d+ix^d;
442 ixbmax^d=ixcmax^d+ix^d;
443 qd(ixc^s)=qd(ixc^s)+gradt(ixb^s,idims)
447 qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
450 qd(ixc^s)=sum(qvec(ixc^s,1:ndim)*bc(ixc^s,1:ndim),dim=ndim+1)
453 gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qd(ixc^s)
454 if(fl%tc_perpendicular) gradt(ixc^s,idims)=gradt(ixc^s,idims)+ke(ixc^s)*qvec(ixc^s,idims)
459 ixb^l=ixo^l-kr(idims,^d);
460 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
462 if({ ix^d==0 .and. ^d==idims | .or.})
then
463 ixbmin^d=ixamin^d-ix^d;
464 ixbmax^d=ixamax^d-ix^d;
465 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
468 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
469 if(fl%tc_saturate)
then
474 if({ ix^d==0 .and. ^d==idims | .or.})
then
475 ixbmin^d=ixamin^d-ix^d;
476 ixbmax^d=ixamax^d-ix^d;
477 bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
481 bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
482 ixb^l=ixa^l+kr(idims,^d);
483 qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bcf(ixa^s,idims))
484 {
do ix^db=ixamin^db,ixamax^db\}
485 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
486 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
495 ixb^l=ixo^l-kr(idims,^d);
496 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
498 ixb^l=ixa^l+kr(idims,^d);
499 bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
504 if({ ix^d==0 .and. ^d==idims | .or.})
then
505 ixbmin^d=ixamin^d-ix^d;
506 ixbmax^d=ixamax^d-ix^d;
507 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)+bc(ixb^s,1:ndim)
508 kaf(ixa^s)=kaf(ixa^s)+ka(ixb^s)
509 if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)+ke(ixb^s)
513 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)*0.5d0**(ndim-1)
515 kaf(ixa^s)=kaf(ixa^s)*0.5d0**(ndim-1)
516 if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)*0.5d0**(ndim-1)
518 minq(ixa^s)=min(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
519 maxq(ixa^s)=max(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
524 if({ ix^d==0 .and. ^d==idims | .or.})
then
525 ixbmin^d=ixcmin^d+ix^d;
526 ixbmax^d=ixcmax^d+ix^d;
527 qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
531 qdd(ixc^s)=qdd(ixc^s)*0.5d0**(ndim-1)
536 if({ ix^d==0 .and. ^d==idims | .or.})
then
537 ixbmin^d=ixamin^d-ix^d;
538 ixbmax^d=ixamax^d-ix^d;
539 where(qd(ixb^s)<=minq(ixa^s))
540 qd(ixb^s)=minq(ixa^s)
541 elsewhere(qd(ixb^s)>=maxq(ixa^s))
542 qd(ixb^s)=maxq(ixa^s)
544 qvec(ixa^s,idims)=qvec(ixa^s,idims)+bc(ixb^s,idims)**2*qd(ixb^s)
545 if(fl%tc_perpendicular) qe(ixa^s)=qe(ixa^s)+qd(ixb^s)
548 qvec(ixa^s,idims)=kaf(ixa^s)*qvec(ixa^s,idims)*0.5d0**(ndim-1)
550 if(fl%tc_perpendicular) qvec(ixa^s,idims)=qvec(ixa^s,idims)+kef(ixa^s)*qe(ixa^s)*0.5d0**(ndim-1)
553 ixbmax^d=ixamax^d+kr(idims,^d);
555 if(idir==idims) cycle
556 qd(ixi^s)=
slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
557 qd(ixi^s)=
slope_limiter(qd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
558 qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qd(ixa^s)
566 if(fl%tc_saturate)
then
569 ixb^l=ixa^l+kr(idims,^d);
570 qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bnorm(ixa^s))
571 {
do ix^db=ixamin^db,ixamax^db\}
572 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
575 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
586 integer,
intent(in) :: ixi^
l, ixo^
l, idims, pm
587 double precision,
intent(in) :: f(ixi^s)
588 double precision :: lf(ixi^s)
589 character(len=*),
intent(in) :: tc_slope_limiter
591 double precision :: signf(ixi^s)
594 ixb^
l=ixo^
l+pm*
kr(idims,^
d);
595 signf(ixo^s)=sign(1.d0,f(ixo^s))
596 select case(tc_slope_limiter)
599 lf(ixo^s)=signf(ixo^s)*max(0.d0,min(abs(f(ixo^s)),signf(ixo^s)*f(ixb^s)))
602 lf(ixo^s)=two*signf(ixo^s)* &
603 max(zero,min(dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s),&
604 signf(ixo^s)*quarter*(f(ixb^s)+f(ixo^s))))
607 lf(ixo^s)=signf(ixo^s)* &
608 max(zero,min(two*dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s)),&
609 min(dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s)))
612 lf(ixo^s)=signf(ixo^s)* &
613 max(zero,min(two*dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s),&
614 (two*f(ixb^s)*signf(ixo^s)+dabs(f(ixo^s)))*third))
616 call mpistop(
"Unknown slope limiter for thermal conduction")
625 integer,
intent(in) :: ixI^L, ixO^L, idir
626 double precision,
intent(in) :: q(ixI^S)
627 double precision,
intent(inout) :: gradq(ixI^S)
630 associate(x=>
block%x)
632 jxo^l=ixo^l+
kr(idir,^
d);
633 select case(coordinate)
635 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/
dxlevel(idir)
636 case(cartesian_stretched)
637 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
641 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,1)-x(ixo^s,1))
644 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/( (x(jxo^s,2)-x(ixo^s,2))*x(ixo^s,1) )
648 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/( (x(jxo^s,3)-x(ixo^s,3))*x(ixo^s,1)*dsin(x(ixo^s,2)) )
653 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/((x(jxo^s,
phi_)-x(ixo^s,
phi_))*x(ixo^s,
r_))
655 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
658 call mpistop(
'Unknown geometry')
668 integer,
intent(in) :: ixi^
l, ixo^
l
669 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
670 double precision,
intent(in) :: w(ixi^s,1:nw)
672 double precision :: dtnew
674 double precision :: dxinv(1:
ndim), tmp(ixi^s), te(ixi^s), rho(ixi^s)
675 double precision :: dtdiff_tcond,dtdiff_tsat
678 ^
d&dxinv(^
d)=one/
dx^
d;
680 call fl%get_temperature_from_conserved(w,x,ixi^
l,ixo^
l,te)
681 call fl%get_rho(w,x,ixi^
l,ixo^
l,rho)
683 tmp(ixo^s)=
tc_gamma_1*fl%tc_k_para*dsqrt((te(ixo^s))**5)/rho(ixo^s)
688 dtdiff_tcond=1d0/maxval(tmp(ixo^s)*dxinv(idim)**2)
689 if(fl%tc_saturate)
then
691 dtdiff_tsat=1d0/maxval(
tc_gamma_1*dsqrt(te(ixo^s))*&
694 dtdiff_tcond=max(dtdiff_tcond,dtdiff_tsat)
697 dtnew=min(dtnew,dtdiff_tcond)
699 dtnew=dtnew/dble(
ndim)
707 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
708 double precision,
intent(in) :: x(ixi^s,1:
ndim)
709 double precision,
intent(in) :: w(ixi^s,1:nw)
710 double precision,
intent(inout) :: wres(ixi^s,1:nw)
711 double precision,
intent(in) :: my_dt
712 logical,
intent(in) :: fix_conserve_at_step
715 double precision :: te(ixi^s),rho(ixi^s)
716 double precision :: qvec(ixi^s,1:
ndim),qd(ixi^s)
717 double precision,
allocatable,
dimension(:^D&,:) :: qvec_equi
718 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
720 double precision :: dxinv(
ndim)
721 integer :: idims,ix^
l,ixb^
l
728 call fl%get_temperature_from_eint(w, x, ixi^
l, ixi^
l, te)
729 call fl%get_rho(w, x, ixi^
l, ixi^
l, rho)
732 allocate(qvec_equi(ixi^s,1:
ndim))
733 call fl%get_temperature_equi(w, x, ixi^
l, ixi^
l, te)
734 call fl%get_rho_equi(w, x, ixi^
l, ixi^
l, rho)
737 qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
739 deallocate(qvec_equi)
746 qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
747 ixb^
l=ixo^
l-
kr(idims,^
d);
748 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
752 qvec(ix^s,idims)=qvec(ix^s,idims)*
block%surfaceC(ix^s,idims)
753 ixb^
l=ixo^
l-
kr(idims,^
d);
754 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
756 qd(ixo^s)=qd(ixo^s)/
block%dvolume(ixo^s)
759 if(fix_conserve_at_step)
then
760 allocate(fluxall(ixi^s,1,1:
ndim))
761 fluxall(ixi^s,1,1:
ndim)=my_dt*qvec(ixi^s,1:
ndim)
765 wres(ixo^s,fl%e_)=qd(ixo^s)
773 integer,
intent(in) :: ixI^L, ixO^L
774 double precision,
intent(in) :: x(ixI^S,1:ndim)
775 double precision,
intent(in) :: w(ixI^S,1:nw)
777 double precision,
intent(in) :: Te(ixI^S),rho(ixI^S)
778 double precision,
intent(out) :: qvec(ixI^S,1:ndim)
781 double precision :: gradT(ixI^S,1:ndim),ke(ixI^S),qd(ixI^S)
783 double precision :: dxinv(ndim)
784 integer :: idims,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixD^L
788 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
795 ixamax^d=ixmax^d; ixamin^d=ixmin^d-1;
797 ixbmin^d=ixamin^d+ix^d;
798 ixbmax^d=ixamax^d+ix^d;
799 ke(ixa^s)=ke(ixa^s)+te(ixb^s)
801 ke(ixa^s)=0.5d0**ndim*ke(ixa^s)
806 ixbmax^d=ixmax^d-kr(idims,^d);
807 call gradient(ke,ixi^l,ixb^l,idims,qd)
808 gradt(ixb^s,idims)=qd(ixb^s)
812 where(ke(ixi^s) < block%wextra(ixi^s,fl%Tcoff_))
813 ke(ixi^s)=block%wextra(ixi^s,fl%Tcoff_)
818 gradt(ixc^s,idims)=gradt(ixc^s,idims)*fl%tc_k_para*sqrt(ke(ixc^s)**5)
821 if(fl%tc_saturate)
then
824 qd(ix^s)=5.d0*rho(ix^s)*dsqrt(te(ix^s)**3)
828 ixbmin^d=ixcmin^d+ix^d;
829 ixbmax^d=ixcmax^d+ix^d;
830 ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
833 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
835 qd(ixc^s)=norm2(gradt(ixc^s,:),dim=ndim+1)
836 {
do ix^db=ixcmin^db,ixcmax^db\}
837 if(qd(ix^d)>ke(ix^d))
then
838 ke(ix^d)=ke(ix^d)/qd(ix^d)
840 gradt(ix^d,idims)=ke(ix^d)*gradt(ix^d,idims)
850 ixb^l=ixo^l-kr(idims,^d);
851 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
853 if({ ix^d==0 .and. ^d==idims | .or.})
then
854 ixbmin^d=ixamin^d-ix^d;
855 ixbmax^d=ixamax^d-ix^d;
856 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
859 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
subroutine mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for flux conservation near refinement boundaries.
subroutine, public store_flux(igrid, fC, idimLIM, nwfluxin)
Module with geometry-related routines (e.g., divergence, curl)
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.
double precision unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision unit_numberdensity
Physical scaling factor for number density.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
integer, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision unit_magneticfield
Physical scaling factor for magnetic field.
double precision unit_velocity
Physical scaling factor for velocity.
logical b0field
split magnetic field as background B0 field
double precision unit_temperature
Physical scaling factor for temperature.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
double precision, dimension(ndim) dxlevel
Thermal conduction for HD and MHD Adaptation of mod_thermal_conduction for the mod_supertimestepping ...
subroutine, public sts_set_source_tc_hd(ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
subroutine, public tc_get_hd_params(fl, read_hd_params)
Read tc module parameters from par file: MHD case.
subroutine set_source_tc_mhd(ixIL, ixOL, w, x, fl, qvec, rho, Te, alpha)
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_mhd(ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
anisotropic thermal conduction with slope limited symmetric scheme Sharma 2007 Journal of Computation...
subroutine gradientc(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q at cell interfaces in direction idir.
double precision function, public get_tc_dt_hd(w, ixIL, ixOL, dxD, x, fl)
subroutine set_source_tc_hd(ixIL, ixOL, w, x, fl, qvec, rho, Te)
double precision tc_gamma_1
The adiabatic index.
subroutine, public tc_get_mhd_params(fl, read_mhd_params)
Init TC coeffiecients: MHD case.
double precision function, dimension(ixi^s) slope_limiter(f, ixIL, ixOL, idims, pm, tc_slope_limiter)
double precision function, public get_tc_dt_mhd(w, ixIL, ixOL, dxD, x, fl)
Get the explicut timestep for the TC (mhd implementation)