56 integer,
intent(in) :: ixI^L, ixO^L
57 double precision,
intent(in) :: w(ixI^S,nw)
58 double precision,
intent(in) :: x(ixI^S,1:ndim)
59 double precision,
intent(out):: res(ixI^S)
67 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_from_eint => null()
68 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_from_conserved => null()
69 procedure(
get_var_subr),
pointer,
nopass :: get_temperature_equi => null()
75 logical :: has_equi=.false.
79 double precision :: tc_k_para
82 double precision :: tc_k_perp
85 integer :: tc_slope_limiter
88 logical :: tc_constant=.false.
91 logical :: tc_perpendicular=.false.
94 logical :: tc_saturate=.false.
113 double precision,
intent(in) :: phys_gamma
123 subroutine read_mhd_params(fl)
128 end subroutine read_mhd_params
132 fl%tc_slope_limiter=1
137 call read_mhd_params(fl)
139 if(fl%tc_k_para==0.d0 .and. fl%tc_k_perp==0.d0)
then
151 if(
mype .eq. 0) print*,
"Spitzer MHD: par: ",fl%tc_k_para, &
152 " ,perp: ",fl%tc_k_perp
154 fl%tc_constant=.true.
164 subroutine read_hd_params(fl)
169 end subroutine read_hd_params
176 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
195 subroutine read_ffhd_params(fl)
199 end subroutine read_ffhd_params
206 call read_ffhd_params(fl)
207 if(fl%tc_k_para==0.d0 )
then
215 if(
mype .eq. 0) print*,
"Spitzer ffHD par: ",fl%tc_k_para
217 fl%tc_constant=.true.
229 integer,
intent(in) :: ixi^
l, ixo^
l
230 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
231 double precision,
intent(in) :: w(ixi^s,1:nw)
232 double precision :: dtnew
234 double precision :: mf(ixo^s,1:
ndim),te(ixi^s),b2(ixi^s),gradt(ixi^s)
235 double precision :: tmp2(ixo^s),tmp(ixo^s),hfs(ixo^s)
236 double precision :: dtdiff_tcond,maxtmp2
240 call fl%get_temperature_from_conserved(w,x,ixi^
l,ixi^
l,te)
246 mf(ixo^s,1:
ndim)=w(ixo^s,iw_mag(1:
ndim))
249 b2(ixo^s)=sum(mf(ixo^s,1:
ndim)**2,dim=
ndim+1)
251 where(b2(ixo^s)/=0.d0)
252 ^
d&mf(ixo^s,^
d)=mf(ixo^s,^
d)**2/b2(ixo^s);
254 ^
d&mf(ixo^s,^
d)=1.d0;
259 call fl%get_rho(w,x,ixi^
l,ixo^
l,b2)
262 if(fl%tc_constant)
then
263 tmp(ixo^s)=fl%tc_k_para
265 if(fl%tc_saturate)
then
273 hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)*sqrt(mf(ixo^s,idims))
276 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))
279 tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)
286 tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*
dxlevel(idims))
287 maxtmp2=maxval(tmp2(ixo^s))
291 dtnew=min(dtnew,dtdiff_tcond)
296 tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*
block%ds(ixo^s,idims))
297 maxtmp2=maxval(tmp2(ixo^s)/
block%ds(ixo^s,idims))
299 dtdiff_tcond=1.d0/(
tc_gamma_1*maxtmp2+smalldouble)
301 dtnew=min(dtnew,dtdiff_tcond)
304 dtnew=dtnew/dble(
ndim)
312 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
313 double precision,
intent(in) :: x(ixi^s,1:
ndim)
314 double precision,
intent(in) :: w(ixi^s,1:nw)
315 double precision,
intent(inout) :: wres(ixi^s,1:nw)
316 double precision,
intent(in) :: my_dt
317 logical,
intent(in) :: fix_conserve_at_step
321 double precision :: qd(ixi^s)
322 double precision :: rho(ixi^s),te(ixi^s)
323 double precision :: qvec(ixi^s,1:
ndim)
324 double precision,
allocatable,
dimension(:^D&,:) :: qvec_equi
326 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
327 double precision :: alpha,dxinv(
ndim)
328 integer :: idims,idir,ix^
d,ix^
l,ixc^
l,ixa^
l,ixb^
l
340 call fl%get_temperature_from_eint(w, x, ixi^
l, ixi^
l, te)
341 call fl%get_rho(w, x, ixi^
l, ixi^
l, rho)
344 allocate(qvec_equi(ixi^s,1:
ndim))
345 call fl%get_temperature_equi(w, x, ixi^
l, ixi^
l, te)
346 call fl%get_rho_equi(w, x, ixi^
l, ixi^
l, rho)
349 qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
351 deallocate(qvec_equi)
357 qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
358 ixb^
l=ixo^
l-
kr(idims,^
d);
359 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
363 qvec(ix^s,idims)=qvec(ix^s,idims)*
block%surfaceC(ix^s,idims)
364 ixb^
l=ixo^
l-
kr(idims,^
d);
365 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
367 qd(ixo^s)=qd(ixo^s)/
block%dvolume(ixo^s)
370 if(fix_conserve_at_step)
then
371 allocate(fluxall(ixi^s,1,1:
ndim))
372 fluxall(ixi^s,1,1:
ndim)=my_dt*qvec(ixi^s,1:
ndim)
377 wres(ixo^s,fl%e_)=qd(ixo^s)
382 integer,
intent(in) :: ixI^L, ixO^L
383 double precision,
intent(in) :: x(ixI^S,1:ndim)
384 double precision,
intent(in) :: w(ixI^S,1:nw)
386 double precision,
intent(in) :: rho(ixI^S),Te(ixI^S)
387 double precision,
intent(in) :: alpha
388 double precision,
intent(out) :: qvec(ixI^S,1:ndim)
390 double precision :: qd(ixI^S)
391 double precision,
dimension(ixI^S,1:ndim) :: mf,Bc,Bcf
392 double precision,
dimension(ixI^S,1:ndim) :: gradT
393 double precision,
dimension(ixI^S) :: ka,kaf,ke,kef,qdd,qe,Binv,minq,maxq,Bnorm
394 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
395 integer,
dimension(ndim) :: lowindex
396 integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L
403 mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))+
block%B0(ixi^s,1:ndim,0)
405 mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))
408 binv(ix^s)=dsqrt(sum(mf(ix^s,1:ndim)**2,dim=ndim+1))
409 where(binv(ix^s)/=0.d0)
410 binv(ix^s)=1.d0/binv(ix^s)
416 mf(ix^s,idims)=mf(ix^s,idims)*binv(ix^s)
419 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
423 ixamin^d=ixcmin^d+ix^d;
424 ixamax^d=ixcmax^d+ix^d;
425 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)+mf(ixa^s,1:ndim)
427 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)*0.5d0**ndim
431 ixbmax^d=ixmax^d-kr(idims,^d);
432 call gradientc(te,ixi^l,ixb^l,idims,minq)
433 gradt(ixb^s,idims)=minq(ixb^s)
435 if(fl%tc_constant)
then
436 if(fl%tc_perpendicular)
then
437 ka(ixc^s)=fl%tc_k_para-fl%tc_k_perp
438 ke(ixc^s)=fl%tc_k_perp
440 ka(ixc^s)=fl%tc_k_para
446 where(minq(ix^s) < block%wextra(ix^s,fl%Tcoff_))
447 minq(ix^s)=block%wextra(ix^s,fl%Tcoff_)
449 minq(ix^s)=fl%tc_k_para*sqrt(minq(ix^s)**5)
451 minq(ix^s)=fl%tc_k_para*sqrt(te(ix^s)**5)
455 ixbmin^d=ixcmin^d+ix^d;
456 ixbmax^d=ixcmax^d+ix^d;
457 ka(ixc^s)=ka(ixc^s)+minq(ixb^s)
460 ka(ixc^s)=0.5d0**ndim*ka(ixc^s)
462 if(fl%tc_perpendicular)
then
463 minq(ix^s)=fl%tc_k_perp*rho(ix^s)**2*binv(ix^s)**2/dsqrt(te(ix^s))
466 ixbmin^d=ixcmin^d+ix^d;
467 ixbmax^d=ixcmax^d+ix^d;
468 ke(ixc^s)=ke(ixc^s)+minq(ixb^s)
471 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
472 where(ke(ixc^s)<ka(ixc^s))
473 ka(ixc^s)=ka(ixc^s)-ke(ixc^s)
480 if(fl%tc_slope_limiter==0)
then
486 if({ ix^d==0 .and. ^d==idims | .or.})
then
487 ixbmin^d=ixcmin^d+ix^d;
488 ixbmax^d=ixcmax^d+ix^d;
489 qd(ixc^s)=qd(ixc^s)+gradt(ixb^s,idims)
493 qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
496 qd(ixc^s)=sum(qvec(ixc^s,1:ndim)*bc(ixc^s,1:ndim),dim=ndim+1)
499 gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qd(ixc^s)
500 if(fl%tc_perpendicular) gradt(ixc^s,idims)=gradt(ixc^s,idims)+ke(ixc^s)*qvec(ixc^s,idims)
505 ixb^l=ixo^l-kr(idims,^d);
506 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
508 if({ ix^d==0 .and. ^d==idims | .or.})
then
509 ixbmin^d=ixamin^d-ix^d;
510 ixbmax^d=ixamax^d-ix^d;
511 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
514 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
515 if(fl%tc_saturate)
then
520 if({ ix^d==0 .and. ^d==idims | .or.})
then
521 ixbmin^d=ixamin^d-ix^d;
522 ixbmax^d=ixamax^d-ix^d;
523 bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
527 bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
528 ixb^l=ixa^l+kr(idims,^d);
529 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))
530 {
do ix^db=ixamin^db,ixamax^db\}
531 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
532 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
541 ixb^l=ixo^l-kr(idims,^d);
542 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
544 ixb^l=ixa^l+kr(idims,^d);
545 bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
550 if({ ix^d==0 .and. ^d==idims | .or.})
then
551 ixbmin^d=ixamin^d-ix^d;
552 ixbmax^d=ixamax^d-ix^d;
553 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)+bc(ixb^s,1:ndim)
554 kaf(ixa^s)=kaf(ixa^s)+ka(ixb^s)
555 if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)+ke(ixb^s)
559 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)*0.5d0**(ndim-1)
561 kaf(ixa^s)=kaf(ixa^s)*0.5d0**(ndim-1)
562 if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)*0.5d0**(ndim-1)
564 minq(ixa^s)=min(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
565 maxq(ixa^s)=max(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
570 if({ ix^d==0 .and. ^d==idims | .or.})
then
571 ixbmin^d=ixcmin^d+ix^d;
572 ixbmax^d=ixcmax^d+ix^d;
573 qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
577 qdd(ixc^s)=qdd(ixc^s)*0.5d0**(ndim-1)
582 if({ ix^d==0 .and. ^d==idims | .or.})
then
583 ixbmin^d=ixamin^d-ix^d;
584 ixbmax^d=ixamax^d-ix^d;
585 where(qd(ixb^s)<=minq(ixa^s))
586 qd(ixb^s)=minq(ixa^s)
587 elsewhere(qd(ixb^s)>=maxq(ixa^s))
588 qd(ixb^s)=maxq(ixa^s)
590 qvec(ixa^s,idims)=qvec(ixa^s,idims)+bc(ixb^s,idims)**2*qd(ixb^s)
591 if(fl%tc_perpendicular) qe(ixa^s)=qe(ixa^s)+qd(ixb^s)
594 qvec(ixa^s,idims)=kaf(ixa^s)*qvec(ixa^s,idims)*0.5d0**(ndim-1)
596 if(fl%tc_perpendicular) qvec(ixa^s,idims)=qvec(ixa^s,idims)+kef(ixa^s)*qe(ixa^s)*0.5d0**(ndim-1)
599 ixbmax^d=ixamax^d+kr(idims,^d);
601 if(idir==idims) cycle
602 qd(ixi^s)=
slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
603 qd(ixi^s)=
slope_limiter(qd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
604 qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qd(ixa^s)
607 if(fl%tc_saturate)
then
610 ixb^l=ixa^l+kr(idims,^d);
611 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))
612 {
do ix^db=ixamin^db,ixamax^db\}
613 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
616 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
626 integer,
intent(in) :: ixi^
l, ixo^
l, idims, pm
627 double precision,
intent(in) :: f(ixi^s)
628 double precision :: lf(ixi^s)
629 integer,
intent(in) :: tc_slope_limiter
631 double precision :: signf(ixi^s)
634 ixb^
l=ixo^
l+pm*
kr(idims,^
d);
635 signf(ixo^s)=sign(1.d0,f(ixo^s))
636 select case(tc_slope_limiter)
639 lf(ixo^s)=two*signf(ixo^s)* &
640 max(zero,min(dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s),&
641 signf(ixo^s)*quarter*(f(ixb^s)+f(ixo^s))))
644 lf(ixo^s)=signf(ixo^s)*max(0.d0,min(abs(f(ixo^s)),signf(ixo^s)*f(ixb^s)))
647 lf(ixo^s)=signf(ixo^s)* &
648 max(zero,min(two*dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s)),&
649 min(dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s)))
652 lf(ixo^s)=signf(ixo^s)* &
653 max(zero,min(two*dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s),&
654 (two*f(ixb^s)*signf(ixo^s)+dabs(f(ixo^s)))*third))
656 call mpistop(
"Unknown slope limiter for thermal conduction")
663 integer,
intent(in) :: ixI^L, ixO^L, idir
664 double precision,
intent(in) :: q(ixI^S)
665 double precision,
intent(inout) :: gradq(ixI^S)
668 associate(x=>
block%x)
670 jxo^l=ixo^l+
kr(idir,^
d);
671 select case(coordinate)
673 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/
dxlevel(idir)
674 case(cartesian_stretched)
675 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
679 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,1)-x(ixo^s,1))
682 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/( (x(jxo^s,2)-x(ixo^s,2))*x(ixo^s,1) )
686 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)) )
691 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/((x(jxo^s,
phi_)-x(ixo^s,
phi_))*x(ixo^s,
r_))
693 gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
696 call mpistop(
'Unknown geometry')
707 integer,
intent(in) :: ixi^
l, ixo^
l
708 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
709 double precision,
intent(in) :: w(ixi^s,1:nw)
711 double precision :: dtnew
713 double precision :: tmp(ixo^s),tmp2(ixo^s),te(ixi^s),rho(ixi^s),hfs(ixo^s),gradt(ixi^s)
714 double precision :: dtdiff_tcond,maxtmp2
717 call fl%get_temperature_from_conserved(w,x,ixi^
l,ixi^
l,te)
718 call fl%get_rho(w,x,ixi^
l,ixo^
l,rho)
720 if(fl%tc_saturate)
then
727 call gradient(te,ixi^
l,ixo^
l,idim,gradt)
728 hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)**2
731 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)))
733 tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/rho(ixo^s)
740 tmp2(ixo^s)=tmp(ixo^s)/
dxlevel(idim)
741 maxtmp2=maxval(tmp2(ixo^s))
745 dtnew=min(dtnew,dtdiff_tcond)
750 tmp2(ixo^s)=tmp(ixo^s)/
block%ds(ixo^s,idim)
751 maxtmp2=maxval(tmp2(ixo^s)/
block%ds(ixo^s,idim))
753 dtdiff_tcond=1.d0/(
tc_gamma_1*maxtmp2+smalldouble)
755 dtnew=min(dtnew,dtdiff_tcond)
758 dtnew=dtnew/dble(
ndim)
765 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
766 double precision,
intent(in) :: x(ixi^s,1:
ndim)
767 double precision,
intent(in) :: w(ixi^s,1:nw)
768 double precision,
intent(inout) :: wres(ixi^s,1:nw)
769 double precision,
intent(in) :: my_dt
770 logical,
intent(in) :: fix_conserve_at_step
773 double precision :: te(ixi^s),rho(ixi^s)
774 double precision :: qvec(ixi^s,1:
ndim),qd(ixi^s)
775 double precision,
allocatable,
dimension(:^D&,:) :: qvec_equi
776 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
778 double precision :: dxinv(
ndim)
779 integer :: idims,ix^
l,ixb^
l
786 call fl%get_temperature_from_eint(w, x, ixi^
l, ixi^
l, te)
787 call fl%get_rho(w, x, ixi^
l, ixi^
l, rho)
790 allocate(qvec_equi(ixi^s,1:
ndim))
791 call fl%get_temperature_equi(w, x, ixi^
l, ixi^
l, te)
792 call fl%get_rho_equi(w, x, ixi^
l, ixi^
l, rho)
795 qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
797 deallocate(qvec_equi)
803 qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
804 ixb^
l=ixo^
l-
kr(idims,^
d);
805 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
809 qvec(ix^s,idims)=qvec(ix^s,idims)*
block%surfaceC(ix^s,idims)
810 ixb^
l=ixo^
l-
kr(idims,^
d);
811 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
813 qd(ixo^s)=qd(ixo^s)/
block%dvolume(ixo^s)
816 if(fix_conserve_at_step)
then
817 allocate(fluxall(ixi^s,1,1:
ndim))
818 fluxall(ixi^s,1,1:
ndim)=my_dt*qvec(ixi^s,1:
ndim)
822 wres(ixo^s,fl%e_)=qd(ixo^s)
827 integer,
intent(in) :: ixI^L, ixO^L
828 double precision,
intent(in) :: x(ixI^S,1:ndim)
829 double precision,
intent(in) :: w(ixI^S,1:nw)
831 double precision,
intent(in) :: Te(ixI^S),rho(ixI^S)
832 double precision,
intent(out) :: qvec(ixI^S,1:ndim)
833 double precision :: gradT(ixI^S,1:ndim),ke(ixI^S),qd(ixI^S)
834 integer :: idims,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixD^L
838 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
844 ixbmax^d=ixmax^d-
kr(idims,^d);
848 if({ix^d==0 .and. ^d==idims |.or. })
then
849 ixbmin^d=ixcmin^d+ix^d;
850 ixbmax^d=ixcmax^d+ix^d;
851 qd(ixc^s)=qd(ixc^s)+ke(ixb^s)
855 qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
860 where(te(ix^s) < block%wextra(ix^s,fl%Tcoff_))
861 qd(ix^s)=fl%tc_k_para*dsqrt(block%wextra(ix^s,fl%Tcoff_))**5
863 qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s))**5
866 qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s))**5
870 ixbmin^d=ixcmin^d+ix^d;
871 ixbmax^d=ixcmax^d+ix^d;
872 ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
875 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
878 gradt(ixc^s,idims)=ke(ixc^s)*qvec(ixc^s,idims)
881 if(fl%tc_saturate)
then
884 qd(ix^s)=5.5d0*rho(ix^s)*dsqrt(te(ix^s)**3)
888 ixbmin^d=ixcmin^d+ix^d;
889 ixbmax^d=ixcmax^d+ix^d;
890 ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
893 ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
895 qd(ixc^s)=norm2(gradt(ixc^s,:),dim=ndim+1)
896 {
do ix^db=ixcmin^db,ixcmax^db\}
897 if(qd(ix^d)>ke(ix^d))
then
898 ke(ix^d)=ke(ix^d)/(qd(ix^d)+smalldouble)
900 gradt(ix^d,idims)=ke(ix^d)*gradt(ix^d,idims)
909 ixb^l=ixo^l-kr(idims,^d);
910 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
912 if({ ix^d==0 .and. ^d==idims | .or.})
then
913 ixbmin^d=ixamin^d-ix^d;
914 ixbmax^d=ixamax^d-ix^d;
915 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
918 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
929 integer,
intent(in) :: ixi^
l, ixo^
l
930 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
931 double precision,
intent(in) :: w(ixi^s,1:nw)
932 double precision :: dtnew
933 double precision :: mf(ixo^s,1:
ndim),te(ixi^s),b2(ixi^s),gradt(ixi^s)
934 double precision :: tmp2(ixo^s),tmp(ixo^s),hfs(ixo^s)
935 double precision :: dtdiff_tcond,maxtmp2
939 call fl%get_temperature_from_conserved(w,x,ixi^
l,ixi^
l,te)
944 call fl%get_rho(w,x,ixi^
l,ixo^
l,b2)
947 if(fl%tc_constant)
then
948 tmp(ixo^s)=fl%tc_k_para
950 if(fl%tc_saturate)
then
957 call gradient(te,ixi^
l,ixo^
l,idims,gradt)
958 hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)*sqrt(mf(ixo^s,idims))
961 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))
964 tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)
971 tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*
dxlevel(idims))
972 maxtmp2=maxval(tmp2(ixo^s))
976 dtnew=min(dtnew,dtdiff_tcond)
981 tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*
block%ds(ixo^s,idims))
982 maxtmp2=maxval(tmp2(ixo^s)/
block%ds(ixo^s,idims))
984 dtdiff_tcond=1.d0/(
tc_gamma_1*maxtmp2+smalldouble)
986 dtnew=min(dtnew,dtdiff_tcond)
989 dtnew=dtnew/dble(
ndim)
997 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
998 double precision,
intent(in) :: x(ixi^s,1:
ndim)
999 double precision,
intent(in) :: w(ixi^s,1:nw)
1000 double precision,
intent(inout) :: wres(ixi^s,1:nw)
1001 double precision,
intent(in) :: my_dt
1002 logical,
intent(in) :: fix_conserve_at_step
1005 double precision :: qd(ixi^s)
1006 double precision :: rho(ixi^s),te(ixi^s)
1007 double precision :: qvec(ixi^s,1:
ndim)
1008 double precision,
allocatable,
dimension(:^D&,:) :: qvec_equi
1009 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
1010 double precision :: alpha,dxinv(
ndim)
1011 integer :: idims,idir,ix^
d,ix^
l,ixc^
l,ixa^
l,ixb^
l
1021 call fl%get_temperature_from_eint(w, x, ixi^
l, ixi^
l, te)
1022 call fl%get_rho(w, x, ixi^
l, ixi^
l, rho)
1028 qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
1029 ixb^
l=ixo^
l-
kr(idims,^
d);
1030 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
1034 qvec(ix^s,idims)=qvec(ix^s,idims)*
block%surfaceC(ix^s,idims)
1035 ixb^
l=ixo^
l-
kr(idims,^
d);
1036 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
1038 qd(ixo^s)=qd(ixo^s)/
block%dvolume(ixo^s)
1041 if(fix_conserve_at_step)
then
1042 allocate(fluxall(ixi^s,1,1:
ndim))
1043 fluxall(ixi^s,1,1:
ndim)=my_dt*qvec(ixi^s,1:
ndim)
1047 wres(ixo^s,fl%e_)=qd(ixo^s)
1052 integer,
intent(in) :: ixI^L, ixO^L
1053 double precision,
intent(in) :: x(ixI^S,1:ndim)
1054 double precision,
intent(in) :: w(ixI^S,1:nw)
1056 double precision,
intent(in) :: rho(ixI^S),Te(ixI^S)
1057 double precision,
intent(in) :: alpha
1058 double precision,
intent(out) :: qvec(ixI^S,1:ndim)
1060 double precision :: qd(ixI^S)
1061 double precision,
dimension(ixI^S,1:ndim) :: mf,Bc,Bcf
1062 double precision,
dimension(ixI^S,1:ndim) :: gradT
1063 double precision,
dimension(ixI^S) :: ka,kaf,ke,kef,qdd,qe,Binv,minq,maxq,Bnorm
1064 double precision,
allocatable,
dimension(:^D&,:,:) :: fluxall
1065 integer,
dimension(ndim) :: lowindex
1066 integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixA^D,ixB^D
1072 mf(ixi^s,1:ndim)=
block%B0(ixi^s,1:ndim,0)
1074 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
1078 ixamin^d=ixcmin^d+ix^d;
1079 ixamax^d=ixcmax^d+ix^d;
1080 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)+mf(ixa^s,1:ndim)
1082 bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)*0.5d0**ndim
1086 ixbmax^d=ixmax^d-kr(idims,^d);
1087 call gradientc(te,ixi^l,ixb^l,idims,minq)
1088 gradt(ixb^s,idims)=minq(ixb^s)
1090 if(fl%tc_constant)
then
1091 ka(ixc^s)=fl%tc_k_para
1096 {
do ix^db=ixmin^db,ixmax^db\}
1097 if(minq(ix^d) < block%wextra(ix^d,fl%Tcoff_))
then
1098 minq(ix^d)=block%wextra(ix^d,fl%Tcoff_)
1101 minq(ix^s)=fl%tc_k_para*sqrt(minq(ix^s)**5)
1103 minq(ix^s)=fl%tc_k_para*sqrt(te(ix^s)**5)
1107 ixbmin^d=ixcmin^d+ix^d;
1108 ixbmax^d=ixcmax^d+ix^d;
1109 ka(ixc^s)=ka(ixc^s)+minq(ixb^s)
1112 ka(ixc^s)=0.5d0**ndim*ka(ixc^s)
1114 if(fl%tc_slope_limiter==0)
then
1120 if({ ix^d==0 .and. ^d==idims | .or.})
then
1121 ixbmin^d=ixcmin^d+ix^d;
1122 ixbmax^d=ixcmax^d+ix^d;
1123 qd(ixc^s)=qd(ixc^s)+gradt(ixb^s,idims)
1127 qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
1132 qd(ixc^s)=qvec(ixc^s,idims)*bc(ixc^s,idims)+qd(ixc^s)
1136 gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qd(ixc^s)
1141 ixb^l=ixo^l-kr(idims,^d);
1142 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
1144 if({ ix^d==0 .and. ^d==idims | .or.})
then
1145 ixbmin^d=ixamin^d-ix^d;
1146 ixbmax^d=ixamax^d-ix^d;
1147 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
1150 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
1151 if(fl%tc_saturate)
then
1156 if({ ix^d==0 .and. ^d==idims | .or.})
then
1157 ixbmin^d=ixamin^d-ix^d;
1158 ixbmax^d=ixamax^d-ix^d;
1159 bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
1163 bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
1164 ixb^l=ixa^l+kr(idims,^d);
1165 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))
1166 {
do ix^db=ixamin^db,ixamax^db\}
1167 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
1168 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
1177 ixb^l=ixo^l-kr(idims,^d);
1178 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
1180 ixb^l=ixa^l+kr(idims,^d);
1181 bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
1186 if({ ix^d==0 .and. ^d==idims | .or.})
then
1187 ixbmin^d=ixamin^d-ix^d;
1188 ixbmax^d=ixamax^d-ix^d;
1189 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)+bc(ixb^s,1:ndim)
1190 kaf(ixa^s)=kaf(ixa^s)+ka(ixb^s)
1194 bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)*0.5d0**(ndim-1)
1196 kaf(ixa^s)=kaf(ixa^s)*0.5d0**(ndim-1)
1198 minq(ixa^s)=min(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
1199 maxq(ixa^s)=max(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
1204 if({ ix^d==0 .and. ^d==idims | .or.})
then
1205 ixbmin^d=ixcmin^d+ix^d;
1206 ixbmax^d=ixcmax^d+ix^d;
1207 qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
1211 qdd(ixc^s)=qdd(ixc^s)*0.5d0**(ndim-1)
1215 qd(ixc^s)=qdd(ixc^s)
1216 if({ ix^d==0 .and. ^d==idims | .or.})
then
1217 ixbmin^d=ixamin^d-ix^d;
1218 ixbmax^d=ixamax^d-ix^d;
1219 {
do ixa^db=ixamin^db,ixamax^db
1220 ixb^db=ixa^db-ix^db\}
1221 if(qd(ixb^d)<=minq(ixa^d))
then
1222 qd(ixb^d)=minq(ixa^d)
1223 else if(qd(ixb^d)>=maxq(ixa^d))
then
1224 qd(ixb^d)=maxq(ixa^d)
1227 qvec(ixa^s,idims)=qvec(ixa^s,idims)+bc(ixb^s,idims)**2*qd(ixb^s)
1230 qvec(ixa^s,idims)=kaf(ixa^s)*qvec(ixa^s,idims)*0.5d0**(ndim-1)
1233 ixbmax^d=ixamax^d+kr(idims,^d);
1235 if(idir==idims) cycle
1236 qd(ixi^s)=
slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
1237 qd(ixi^s)=
slope_limiter(qd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
1238 qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qd(ixa^s)
1240 if(fl%tc_saturate)
then
1243 ixb^l=ixa^l+kr(idims,^d);
1244 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))
1245 {
do ix^db=ixamin^db,ixamax^db\}
1246 if(dabs(qvec(ix^d,idims))>qd(ix^d))
then
1249 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
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 or RHD and RMHD or twofl (plasma-neutral) module Adaptation of mod_...
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)
Init TC coefficients: HD case.
subroutine, public sts_set_source_tc_ffhd(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 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)
Get the explicit timestep for the TC (hd implementation)
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 coefficients: MHD case.
subroutine set_source_tc_ffhd(ixIL, ixOL, w, x, fl, qvec, rho, Te, alpha)
double precision function, public get_tc_dt_ffhd(w, ixIL, ixOL, dxD, x, fl)
Get the explicut timestep for the TC (ffhd implementation)
double precision function, dimension(ixi^s) slope_limiter(f, ixIL, ixOL, idims, pm, tc_slope_limiter)
subroutine, public tc_get_ffhd_params(fl, read_ffhd_params)
double precision function, public get_tc_dt_mhd(w, ixIL, ixOL, dxD, x, fl)
Get the explicut timestep for the TC (mhd implementation)