24 logical,
public,
protected ::
hd_dust = .false.
48 integer,
public,
protected ::
rho_
51 integer,
allocatable,
public,
protected ::
mom(:)
54 integer,
allocatable,
public,
protected ::
tracer(:)
57 integer,
public,
protected ::
e_
60 integer,
public,
protected ::
p_
63 integer,
public,
protected ::
te_
69 double precision,
public ::
hd_gamma = 5.d0/3.0d0
75 double precision,
protected :: small_e
78 logical,
public,
protected ::
hd_trac = .false.
88 double precision,
public,
protected ::
h_ion_fr=1d0
91 double precision,
public,
protected ::
he_ion_fr=1d0
98 double precision,
public,
protected ::
rr=1d0
118 subroutine hd_read_params(files)
120 character(len=*),
intent(in) :: files(:)
129 do n = 1,
size(files)
130 open(
unitpar, file=trim(files(n)), status=
"old")
131 read(
unitpar, hd_list,
end=111)
135 end subroutine hd_read_params
140 integer,
intent(in) :: fh
141 integer,
parameter :: n_par = 1
142 double precision :: values(n_par)
143 character(len=name_len) :: names(n_par)
144 integer,
dimension(MPI_STATUS_SIZE) :: st
147 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
151 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
152 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
178 phys_internal_e = .false.
187 if(
mype==0)
write(*,*)
'WARNING: set hd_trac_type=1'
192 if(
mype==0)
write(*,*)
'WARNING: set hd_trac=F when ndim>=2'
200 if(
mype==0)
write(*,*)
'WARNING: set hd_thermal_conduction=F when hd_energy=F'
204 if(
mype==0)
write(*,*)
'WARNING: set hd_radiative_cooling=F when hd_energy=F'
210 if(
mype==0)
write(*,*)
'WARNING: set hd_partial_ionization=F when eq_state_units=F'
215 allocate(start_indices(number_species),stop_indices(number_species))
223 mom(:) = var_set_momentum(
ndir)
227 e_ = var_set_energy()
252 phys_req_diagonal = .false.
264 phys_req_diagonal = .true.
271 tracer(itr) = var_set_fluxvar(
"trc",
"trp", itr, need_bc=.false.)
278 stop_indices(1)=nwflux
282 te_ = var_set_auxvar(
'Te',
'Te')
307 call mpistop(
"thermal conduction needs hd_energy=T")
308 phys_req_diagonal = .true.
327 call mpistop(
"radiative cooling needs hd_energy=T")
333 rc_fl%get_var_Rfactor => hd_get_rfactor
340 te_fl_hd%get_var_Rfactor => hd_get_rfactor
359 phys_req_diagonal = .true.
363 if (.not.
allocated(flux_type))
then
364 allocate(flux_type(
ndir, nw))
365 flux_type = flux_default
366 else if (any(shape(flux_type) /= [
ndir, nw]))
then
367 call mpistop(
"phys_check error: flux_type has wrong shape")
371 allocate(iw_vector(nvector))
372 iw_vector(1) =
mom(1) - 1
383 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
385 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
387 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
389 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
392 call mpistop(
"Error in synthesize emission: Unknown convert_type")
403 integer,
intent(in) :: ixI^L, ixO^L, igrid, nflux
404 double precision,
intent(in) :: x(ixI^S,1:ndim)
405 double precision,
intent(inout) :: wres(ixI^S,1:nw), w(ixI^S,1:nw)
406 double precision,
intent(in) :: my_dt
407 logical,
intent(in) :: fix_conserve_at_step
418 integer,
intent(in) :: ixi^
l, ixo^
l
419 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
420 double precision,
intent(in) :: w(ixi^s,1:nw)
421 double precision :: dtnew
431 integer,
intent(in) :: ixI^L,ixO^L
432 double precision,
intent(inout) :: w(ixI^S,1:nw)
433 double precision,
intent(in) :: x(ixI^S,1:ndim)
434 integer,
intent(in) :: step
437 logical :: flag(ixI^S,1:nw)
438 character(len=140) :: error_msg
441 where(w(ixo^s,
e_)<small_e) flag(ixo^s,
e_)=.true.
442 if(any(flag(ixo^s,
e_)))
then
445 where(flag(ixo^s,
e_)) w(ixo^s,
e_)=small_e
452 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,
rho_)
454 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
463 type(tc_fluid),
intent(inout) :: fl
465 logical :: tc_saturate=.false.
466 double precision :: tc_k_para=0d0
468 namelist /tc_list/ tc_saturate, tc_k_para
472 read(
unitpar, tc_list,
end=111)
475 fl%tc_saturate = tc_saturate
476 fl%tc_k_para = tc_k_para
482 integer,
intent(in) :: ixI^L, ixO^L
483 double precision,
intent(in) :: w(ixI^S,1:nw),x(ixI^S,1:ndim)
484 double precision,
intent(out) :: rho(ixI^S)
486 rho(ixo^s) = w(ixo^s,
rho_)
496 type(rc_fluid),
intent(inout) :: fl
499 integer :: ncool = 4000
500 double precision :: cfrac=0.1d0
503 character(len=std_len) :: coolcurve=
'JCcorona'
506 character(len=std_len) :: coolmethod=
'exact'
509 logical :: Tfix=.false.
515 logical :: rc_split=.false.
518 namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
522 read(
unitpar, rc_list,
end=111)
527 fl%coolcurve=coolcurve
528 fl%coolmethod=coolmethod
541 if (
hd_gamma <= 0.0d0)
call mpistop (
"Error: hd_gamma <= 0")
542 if (
hd_adiab < 0.0d0)
call mpistop (
"Error: hd_adiab < 0")
546 call mpistop (
"Error: hd_gamma <= 0 or hd_gamma == 1.0")
561 double precision :: mp,kB
562 double precision :: a,b
616 logical,
intent(in) :: primitive
617 integer,
intent(in) :: ixi^
l, ixo^
l
618 double precision,
intent(in) :: w(ixi^s, nw)
619 logical,
intent(inout) :: flag(ixi^s,1:nw)
620 double precision :: tmp(ixi^s)
628 tmp(ixo^s) = (
hd_gamma - 1.0d0)*(w(ixo^s,
e_) - &
644 integer,
intent(in) :: ixi^
l, ixo^
l
645 double precision,
intent(inout) :: w(ixi^s, nw)
646 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
647 double precision :: invgam
657 w(ixo^s,
e_) = w(ixo^s,
e_) * invgam + &
658 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * w(ixo^s,
rho_)
663 w(ixo^s,
mom(idir)) = w(ixo^s,
rho_) * w(ixo^s,
mom(idir))
676 integer,
intent(in) :: ixi^
l, ixo^
l
677 double precision,
intent(inout) :: w(ixi^s, nw)
678 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
680 double precision :: inv_rho(ixo^s)
686 inv_rho = 1.0d0 / w(ixo^s,
rho_)
696 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir)) * inv_rho
709 integer,
intent(in) :: ixI^L, ixO^L
710 double precision,
intent(inout) :: w(ixI^S, nw)
711 double precision,
intent(in) :: x(ixI^S, 1:ndim)
714 w(ixo^s,
e_)=w(ixo^s,
e_)&
722 integer,
intent(in) :: ixI^L, ixO^L
723 double precision,
intent(inout) :: w(ixI^S, nw)
724 double precision,
intent(in) :: x(ixI^S, 1:ndim)
727 w(ixo^s,
e_)=w(ixo^s,
e_)&
735 integer,
intent(in) :: ixI^L, ixO^L
736 double precision :: w(ixI^S, nw)
737 double precision,
intent(in) :: x(ixI^S, 1:ndim)
743 call mpistop(
"energy from entropy can not be used with -eos = iso !")
750 integer,
intent(in) :: ixI^L, ixO^L
751 double precision :: w(ixI^S, nw)
752 double precision,
intent(in) :: x(ixI^S, 1:ndim)
758 call mpistop(
"entropy from energy can not be used with -eos = iso !")
765 integer,
intent(in) :: ixI^L, ixO^L, idim
766 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
767 double precision,
intent(out) :: v(ixI^S)
769 v(ixo^s) = w(ixo^s,
mom(idim)) / w(ixo^s,
rho_)
776 integer,
intent(in) :: ixI^L, ixO^L
777 double precision,
intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
778 double precision,
intent(out) :: v(ixI^S,1:ndir)
783 v(ixo^s,idir) = w(ixo^s,
mom(idir)) / w(ixo^s,
rho_)
793 integer,
intent(in) :: ixI^L, ixO^L, idim
794 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
795 double precision,
intent(inout) :: cmax(ixI^S)
796 double precision :: csound(ixI^S)
797 double precision :: v(ixI^S)
801 csound(ixo^s) = dsqrt(csound(ixo^s))
803 cmax(ixo^s) = dabs(v(ixo^s))+csound(ixo^s)
813 integer,
intent(in) :: ixI^L, ixO^L
814 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
815 double precision,
intent(inout) :: a2max(ndim)
816 double precision :: a2(ixI^S,ndim,nw)
817 integer :: gxO^L,hxO^L,jxO^L,kxO^L,i,j
822 hxo^l=ixo^l-
kr(i,^
d);
823 gxo^l=hxo^l-
kr(i,^
d);
824 jxo^l=ixo^l+
kr(i,^
d);
825 kxo^l=jxo^l+
kr(i,^
d);
826 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
827 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
828 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/
dxlevel(i)**2
835 integer,
intent(in) :: ixI^L,ixO^L
836 double precision,
intent(in) :: x(ixI^S,1:ndim)
837 double precision,
intent(inout) :: w(ixI^S,1:nw)
838 double precision,
intent(out) :: tco_local, Tmax_local
840 double precision,
parameter :: trac_delta=0.25d0
841 double precision :: tmp1(ixI^S),Te(ixI^S),lts(ixI^S)
842 double precision :: ltr(ixI^S),ltrc,ltrp,tcoff(ixI^S)
843 integer :: jxO^L,hxO^L
844 integer :: jxP^L,hxP^L,ixP^L
845 logical :: lrlt(ixI^S)
851 tmax_local=maxval(te(ixo^s))
858 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
860 where(lts(ixo^s) > trac_delta)
863 if(any(lrlt(ixo^s)))
then
864 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
875 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
876 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
877 w(ixo^s,
tcoff_)=te(ixo^s)*&
878 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
880 call mpistop(
"mhd_trac_type not allowed for 1D simulation")
886 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
891 integer,
intent(in) :: ixI^L, ixO^L, idim
893 double precision,
intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
895 double precision,
intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
896 double precision,
intent(in) :: x(ixI^S, 1:ndim)
897 double precision,
intent(inout) :: cmax(ixI^S,1:number_species)
898 double precision,
intent(inout),
optional :: cmin(ixI^S,1:number_species)
899 double precision,
intent(in) :: Hspeed(ixI^S,1:number_species)
901 double precision :: wmean(ixI^S,nw)
902 double precision,
dimension(ixI^S) :: umean, dmean, csoundL, csoundR, tmp1,tmp2,tmp3
910 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
911 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
912 tmp3(ixo^s)=1.d0/(dsqrt(wlp(ixo^s,
rho_))+dsqrt(wrp(ixo^s,
rho_)))
913 umean(ixo^s)=(wlp(ixo^s,
mom(idim))*tmp1(ixo^s)+wrp(ixo^s,
mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
923 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
924 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
925 (wrp(ixo^s,
mom(idim))-wlp(ixo^s,
mom(idim)))**2
927 dmean(ixo^s)=dsqrt(dmean(ixo^s))
928 if(
present(cmin))
then
929 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
930 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
932 {
do ix^db=ixomin^db,ixomax^db\}
933 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
934 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
938 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
942 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
943 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
947 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
948 tmp1(ixo^s)=wmean(ixo^s,
mom(idim))/wmean(ixo^s,
rho_)
950 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
952 if(
present(cmin))
then
953 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
954 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
955 if(h_correction)
then
956 {
do ix^db=ixomin^db,ixomax^db\}
957 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
958 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
962 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
966 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
977 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
978 if(
present(cmin))
then
979 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))-csoundl(ixo^s)
980 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
981 if(h_correction)
then
982 {
do ix^db=ixomin^db,ixomax^db\}
983 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
984 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
988 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
991 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
992 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1002 integer,
intent(in) :: ixi^
l, ixo^
l
1003 double precision,
intent(in) :: w(ixi^s,nw)
1004 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1005 double precision,
intent(out) :: csound2(ixi^s)
1018 integer,
intent(in) :: ixi^
l, ixo^
l
1019 double precision,
intent(in) :: w(ixi^s, 1:nw)
1020 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1021 double precision,
intent(out):: pth(ixi^s)
1025 pth(ixo^s) = (
hd_gamma - 1.0d0) * (w(ixo^s,
e_) - &
1036 {
do ix^db= ixo^lim^db\}
1042 {
do ix^db= ixo^lim^db\}
1044 write(*,*)
"Error: small value of gas pressure",pth(ix^
d),&
1045 " encountered when call hd_get_pthermal"
1047 write(*,*)
"Location: ", x(ix^
d,:)
1048 write(*,*)
"Cell number: ", ix^
d
1050 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^
d,iw)
1054 write(*,*)
"Saving status at the previous time step"
1065 integer,
intent(in) :: ixI^L, ixO^L
1066 double precision,
intent(in) :: w(ixI^S, 1:nw)
1067 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1068 double precision,
intent(out):: res(ixI^S)
1070 double precision :: R(ixI^S)
1072 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1074 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1080 integer,
intent(in) :: ixI^L, ixO^L
1081 double precision,
intent(in) :: w(ixI^S, 1:nw)
1082 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1083 double precision,
intent(out):: res(ixI^S)
1085 double precision :: R(ixI^S)
1087 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1088 res(ixo^s) = (
hd_gamma - 1.0d0) * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1096 integer,
intent(in) :: ixI^L, ixO^L, idim
1097 double precision,
intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
1098 double precision,
intent(out) :: f(ixI^S, nwflux)
1099 double precision :: pth(ixI^S), v(ixI^S),frame_vel(ixI^S)
1100 integer :: idir, itr
1105 f(ixo^s,
rho_) = v(ixo^s) * w(ixo^s,
rho_)
1109 f(ixo^s,
mom(idir)) = v(ixo^s) * w(ixo^s,
mom(idir))
1112 f(ixo^s,
mom(idim)) = f(ixo^s,
mom(idim)) + pth(ixo^s)
1116 f(ixo^s,
e_) = v(ixo^s) * (w(ixo^s,
e_) + pth(ixo^s))
1120 f(ixo^s,
tracer(itr)) = v(ixo^s) * w(ixo^s,
tracer(itr))
1136 integer,
intent(in) :: ixI^L, ixO^L, idim
1138 double precision,
intent(in) :: wC(ixI^S, 1:nw)
1140 double precision,
intent(in) :: w(ixI^S, 1:nw)
1141 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1142 double precision,
intent(out) :: f(ixI^S, nwflux)
1143 double precision :: pth(ixI^S),frame_vel(ixI^S)
1144 integer :: idir, itr
1147 pth(ixo^s) = w(ixo^s,
p_)
1152 f(ixo^s,
rho_) = w(ixo^s,
mom(idim)) * w(ixo^s,
rho_)
1156 f(ixo^s,
mom(idir)) = w(ixo^s,
mom(idim)) * wc(ixo^s,
mom(idir))
1159 f(ixo^s,
mom(idim)) = f(ixo^s,
mom(idim)) + pth(ixo^s)
1163 f(ixo^s,
e_) = w(ixo^s,
mom(idim)) * (wc(ixo^s,
e_) + w(ixo^s,
p_))
1196 integer,
intent(in) :: ixI^L, ixO^L
1197 double precision,
intent(in) :: qdt, dtfactor, x(ixI^S, 1:ndim)
1198 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
1202 double precision :: pth(ixI^S), source(ixI^S), minrho
1203 integer :: iw,idir, h1x^L{^NOONED, h2x^L}
1204 integer :: mr_,mphi_
1205 integer :: irho, ifluid, n_fluids
1206 double precision :: exp_factor(ixI^S), del_exp_factor(ixI^S), exp_factor_primitive(ixI^S)
1220 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1224 do ifluid = 0, n_fluids-1
1226 if (ifluid == 0)
then
1242 where (wct(ixo^s, irho) > minrho)
1243 source(ixo^s) =
source(ixo^s) + wct(ixo^s, mphi_)**2 / wct(ixo^s, irho)
1244 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1247 where (wct(ixo^s, irho) > minrho)
1248 source(ixo^s) = -wct(ixo^s, mphi_) * wct(ixo^s, mr_) / wct(ixo^s, irho)
1249 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1253 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1258 call mpistop(
"Dust geom source terms not implemented yet with spherical geometries")
1262 h1x^l=ixo^l-
kr(1,^
d); {^nooned h2x^l=ixo^l-
kr(2,^
d);}
1265 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1266 *(
block%surfaceC(ixo^s, 1) -
block%surfaceC(h1x^s, 1)) &
1267 /
block%dvolume(ixo^s)
1273 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
1277 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1278 * (
block%surfaceC(ixo^s, 2) -
block%surfaceC(h2x^s, 2)) &
1279 /
block%dvolume(ixo^s)
1281 source(ixo^s) =
source(ixo^s) + (wct(ixo^s,
mom(3))**2 / wct(ixo^s,
rho_)) / tan(x(ixo^s, 2))
1284 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1288 source(ixo^s) = -(wct(ixo^s,
mom(3)) * wct(ixo^s, mr_)) / wct(ixo^s,
rho_)&
1289 - (wct(ixo^s,
mom(2)) * wct(ixo^s,
mom(3))) / wct(ixo^s,
rho_) / tan(x(ixo^s, 2))
1290 w(ixo^s,
mom(3)) = w(ixo^s,
mom(3)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1299 call mpistop(
"Rotating frame not implemented yet with dust")
1308 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1317 integer,
intent(in) :: ixI^L, ixO^L
1318 double precision,
intent(in) :: qdt, dtfactor
1319 double precision,
intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
1320 double precision,
intent(inout) :: w(ixI^S, 1:nw)
1321 logical,
intent(in) :: qsourcesplit
1322 logical,
intent(inout) :: active
1324 double precision :: gravity_field(ixI^S, 1:ndim)
1325 integer :: idust, idim
1333 qsourcesplit,active,
rc_fl)
1348 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1352 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s,
dust_rho(idust))
1363 if(.not.qsourcesplit)
then
1379 integer,
intent(in) :: ixI^L, ixO^L
1380 double precision,
intent(in) :: dx^D, x(ixI^S, 1:^ND)
1381 double precision,
intent(in) :: w(ixI^S, 1:nw)
1382 double precision,
intent(inout) :: dtnew
1410 integer,
intent(in) :: ixi^l, ixo^l
1411 double precision,
intent(in) :: w(ixi^s, nw)
1412 double precision :: ke(ixo^s)
1413 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1415 if (
present(inv_rho))
then
1416 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * inv_rho
1418 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) / w(ixo^s,
rho_)
1424 integer,
intent(in) :: ixi^l, ixo^l
1425 double precision,
intent(in) :: w(ixi^s, nw)
1426 double precision :: inv_rho(ixo^s)
1429 inv_rho = 1.0d0 / w(ixo^s,
rho_)
1439 logical,
intent(in) :: primitive
1440 integer,
intent(in) :: ixI^L,ixO^L
1441 double precision,
intent(inout) :: w(ixI^S,1:nw)
1442 double precision,
intent(in) :: x(ixI^S,1:ndim)
1443 character(len=*),
intent(in) :: subname
1446 logical :: flag(ixI^S,1:nw)
1448 call hd_check_w(primitive, ixi^l, ixo^l, w, flag)
1456 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(idir)) = 0.0d0
1464 where(flag(ixo^s,
rho_)) w(ixo^s,
e_) = small_e +
hd_kin_en(w,ixi^l,ixo^l)
1473 where(flag(ixo^s,
e_))
1498 -0.5d0*sum(w(ixi^s,
mom(:))**2, dim=ndim+1)/w(ixi^s,
rho_))
1501 +0.5d0*sum(w(ixi^s,
mom(:))**2, dim=ndim+1)/w(ixi^s,
rho_)
1513 if(.not.primitive)
then
1521 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir))/w(ixo^s,
rho_)
1533 integer,
intent(in) :: ixI^L, ixO^L
1534 double precision,
intent(in) :: w(ixI^S,1:nw)
1535 double precision,
intent(in) :: x(ixI^S,1:ndim)
1536 double precision,
intent(out):: Rfactor(ixI^S)
1538 double precision :: iz_H(ixO^S),iz_He(ixO^S)
1542 rfactor(ixo^s)=(1.d0+iz_h(ixo^s)+0.1d0*(1.d0+iz_he(ixo^s)*(1.d0+iz_he(ixo^s))))/2.3d0
1548 integer,
intent(in) :: ixI^L, ixO^L
1549 double precision,
intent(in) :: w(ixI^S,1:nw)
1550 double precision,
intent(in) :: x(ixI^S,1:ndim)
1551 double precision,
intent(out):: Rfactor(ixI^S)
1561 integer,
intent(in) :: ixI^L, ixO^L
1562 double precision,
intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
1563 double precision,
intent(inout) :: w(ixI^S,1:nw)
1565 double precision :: iz_H(ixO^S),iz_He(ixO^S), pth(ixI^S)
Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices iw=iwmin......
Module with basic data types used in amrvac.
integer, parameter std_len
Default length for strings.
Module to include CAK radiation line force in (magneto)hydrodynamic models Computes both the force fr...
subroutine cak_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
Check time step for total radiation contribution.
subroutine cak_init(phys_gamma)
Initialize the module.
subroutine cak_add_source(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter bigdouble
A very large real number.
Module for including dust species, which interact with the gas through a drag force.
subroutine, public dust_check_w(ixIL, ixOL, w, flag)
subroutine, public dust_get_flux(w, x, ixIL, ixOL, idim, f)
subroutine, public dust_implicit_update(dtfactor, qdt, qtC, psb, psa)
Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
subroutine, public dust_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
Get dt related to dust and gas stopping time (Laibe 2011)
subroutine, public dust_add_source(qdt, ixIL, ixOL, wCT, w, x, qsourcesplit, active)
w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
subroutine, public dust_to_primitive(ixIL, ixOL, w, x)
integer, dimension(:, :), allocatable, public, protected dust_mom
Indices of the dust momentum densities.
subroutine, public dust_get_cmax(w, x, ixIL, ixOL, idim, cmax, cmin)
integer, public, protected dust_n_species
The number of dust species.
integer, dimension(:), allocatable, public, protected dust_rho
Indices of the dust densities.
subroutine, public dust_get_flux_prim(w, x, ixIL, ixOL, idim, f)
subroutine, public dust_check_params()
subroutine, public dust_evaluate_implicit(qtC, psa)
inplace update of psa==>F_im(psa)
subroutine, public dust_init(g_rho, g_mom, g_energy)
subroutine, public dust_to_conserved(ixIL, ixOL, w, x)
Module for flux conservation near refinement boundaries.
Module with geometry-related routines (e.g., divergence, curl)
integer, parameter spherical
integer, parameter cylindrical
integer, parameter cartesian_expansion
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.
double precision small_pressure
double precision unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
double precision global_time
The global simulation time.
double precision unit_mass
Physical scaling factor for mass.
logical use_imex_scheme
whether IMEX in use or not
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer it
Number of time steps taken.
double precision unit_numberdensity
Physical scaling factor for number density.
character(len=std_len) convert_type
Which format to use when converting.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
logical use_particles
Use particles module or not.
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_velocity
Physical scaling factor for velocity.
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 phys_trac
Use TRAC (Johnston 2019 ApJL, 873, L22) for MHD or 1D HD.
logical crash
Save a snapshot before crash a run met unphysical values.
double precision small_density
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
double precision, dimension(ndim) dxlevel
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module for including gravity in (magneto)hydrodynamics simulations.
logical grav_split
source split or not
subroutine gravity_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine gravity_add_source(qdt, ixIL, ixOL, wCT, wCTprim, w, x, energy, rhov, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine gravity_init()
Initialize the module.
Hydrodynamics physics module.
subroutine, public hd_check_params
logical, public, protected hd_energy
Whether an energy equation is used.
logical, public, protected hd_dust
Whether dust is added.
subroutine rfactor_from_temperature_ionization(w, x, ixIL, ixOL, Rfactor)
integer, public, protected e_
Index of the energy density (-1 if not present)
logical, public, protected hd_radiative_cooling
Whether radiative cooling is added.
double precision, public, protected rr
subroutine hd_get_flux(wC, w, x, ixIL, ixOL, idim, f)
double precision, public hd_gamma
The adiabatic index.
subroutine rhos_to_e(ixIL, ixOL, w, x)
subroutine hd_get_a2max(w, x, ixIL, ixOL, a2max)
subroutine hd_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
Add geometrical source terms to w.
integer, public, protected hd_trac_type
subroutine hd_get_flux_cons(w, x, ixIL, ixOL, idim, f)
logical, public, protected hd_particles
Whether particles module is added.
type(tc_fluid), allocatable, public tc_fl
logical, public, protected hd_viscosity
Whether viscosity is added.
subroutine rc_params_read(fl)
subroutine hd_get_rho(w, x, ixIL, ixOL, rho)
subroutine, public hd_check_w(primitive, ixIL, ixOL, w, flag)
Returns logical argument flag where values are ok.
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
subroutine hd_e_to_ei(ixIL, ixOL, w, x)
Transform total energy to internal energy.
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
subroutine hd_sts_set_source_tc_hd(ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux)
subroutine hd_ei_to_e(ixIL, ixOL, w, x)
Transform internal energy to total energy.
subroutine hd_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
subroutine hd_write_info(fh)
Write this module's parameters to a snapsoht.
integer, public, protected te_
Indices of temperature.
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
subroutine hd_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
double precision function hd_get_tc_dt_hd(w, ixIL, ixOL, dxD, x)
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
logical, public, protected hd_cak_force
Whether CAK radiation line force is activated.
subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
Calculate cmax_idim = csound + abs(v_idim) within ixO^L.
subroutine, public hd_phys_init()
Initialize the module.
subroutine hd_get_temperature_from_eint(w, x, ixIL, ixOL, res)
Calculate temperature=p/rho when in e_ the internal energy is stored.
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
logical, public, protected hd_thermal_conduction
Whether thermal conduction is added.
logical, public, protected eq_state_units
logical, public, protected hd_force_diagonal
Allows overruling default corner filling (for debug mode, since otherwise corner primitives fail)
subroutine, public hd_get_pthermal(w, x, ixIL, ixOL, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
double precision, public hd_adiab
The adiabatic constant.
subroutine hd_get_tcutoff(ixIL, ixOL, w, x, tco_local, Tmax_local)
get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
integer, public, protected rho_
Index of the density (in the w array)
subroutine tc_params_read_hd(fl)
subroutine hd_get_v_idim(w, x, ixIL, ixOL, idim, v)
Calculate v_i = m_i / rho within ixO^L.
logical, public, protected hd_partial_ionization
Whether plasma is partially ionized.
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
logical, public, protected hd_gravity
Whether gravity is added.
subroutine hd_physical_units
subroutine rfactor_from_constant_ionization(w, x, ixIL, ixOL, Rfactor)
subroutine, public hd_to_conserved(ixIL, ixOL, w, x)
Transform primitive variables into conservative ones.
type(rc_fluid), allocatable, public rc_fl
subroutine, public hd_to_primitive(ixIL, ixOL, w, x)
Transform conservative variables into primitive ones.
logical, public, protected hd_trac
Whether TRAC method is used.
integer, public, protected hd_n_tracer
Number of tracer species.
type(te_fluid), allocatable, public te_fl_hd
subroutine e_to_rhos(ixIL, ixOL, w, x)
subroutine hd_get_cmax(w, x, ixIL, ixOL, idim, cmax)
Calculate cmax_idim = csound + abs(v_idim) within ixO^L.
logical, public, protected hd_rotating_frame
Whether rotating frame is activated.
subroutine hd_get_v(w, x, ixIL, ixOL, v)
Calculate velocity vector v_i = m_i / rho within ixO^L.
double precision function, dimension(ixo^s), public hd_kin_en(w, ixIL, ixOL, inv_rho)
subroutine hd_update_temperature(ixIL, ixOL, wCT, w, x)
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
subroutine hd_tc_handle_small_e(w, x, ixIL, ixOL, step)
subroutine hd_get_temperature_from_etot(w, x, ixIL, ixOL, res)
Calculate temperature=p/rho when in e_ the total energy is stored.
subroutine hd_handle_small_values(primitive, w, x, ixIL, ixOL, subname)
subroutine, public hd_get_csound2(w, x, ixIL, ixOL, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. csound2=gamma*p/rho.
double precision function, dimension(ixo^s) hd_inv_rho(w, ixIL, ixOL)
module ionization degree - get ionization degree for given temperature
subroutine ionization_degree_init()
subroutine ionization_degree_from_temperature(ixIL, ixOL, Te, iz_H, iz_He)
Module containing all the particle routines.
subroutine particles_init()
Initialize particle data and parameters.
This module defines the procedures of a physics module. It contains function pointers for the various...
module radiative cooling – add optically thin radiative cooling for HD and MHD
subroutine radiative_cooling_init(fl, read_params)
subroutine cooling_get_dt(w, ixIL, ixOL, dtnew, dxD, x, fl)
subroutine radiative_cooling_init_params(phys_gamma, He_abund)
Radiative cooling initialization.
subroutine radiative_cooling_add_source(qdt, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active, fl)
Module for including rotating frame in (magneto)hydrodynamics simulations The rotation vector is assu...
subroutine rotating_frame_add_source(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine rotating_frame_init()
Initialize the module.
Module for handling problematic values in simulations, such as negative pressures.
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
subroutine, public small_values_average(ixIL, ixOL, w, x, w_flag, windex)
logical, dimension(:), allocatable, public small_values_fix_iw
Whether to apply small value fixes to certain variables.
subroutine, public small_values_error(wprim, x, ixIL, ixOL, w_flag, subname)
character(len=20), public small_values_method
How to handle small values.
Generic supertimestepping method 1) in amrvac.par in sts_list set the following parameters which have...
subroutine, public add_sts_method(sts_getdt, sts_set_sources, startVar, nflux, startwbc, nwbc, evolve_B)
subroutine which added programatically a term to be calculated using STS Params: sts_getdt function c...
subroutine, public set_conversion_methods_to_head(sts_before_first_cycle, sts_after_last_cycle)
Set the hooks called before the first cycle and after the last cycle in the STS update This method sh...
subroutine, public set_error_handling_to_head(sts_error_handling)
Set the hook of error handling in the STS update. This method is called before updating the BC....
subroutine, public sts_init()
Initialize sts module.
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 tc_init_params(phys_gamma)
double precision function, public get_tc_dt_hd(w, ixIL, ixOL, dxD, x, fl)
subroutine get_euv_image(qunit, fl)
subroutine get_sxr_image(qunit, fl)
subroutine get_euv_spectrum(qunit, fl)
subroutine get_whitelight_image(qunit, fl)
Module with all the methods that users can customize in AMRVAC.
procedure(rfactor), pointer usr_rfactor
procedure(set_surface), pointer usr_set_surface
procedure(phys_gravity), pointer usr_gravity
procedure(hd_pthermal), pointer usr_set_pthermal
The module add viscous source terms and check time step.
subroutine viscosity_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine viscosity_add_source(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
subroutine viscosity_init(phys_wider_stencil, phys_req_diagonal)
Initialize the module.
subroutine visc_add_source_geom(qdt, ixIL, ixOL, wCT, w, x)
subroutine, public visc_get_flux_prim(w, x, ixIL, ixOL, idim, f, energy)