52 integer,
intent(in) :: ixI^L, ixO^L
53 double precision,
intent(in) :: w(ixI^S,nw)
54 double precision,
intent(in) :: x(ixI^S,1:ndim)
55 double precision,
intent(out):: res(ixI^S)
65 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_from_eint => null()
66 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_from_conserved => null()
67 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_equi => null()
73 logical :: has_equi=.false.
77 double precision :: tc_k_para
80 double precision :: tc_k_perp
83 integer :: tc_slope_limiter
86 logical :: tc_constant=.false.
89 logical :: tc_perpendicular=.false.
92 logical :: tc_saturate=.false.
108 double precision,
intent(in) :: phys_gamma
118 subroutine read_mhd_params(fl)
123 end subroutine read_mhd_params
128 fl%tc_slope_limiter=1
134 call read_mhd_params(fl)
136 if(fl%tc_k_para==0.d0 .and. fl%tc_k_perp==0.d0)
then
148 if(
mype .eq. 0) print*,
"Spitzer MHD: par: ",fl%tc_k_para, &
149 " ,perp: ",fl%tc_k_perp
151 fl%tc_constant=.true.
164 subroutine read_hd_params(fl)
169 end subroutine read_hd_params
175 call read_hd_params(fl)
176 if(fl%tc_k_para==0.d0 )
then
184 if(
mype .eq. 0) print*,
"Spitzer HD par: ",fl%tc_k_para
196 integer,
intent(in) :: ixi^
l, ixo^
l
197 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
198 double precision,
intent(in) :: w(ixi^s,1:nw)
199 double precision :: dtnew
201 double precision :: mf(ixo^s,1:
ndim),te(ixi^s),b2(ixi^s),gradt(ixi^s)
202 double precision :: tmp2(ixo^s),tmp(ixo^s),hfs(ixo^s)
203 double precision :: dtdiff_tcond,maxtmp2
207 call fl%get_temperature_from_conserved(w,x,ixi^
l,ixi^
l,te)
213 mf(ixo^s,1:
ndim)=w(ixo^s,iw_mag(1:
ndim))
216 tmp=sum(mf**2,dim=
ndim+1)
218 where(tmp(ixo^s)/=0.d0)
219 ^
d&mf(ixo^s,^
d)=mf(ixo^s,^
d)**2/tmp(ixo^s);
221 ^
d&mf(ixo^s,^
d)=1.d0;
226 call fl%get_rho(w,x,ixi^
l,ixo^
l,b2)
229 if(fl%tc_constant)
then
230 tmp(ixo^s)=fl%tc_k_para
232 if(fl%tc_saturate)
then
240 hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)*sqrt(mf(ixo^s,idims))
243 tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)/(1.d0+4.2d0*tmp2(ixo^s)*dabs(hfs(ixo^s))/te(ixo^s))
246 tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)
253 tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*
dxlevel(idims))
254 maxtmp2=maxval(tmp2(ixo^s))
258 dtnew=min(dtnew,dtdiff_tcond)
263 tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*
block%ds(ixo^s,idims))
264 maxtmp2=maxval(tmp2(ixo^s)/
block%ds(ixo^s,idims))
266 dtdiff_tcond=1.d0/(
tc_gamma_1*maxtmp2+smalldouble)
268 dtnew=min(dtnew,dtdiff_tcond)
271 dtnew=dtnew/dble(
ndim)
280 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
281 double precision,
intent(in) :: x(ixi^s,1:
ndim)
282 double precision,
intent(in) :: w(ixi^s,1:nw)
283 double precision,
intent(inout) :: wres(ixi^s,1:nw)
284 double precision,
intent(in) :: my_dt
285 logical,
intent(in) :: fix_conserve_at_step
289 double precision :: qd(ixi^s)
290 double precision :: rho(ixi^s),te(ixi^s)
291 double precision :: qvec(ixi^s,1:
ndim)
292 double precision,
allocatable,
dimension(:^D&,:) :: qvec_equi
294 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
295 double precision :: alpha,dxinv(
ndim)
296 integer :: idims,idir,ix^
d,ix^
l,ixc^
l,ixa^
l,ixb^
l
308 call fl%get_temperature_from_eint(w, x, ixi^
l, ixi^
l, te)
309 call fl%get_rho(w, x, ixi^
l, ixi^
l, rho)
312 allocate(qvec_equi(ixi^s,1:
ndim))
313 call fl%get_temperature_equi(w, x, ixi^
l, ixi^
l, te)
314 call fl%get_rho_equi(w, x, ixi^
l, ixi^
l, rho)
317 qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
319 deallocate(qvec_equi)
325 qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
326 ixb^
l=ixo^
l-
kr(idims,^
d);
327 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
331 qvec(ix^s,idims)=qvec(ix^s,idims)*
block%surfaceC(ix^s,idims)
332 ixb^
l=ixo^
l-
kr(idims,^
d);
333 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
335 qd(ixo^s)=qd(ixo^s)/
block%dvolume(ixo^s)
338 if(fix_conserve_at_step)
then
339 allocate(fluxall(ixi^s,1,1:
ndim))
340 fluxall(ixi^s,1,1:
ndim)=my_dt*qvec(ixi^s,1:
ndim)
345 wres(ixo^s,fl%e_)=qd(ixo^s)
352 integer,
intent(in) :: ixI^L, ixO^L
353 double precision,
intent(in) :: x(ixI^S,1:ndim)
354 double precision,
intent(in) :: w(ixI^S,1:nw)
356 double precision,
intent(in) :: rho(ixI^S),Te(ixI^S)
357 double precision,
intent(in) :: alpha
358 double precision,
intent(out) :: qvec(ixI^S,1:ndim)
361 double precision :: qd(ixI^S)
363 double precision,
dimension(ixI^S,1:ndim) :: mf,Bc,Bcf
364 double precision,
dimension(ixI^S,1:ndim) :: gradT
365 double precision,
dimension(ixI^S) :: ka,kaf,ke,kef,qdd,qe,Binv,minq,maxq,Bnorm
366 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
367 integer,
dimension(ndim) :: lowindex
368 integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixA^D,ixB^D
375 mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))+
block%B0(ixi^s,1:ndim,0)
377 mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))
380 binv=dsqrt(sum(mf**2,dim=ndim+1))
381 {
do ix^db=ixmin^db,ixmax^db\}
382 if(binv(ix^d)/=0.d0)
then
383 binv(ix^d)=1.d0/binv(ix^d)
390 mf(ix^s,idims)=mf(ix^s,idims)*binv(ix^s)
393 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
397 ixamin^d=ixcmin^d+ix^d;
398 ixamax^d=ixcmax^d+ix^d;
399 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)+mf(ixa^s,1:ndim)
401 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)*0.5d0**ndim
405 ixbmax^d=ixmax^d-kr(idims,^d);
406 call gradientc(te,ixi^l,ixb^l,idims,minq)
407 gradt(ixb^s,idims)=minq(ixb^s)
409 if(fl%tc_constant)
then
410 if(fl%tc_perpendicular)
then
411 ka(ixc^s)=fl%tc_k_para-fl%tc_k_perp
412 ke(ixc^s)=fl%tc_k_perp
414 ka(ixc^s)=fl%tc_k_para
420 {
do ix^db=ixmin^db,ixmax^db\}
421 if(minq(ix^d) < block%wextra(ix^d,fl%Tcoff_))
then
422 minq(ix^d)=block%wextra(ix^d,fl%Tcoff_)
425 minq(ix^s)=fl%tc_k_para*sqrt(minq(ix^s)**5)
427 minq(ix^s)=fl%tc_k_para*sqrt(te(ix^s)**5)
431 ixbmin^d=ixcmin^d+ix^d;
432 ixbmax^d=ixcmax^d+ix^d;
433 ka(ixc^s)=ka(ixc^s)+minq(ixb^s)
436 ka(ixc^s)=0.5d0**ndim*ka(ixc^s)
438 if(fl%tc_perpendicular)
then
439 minq(ix^s)=fl%tc_k_perp*rho(ix^s)**2*binv(ix^s)**2/dsqrt(te(ix^s))
442 ixbmin^d=ixcmin^d+ix^d;
443 ixbmax^d=ixcmax^d+ix^d;
444 ke(ixc^s)=ke(ixc^s)+minq(ixb^s)
447 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
448 {
do ix^db=ixcmin^db,ixcmax^db\}
449 if(ke(ix^d)<ka(ix^d))
then
450 ka(ix^d)=ka(ix^d)-ke(ix^d)
458 if(fl%tc_slope_limiter==0)
then
464 if({ ix^d==0 .and. ^d==idims | .or.})
then
465 ixbmin^d=ixcmin^d+ix^d;
466 ixbmax^d=ixcmax^d+ix^d;
467 qd(ixc^s)=qd(ixc^s)+gradt(ixb^s,idims)
471 qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
476 qd(ixc^s)=qvec(ixc^s,idims)*bc(ixc^s,idims)+qd(ixc^s)
480 gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qd(ixc^s)
481 if(fl%tc_perpendicular) gradt(ixc^s,idims)=gradt(ixc^s,idims)+ke(ixc^s)*qvec(ixc^s,idims)
486 ixb^l=ixo^l-kr(idims,^d);
487 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
489 if({ ix^d==0 .and. ^d==idims | .or.})
then
490 ixbmin^d=ixamin^d-ix^d;
491 ixbmax^d=ixamax^d-ix^d;
492 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
495 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
496 if(fl%tc_saturate)
then
501 if({ ix^d==0 .and. ^d==idims | .or.})
then
502 ixbmin^d=ixamin^d-ix^d;
503 ixbmax^d=ixamax^d-ix^d;
504 bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
508 bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
509 ixb^l=ixa^l+kr(idims,^d);
510 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))
511 {
do ix^db=ixamin^db,ixamax^db\}
512 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
513 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
522 ixb^l=ixo^l-kr(idims,^d);
523 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
525 ixb^l=ixa^l+kr(idims,^d);
526 bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
531 if({ ix^d==0 .and. ^d==idims | .or.})
then
532 ixbmin^d=ixamin^d-ix^d;
533 ixbmax^d=ixamax^d-ix^d;
534 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)+bc(ixb^s,1:ndim)
535 kaf(ixa^s)=kaf(ixa^s)+ka(ixb^s)
536 if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)+ke(ixb^s)
540 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)*0.5d0**(ndim-1)
542 kaf(ixa^s)=kaf(ixa^s)*0.5d0**(ndim-1)
543 if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)*0.5d0**(ndim-1)
545 minq(ixa^s)=min(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
546 maxq(ixa^s)=max(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
551 if({ ix^d==0 .and. ^d==idims | .or.})
then
552 ixbmin^d=ixcmin^d+ix^d;
553 ixbmax^d=ixcmax^d+ix^d;
554 qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
558 qdd(ixc^s)=qdd(ixc^s)*0.5d0**(ndim-1)
563 if({ ix^d==0 .and. ^d==idims | .or.})
then
564 ixbmin^d=ixamin^d-ix^d;
565 ixbmax^d=ixamax^d-ix^d;
566 {
do ixa^db=ixamin^db,ixamax^db
567 ixb^db=ixa^db-ix^db\}
568 if(qd(ixb^d)<=minq(ixa^d))
then
569 qd(ixb^d)=minq(ixa^d)
570 else if(qd(ixb^d)>=maxq(ixa^d))
then
571 qd(ixb^d)=maxq(ixa^d)
574 qvec(ixa^s,idims)=qvec(ixa^s,idims)+bc(ixb^s,idims)**2*qd(ixb^s)
575 if(fl%tc_perpendicular) qe(ixa^s)=qe(ixa^s)+qd(ixb^s)
578 qvec(ixa^s,idims)=kaf(ixa^s)*qvec(ixa^s,idims)*0.5d0**(ndim-1)
580 if(fl%tc_perpendicular) qvec(ixa^s,idims)=qvec(ixa^s,idims)+kef(ixa^s)*qe(ixa^s)*0.5d0**(ndim-1)
583 ixbmax^d=ixamax^d+kr(idims,^d);
585 if(idir==idims) cycle
586 qd(ixi^s)=
slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
587 qd(ixi^s)=
slope_limiter(qd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
588 qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qd(ixa^s)
591 if(fl%tc_saturate)
then
594 ixb^l=ixa^l+kr(idims,^d);
595 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))
596 {
do ix^db=ixamin^db,ixamax^db\}
597 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
600 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
611 integer,
intent(in) :: ixi^
l, ixo^
l, idims, pm
612 double precision,
intent(in) :: f(ixi^s)
613 double precision :: lf(ixi^s)
614 integer,
intent(in) :: tc_slope_limiter
616 double precision :: signf(ixi^s)
619 ixb^
l=ixo^
l+pm*
kr(idims,^
d);
620 signf(ixo^s)=sign(1.d0,f(ixo^s))
621 select case(tc_slope_limiter)
624 lf(ixo^s)=two*signf(ixo^s)* &
625 max(zero,min(dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s),&
626 signf(ixo^s)*quarter*(f(ixb^s)+f(ixo^s))))
629 lf(ixo^s)=signf(ixo^s)*max(0.d0,min(abs(f(ixo^s)),signf(ixo^s)*f(ixb^s)))
632 lf(ixo^s)=signf(ixo^s)* &
633 max(zero,min(two*dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s)),&
634 min(dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s)))
637 lf(ixo^s)=signf(ixo^s)* &
638 max(zero,min(two*dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s),&
639 (two*f(ixb^s)*signf(ixo^s)+dabs(f(ixo^s)))*third))
641 call mpistop(
"Unknown slope limiter for thermal conduction")
650 integer,
intent(in) :: ixI^L, ixO^L, idir
651 double precision,
intent(in) :: q(ixI^S)
652 double precision,
intent(inout) :: gradq(ixI^S)
655 associate(x=>
block%x)
657 jxo^l=ixo^l+
kr(idir,^
d);
658 select case(coordinate)
660 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/
dxlevel(idir)
661 case(cartesian_stretched)
662 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
666 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,1)-x(ixo^s,1))
669 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/( (x(jxo^s,2)-x(ixo^s,2))*x(ixo^s,1) )
673 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)) )
678 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/((x(jxo^s,
phi_)-x(ixo^s,
phi_))*x(ixo^s,
r_))
680 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
683 call mpistop(
'Unknown geometry')
693 integer,
intent(in) :: ixi^
l, ixo^
l
694 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
695 double precision,
intent(in) :: w(ixi^s,1:nw)
697 double precision :: dtnew
699 double precision :: tmp(ixo^s),tmp2(ixo^s),te(ixi^s),rho(ixi^s),hfs(ixo^s),gradt(ixi^s)
700 double precision :: dtdiff_tcond,maxtmp2
703 call fl%get_temperature_from_conserved(w,x,ixi^
l,ixi^
l,te)
704 call fl%get_rho(w,x,ixi^
l,ixo^
l,rho)
706 if(fl%tc_saturate)
then
713 call gradient(te,ixi^
l,ixo^
l,idim,gradt)
714 hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)**2
717 tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/(rho(ixo^s)*(1.d0+4.2d0*tmp2(ixo^s)*sqrt(hfs(ixo^s))/te(ixo^s)))
719 tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/rho(ixo^s)
726 tmp2(ixo^s)=tmp(ixo^s)/
dxlevel(idim)
727 maxtmp2=maxval(tmp2(ixo^s))
731 dtnew=min(dtnew,dtdiff_tcond)
736 tmp2(ixo^s)=tmp(ixo^s)/
block%ds(ixo^s,idim)
737 maxtmp2=maxval(tmp2(ixo^s)/
block%ds(ixo^s,idim))
739 dtdiff_tcond=1.d0/(
tc_gamma_1*maxtmp2+smalldouble)
741 dtnew=min(dtnew,dtdiff_tcond)
744 dtnew=dtnew/dble(
ndim)
752 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
753 double precision,
intent(in) :: x(ixi^s,1:
ndim)
754 double precision,
intent(in) :: w(ixi^s,1:nw)
755 double precision,
intent(inout) :: wres(ixi^s,1:nw)
756 double precision,
intent(in) :: my_dt
757 logical,
intent(in) :: fix_conserve_at_step
760 double precision :: te(ixi^s),rho(ixi^s)
761 double precision :: qvec(ixi^s,1:
ndim),qd(ixi^s)
762 double precision,
allocatable,
dimension(:^D&,:) :: qvec_equi
763 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
765 double precision :: dxinv(
ndim)
766 integer :: idims,ix^
l,ixb^
l
773 call fl%get_temperature_from_eint(w, x, ixi^
l, ixi^
l, te)
774 call fl%get_rho(w, x, ixi^
l, ixi^
l, rho)
777 allocate(qvec_equi(ixi^s,1:
ndim))
778 call fl%get_temperature_equi(w, x, ixi^
l, ixi^
l, te)
779 call fl%get_rho_equi(w, x, ixi^
l, ixi^
l, rho)
782 qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
784 deallocate(qvec_equi)
791 qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
792 ixb^
l=ixo^
l-
kr(idims,^
d);
793 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
797 qvec(ix^s,idims)=qvec(ix^s,idims)*
block%surfaceC(ix^s,idims)
798 ixb^
l=ixo^
l-
kr(idims,^
d);
799 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
801 qd(ixo^s)=qd(ixo^s)/
block%dvolume(ixo^s)
804 if(fix_conserve_at_step)
then
805 allocate(fluxall(ixi^s,1,1:
ndim))
806 fluxall(ixi^s,1,1:
ndim)=my_dt*qvec(ixi^s,1:
ndim)
810 wres(ixo^s,fl%e_)=qd(ixo^s)
817 integer,
intent(in) :: ixI^L, ixO^L
818 double precision,
intent(in) :: x(ixI^S,1:ndim)
819 double precision,
intent(in) :: w(ixI^S,1:nw)
821 double precision,
intent(in) :: Te(ixI^S),rho(ixI^S)
822 double precision,
intent(out) :: qvec(ixI^S,1:ndim)
825 double precision :: gradT(ixI^S,1:ndim),ke(ixI^S),qd(ixI^S)
827 integer :: idims,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixD^L
831 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
865 ixbmax^d=ixmax^d-
kr(idims,^d);
869 if({ix^d==0 .and. ^d==idims |.or. })
then
870 ixbmin^d=ixcmin^d+ix^d;
871 ixbmax^d=ixcmax^d+ix^d;
872 qd(ixc^s)=qd(ixc^s)+ke(ixb^s)
876 qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
881 {
do ix^db=ixmin^db,ixmax^db\}
882 if(te(ix^d) < block%wextra(ix^d,fl%Tcoff_))
then
883 qd(ix^d)=fl%tc_k_para*dsqrt(block%wextra(ix^d,fl%Tcoff_))**5
885 qd(ix^d)=fl%tc_k_para*dsqrt(te(ix^d))**5
889 qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s))**5
893 ixbmin^d=ixcmin^d+ix^d;
894 ixbmax^d=ixcmax^d+ix^d;
895 ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
898 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
901 gradt(ixc^s,idims)=ke(ixc^s)*qvec(ixc^s,idims)
904 if(fl%tc_saturate)
then
907 qd(ix^s)=5.5d0*rho(ix^s)*dsqrt(te(ix^s)**3)
911 ixbmin^d=ixcmin^d+ix^d;
912 ixbmax^d=ixcmax^d+ix^d;
913 ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
916 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
918 qd(ixc^s)=norm2(gradt(ixc^s,:),dim=ndim+1)
919 {
do ix^db=ixcmin^db,ixcmax^db\}
920 if(qd(ix^d)>ke(ix^d))
then
921 ke(ix^d)=ke(ix^d)/(qd(ix^d)+smalldouble)
923 gradt(ix^d,idims)=ke(ix^d)*gradt(ix^d,idims)
932 ixb^l=ixo^l-kr(idims,^d);
933 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
935 if({ ix^d==0 .and. ^d==idims | .or.})
then
936 ixbmin^d=ixamin^d-ix^d;
937 ixbmax^d=ixamax^d-ix^d;
938 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
941 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
subroutine, public 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)
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
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
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)