33 logical,
public,
protected ::
rhd_dust = .false.
51 integer,
public,
protected ::
rho_
54 integer,
allocatable,
public,
protected ::
mom(:)
57 integer,
allocatable,
public,
protected ::
tracer(:)
60 integer,
public,
protected ::
e_
63 integer,
public,
protected ::
p_
66 integer,
public,
protected ::
r_e
69 integer,
public,
protected ::
te_
81 double precision,
protected :: small_e
84 double precision,
public,
protected ::
small_r_e = 0.d0
87 logical,
public,
protected ::
rhd_trac = .false.
118 logical,
protected :: radio_acoustic_filter = .false.
119 integer,
protected :: size_ra_filter = 1
125 logical :: dt_c = .false.
129 double precision,
public,
protected ::
h_ion_fr=1d0
139 double precision,
public,
protected ::
rr=1d0
164 subroutine rhd_read_params(files)
166 character(len=*),
intent(in) :: files(:)
177 do n = 1,
size(files)
178 open(
unitpar, file=trim(files(n)), status=
"old")
179 read(
unitpar, rhd_list,
end=111)
183 end subroutine rhd_read_params
188 integer,
intent(in) :: fh
189 integer,
parameter :: n_par = 1
190 double precision :: values(n_par)
191 character(len=name_len) :: names(n_par)
192 integer,
dimension(MPI_STATUS_SIZE) :: st
195 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
199 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
200 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
234 if(
mype==0)
write(*,*)
'WARNING: set rhd_trac_type=1'
239 if(
mype==0)
write(*,*)
'WARNING: set rhd_trac=F when ndim>=2'
247 if(
mype==0)
write(*,*)
'WARNING: set rhd_thermal_conduction=F when rhd_energy=F'
251 if(
mype==0)
write(*,*)
'WARNING: set rhd_radiative_cooling=F when rhd_energy=F'
257 if(
mype==0)
write(*,*)
'WARNING: set rhd_partial_ionization=F when eq_state_units=F'
262 allocate(start_indices(number_species),stop_indices(number_species))
270 mom(:) = var_set_momentum(
ndir)
274 e_ = var_set_energy()
282 r_e = var_set_radiation_energy()
286 te_ = var_set_auxvar(
'Te',
'Te')
313 phys_req_diagonal = .false.
329 call mpistop(
'Radiation formalism unknown')
335 phys_req_diagonal = .true.
342 tracer(itr) = var_set_fluxvar(
"trc",
"trp", itr, need_bc=.false.)
349 stop_indices(1)=nwflux
371 call mpistop(
"thermal conduction needs rhd_energy=T")
372 phys_req_diagonal = .true.
391 call mpistop(
"radiative cooling needs rhd_energy=T")
397 rc_fl%get_var_Rfactor => rhd_get_rfactor
404 te_fl_rhd%get_var_Rfactor => rhd_get_rfactor
420 phys_req_diagonal = .true.
424 if (.not.
allocated(flux_type))
then
425 allocate(flux_type(
ndir, nw))
426 flux_type = flux_default
427 else if (any(shape(flux_type) /= [
ndir, nw]))
then
428 call mpistop(
"phys_check error: flux_type has wrong shape")
432 allocate(iw_vector(nvector))
433 iw_vector(1) =
mom(1) - 1
447 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
449 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
451 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
454 call mpistop(
"Error in synthesize emission: Unknown convert_type")
465 integer,
intent(in) :: ixI^L, ixO^L, igrid, nflux
466 double precision,
intent(in) :: x(ixI^S,1:ndim)
467 double precision,
intent(inout) :: wres(ixI^S,1:nw), w(ixI^S,1:nw)
468 double precision,
intent(in) :: my_dt
469 logical,
intent(in) :: fix_conserve_at_step
481 integer,
intent(in) :: ixi^
l, ixo^
l
482 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
483 double precision,
intent(in) :: w(ixi^s,1:nw)
484 double precision :: dtnew
495 integer,
intent(in) :: ixI^L,ixO^L
496 double precision,
intent(inout) :: w(ixI^S,1:nw)
497 double precision,
intent(in) :: x(ixI^S,1:ndim)
498 integer,
intent(in) :: step
501 logical :: flag(ixI^S,1:nw)
502 character(len=140) :: error_msg
505 where(w(ixo^s,
e_)<small_e) flag(ixo^s,
e_)=.true.
506 if(any(flag(ixo^s,
e_)))
then
509 where(flag(ixo^s,
e_)) w(ixo^s,
e_)=small_e
516 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,
rho_)
518 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
528 type(tc_fluid),
intent(inout) :: fl
530 logical :: tc_saturate=.false.
531 double precision :: tc_k_para=0d0
533 namelist /tc_list/ tc_saturate, tc_k_para
535 do n = 1,
size(par_files)
536 open(
unitpar, file=trim(par_files(n)), status=
"old")
537 read(
unitpar, tc_list,
end=111)
540 fl%tc_saturate = tc_saturate
541 fl%tc_k_para = tc_k_para
547 integer,
intent(in) :: ixI^L, ixO^L
548 double precision,
intent(in) :: w(ixI^S,1:nw),x(ixI^S,1:ndim)
549 double precision,
intent(out) :: rho(ixI^S)
551 rho(ixo^s) = w(ixo^s,
rho_)
561 type(rc_fluid),
intent(inout) :: fl
564 integer :: ncool = 4000
565 double precision :: cfrac=0.1d0
568 character(len=std_len) :: coolcurve=
'JCcorona'
571 character(len=std_len) :: coolmethod=
'exact'
574 logical :: Tfix=.false.
580 logical :: rc_split=.false.
583 namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
587 read(
unitpar, rc_list,
end=111)
592 fl%coolcurve=coolcurve
593 fl%coolmethod=coolmethod
606 if (
rhd_gamma <= 0.0d0)
call mpistop (
"Error: rhd_gamma <= 0")
607 if (
rhd_adiab < 0.0d0)
call mpistop (
"Error: rhd_adiab < 0")
613 call mpistop (
"Error: rhd_gamma <= 0 or rhd_gamma == 1.0")
623 call mpistop(
"Use an IMEX scheme when doing FLD")
642 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
643 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
646 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
647 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
651 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
652 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
660 call mpistop(
"divE_multigrid warning: unknown b.c. ")
667 double precision :: mp,kB
668 double precision :: a,b
725 logical,
intent(in) :: primitive
726 integer,
intent(in) :: ixi^
l, ixo^
l
727 double precision,
intent(in) :: w(ixi^s, nw)
728 logical,
intent(inout) :: flag(ixi^s,1:nw)
729 double precision :: tmp(ixi^s)
755 integer,
intent(in) :: ixi^
l, ixo^
l
756 double precision,
intent(inout) :: w(ixi^s, nw)
757 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
758 double precision :: invgam
768 w(ixo^s,
e_) = w(ixo^s,
e_) * invgam + &
769 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * w(ixo^s,
rho_)
774 w(ixo^s,
mom(idir)) = w(ixo^s,
rho_) * w(ixo^s,
mom(idir))
787 integer,
intent(in) :: ixi^
l, ixo^
l
788 double precision,
intent(inout) :: w(ixi^s, nw)
789 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
791 double precision :: inv_rho(ixo^s)
797 inv_rho = 1.0d0 / w(ixo^s,
rho_)
807 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir)) * inv_rho
820 integer,
intent(in) :: ixI^L, ixO^L
821 double precision,
intent(inout) :: w(ixI^S, nw)
822 double precision,
intent(in) :: x(ixI^S, 1:ndim)
825 w(ixo^s,
e_)=w(ixo^s,
e_)&
833 integer,
intent(in) :: ixI^L, ixO^L
834 double precision,
intent(inout) :: w(ixI^S, nw)
835 double precision,
intent(in) :: x(ixI^S, 1:ndim)
838 w(ixo^s,
e_)=w(ixo^s,
e_)&
846 integer,
intent(in) :: ixI^L, ixO^L
847 double precision :: w(ixI^S, nw)
848 double precision,
intent(in) :: x(ixI^S, 1:ndim)
854 call mpistop(
"energy from entropy can not be used with -eos = iso !")
861 integer,
intent(in) :: ixI^L, ixO^L
862 double precision :: w(ixI^S, nw)
863 double precision,
intent(in) :: x(ixI^S, 1:ndim)
869 call mpistop(
"entropy from energy can not be used with -eos = iso !")
876 integer,
intent(in) :: ixI^L, ixO^L, idim
877 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
878 double precision,
intent(out) :: v(ixI^S)
880 v(ixo^s) = w(ixo^s,
mom(idim)) / w(ixo^s,
rho_)
888 integer,
intent(in) :: ixI^L, ixO^L, idim
889 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
890 double precision,
intent(inout) :: cmax(ixI^S)
891 double precision :: csound(ixI^S)
892 double precision :: v(ixI^S)
894 call rhd_get_v(w, x, ixi^l, ixo^l, idim, v)
896 csound(ixo^s) = dsqrt(csound(ixo^s))
898 cmax(ixo^s) = dabs(v(ixo^s))+csound(ixo^s)
908 integer,
intent(in) :: ixI^L, ixO^L
909 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
910 double precision,
intent(inout) :: a2max(ndim)
911 double precision :: a2(ixI^S,ndim,nw)
912 integer :: gxO^L,hxO^L,jxO^L,kxO^L,i,j
917 hxo^l=ixo^l-
kr(i,^
d);
918 gxo^l=hxo^l-
kr(i,^
d);
919 jxo^l=ixo^l+
kr(i,^
d);
920 kxo^l=jxo^l+
kr(i,^
d);
921 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
922 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
923 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/
dxlevel(i)**2
930 integer,
intent(in) :: ixI^L,ixO^L
931 double precision,
intent(in) :: x(ixI^S,1:ndim)
932 double precision,
intent(inout) :: w(ixI^S,1:nw)
933 double precision,
intent(out) :: tco_local, Tmax_local
935 double precision,
parameter :: trac_delta=0.25d0
936 double precision :: tmp1(ixI^S),Te(ixI^S),lts(ixI^S)
937 double precision :: ltr(ixI^S),ltrc,ltrp,tcoff(ixI^S)
938 integer :: jxO^L,hxO^L
939 integer :: jxP^L,hxP^L,ixP^L
940 logical :: lrlt(ixI^S)
946 tmax_local=maxval(te(ixo^s))
953 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
955 where(lts(ixo^s) > trac_delta)
958 if(any(lrlt(ixo^s)))
then
959 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
970 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
971 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
972 w(ixo^s,
tcoff_)=te(ixo^s)*&
973 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
975 call mpistop(
"mrhd_trac_type not allowed for 1D simulation")
981 subroutine rhd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
986 integer,
intent(in) :: ixI^L, ixO^L, idim
988 double precision,
intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
990 double precision,
intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
991 double precision,
intent(in) :: x(ixI^S, 1:ndim)
992 double precision,
intent(inout) :: cmax(ixI^S,1:number_species)
993 double precision,
intent(inout),
optional :: cmin(ixI^S,1:number_species)
994 double precision,
intent(in) :: Hspeed(ixI^S,1:number_species)
996 double precision :: wmean(ixI^S,nw)
997 double precision,
dimension(ixI^S) :: umean, dmean, csoundL, csoundR, tmp1,tmp2,tmp3
1005 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1006 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1007 tmp3(ixo^s)=1.d0/(dsqrt(wlp(ixo^s,
rho_))+dsqrt(wrp(ixo^s,
rho_)))
1008 umean(ixo^s)=(wlp(ixo^s,
mom(idim))*tmp1(ixo^s)+wrp(ixo^s,
mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1024 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1025 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1026 (wrp(ixo^s,
mom(idim))-wlp(ixo^s,
mom(idim)))**2
1028 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1029 if(
present(cmin))
then
1030 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1031 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1033 {
do ix^db=ixomin^db,ixomax^db\}
1034 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1035 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1039 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1043 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1044 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1048 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1049 tmp1(ixo^s)=wmean(ixo^s,
mom(idim))/wmean(ixo^s,
rho_)
1051 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1053 if(
present(cmin))
then
1054 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1055 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1056 if(h_correction)
then
1057 {
do ix^db=ixomin^db,ixomax^db\}
1058 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1059 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1063 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1067 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1084 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1085 if(
present(cmin))
then
1086 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))-csoundl(ixo^s)
1087 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1088 if(h_correction)
then
1089 {
do ix^db=ixomin^db,ixomax^db\}
1090 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1091 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1095 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1098 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1099 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1109 integer,
intent(in) :: ixi^
l, ixo^
l
1110 double precision,
intent(in) :: w(ixi^s,nw)
1111 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1112 double precision,
intent(out) :: csound2(ixi^s)
1115 csound2(ixo^s)=max(
rhd_gamma,4.d0/3.d0)*csound2(ixo^s)/w(ixo^s,
rho_)
1126 integer,
intent(in) :: ixi^
l, ixo^
l
1127 double precision,
intent(in) :: w(ixi^s, 1:nw)
1128 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1129 double precision,
intent(out):: pth(ixi^s)
1133 pth(ixi^s) = (
rhd_gamma - 1.d0) * (w(ixi^s,
e_) - &
1134 0.5d0 * sum(w(ixi^s,
mom(:))**2, dim=
ndim+1) / w(ixi^s,
rho_))
1146 call mpistop(
'rhd_pressure unknown, use Trad or adiabatic')
1154 {
do ix^db= ixo^lim^db\}
1160 {
do ix^db= ixo^lim^db\}
1162 write(*,*)
"Error: small value of gas pressure",pth(ix^
d),&
1163 " encountered when call rhd_get_pthermal"
1165 write(*,*)
"Location: ", x(ix^
d,:)
1166 write(*,*)
"Cell number: ", ix^
d
1168 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^
d,iw)
1172 write(*,*)
"Saving status at the previous time step"
1186 integer,
intent(in) :: ixi^
l, ixo^
l
1187 double precision,
intent(in) :: w(ixi^s, 1:nw)
1188 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1189 double precision,
intent(out):: prad(ixo^s, 1:
ndim, 1:
ndim)
1197 call mpistop(
'Radiation formalism unknown')
1205 integer,
intent(in) :: ixi^
l, ixo^
l
1206 double precision,
intent(in) :: w(ixi^s, 1:nw)
1207 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1208 double precision :: pth(ixi^s)
1209 double precision :: prad_tensor(ixo^s, 1:
ndim, 1:
ndim)
1210 double precision :: prad_max(ixo^s)
1211 double precision,
intent(out):: ptot(ixi^s)
1217 {
do ix^
d = ixomin^
d,ixomax^
d\}
1218 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1222 if (radio_acoustic_filter)
then
1226 ptot(ixo^s) = pth(ixo^s) + prad_max(ixo^s)
1234 integer,
intent(in) :: ixI^L, ixO^L
1235 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1236 double precision,
intent(inout) :: prad_max(ixO^S)
1238 double precision :: tmp_prad(ixI^S)
1239 integer :: ix^D, filter, idim
1241 if (size_ra_filter .lt. 1)
call mpistop(
"ra filter of size < 1 makes no sense")
1242 if (size_ra_filter .gt.
nghostcells)
call mpistop(
"ra filter of size < nghostcells makes no sense")
1244 tmp_prad(ixi^s) = zero
1245 tmp_prad(ixo^s) = prad_max(ixo^s)
1247 do filter = 1,size_ra_filter
1250 {
do ix^d = ixomin^d,ixomax^d\}
1251 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d+filter*
kr(idim,^d)))
1252 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d-filter*
kr(idim,^d)))
1261 integer,
intent(in) :: ixI^L, ixO^L
1262 double precision,
intent(in) :: w(ixI^S, 1:nw)
1263 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1264 double precision,
intent(out):: res(ixI^S)
1266 double precision :: R(ixI^S)
1268 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1270 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1277 integer,
intent(in) :: ixI^L, ixO^L
1278 double precision,
intent(in) :: w(ixI^S, 1:nw)
1279 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1280 double precision,
intent(out):: res(ixI^S)
1281 double precision :: R(ixI^S)
1283 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1284 res(ixo^s) = (
rhd_gamma - 1.0d0) * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1291 integer,
intent(in) :: ixi^
l, ixo^
l
1292 double precision,
intent(in) :: w(ixi^s, 1:nw)
1293 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1294 double precision :: pth(ixi^s)
1295 double precision,
intent(out):: tgas(ixi^s)
1298 tgas(ixi^s) = pth(ixi^s)/w(ixi^s,
rho_)
1307 integer,
intent(in) :: ixi^
l, ixo^
l
1308 double precision,
intent(in) :: w(ixi^s, 1:nw)
1309 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1310 double precision,
intent(out):: trad(ixi^s)
1322 integer,
intent(in) :: ixI^L, ixO^L
1323 double precision,
intent(inout) :: w(ixI^S, nw)
1324 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1327 w(ixo^s,
e_)=w(ixo^s,
e_)&
1336 integer,
intent(in) :: ixI^L, ixO^L
1337 double precision,
intent(inout) :: w(ixI^S, nw)
1338 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1341 w(ixo^s,
e_)=w(ixo^s,
e_)&
1351 integer,
intent(in) :: ixI^L, ixO^L, idim
1352 double precision,
intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
1353 double precision,
intent(out) :: f(ixI^S, nwflux)
1354 double precision :: pth(ixI^S), v(ixI^S),frame_vel(ixI^S)
1355 integer :: idir, itr
1358 call rhd_get_v(w, x, ixi^l, ixo^l, idim, v)
1360 f(ixo^s,
rho_) = v(ixo^s) * w(ixo^s,
rho_)
1364 f(ixo^s,
mom(idir)) = v(ixo^s) * w(ixo^s,
mom(idir))
1367 f(ixo^s,
mom(idim)) = f(ixo^s,
mom(idim)) + pth(ixo^s)
1371 f(ixo^s,
e_) = v(ixo^s) * (w(ixo^s,
e_) + pth(ixo^s))
1375 f(ixo^s,
r_e) = v(ixo^s) * w(ixo^s,
r_e)
1377 f(ixo^s,
r_e) = zero
1381 f(ixo^s,
tracer(itr)) = v(ixo^s) * w(ixo^s,
tracer(itr))
1397 integer,
intent(in) :: ixI^L, ixO^L, idim
1399 double precision,
intent(in) :: wC(ixI^S, 1:nw)
1401 double precision,
intent(in) :: w(ixI^S, 1:nw)
1402 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1403 double precision,
intent(out) :: f(ixI^S, nwflux)
1404 double precision :: pth(ixI^S),frame_vel(ixI^S)
1405 integer :: idir, itr
1408 pth(ixo^s) = w(ixo^s,
p_)
1413 f(ixo^s,
rho_) = w(ixo^s,
mom(idim)) * w(ixo^s,
rho_)
1417 f(ixo^s,
mom(idir)) = w(ixo^s,
mom(idim)) * wc(ixo^s,
mom(idir))
1420 f(ixo^s,
mom(idim)) = f(ixo^s,
mom(idim)) + pth(ixo^s)
1424 f(ixo^s,
e_) = w(ixo^s,
mom(idim)) * (wc(ixo^s,
e_) + w(ixo^s,
p_))
1428 f(ixo^s,
r_e) = w(ixo^s,
mom(idim)) * wc(ixo^s,
r_e)
1430 f(ixo^s,
r_e) = zero
1461 integer,
intent(in) :: ixI^L, ixO^L
1462 double precision,
intent(in) :: qdt, dtfactor, x(ixI^S, 1:ndim)
1463 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
1467 double precision :: pth(ixI^S), source(ixI^S), minrho
1468 integer :: iw,idir, h1x^L{^NOONED, h2x^L}
1469 integer :: mr_,mphi_
1470 integer :: irho, ifluid, n_fluids
1471 double precision :: exp_factor(ixI^S), del_exp_factor(ixI^S), exp_factor_primitive(ixI^S)
1485 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1490 call mpistop(
"Diffusion term not implemented yet with cylkindrical geometries")
1493 do ifluid = 0, n_fluids-1
1495 if (ifluid == 0)
then
1511 where (wct(ixo^s, irho) > minrho)
1512 source(ixo^s) =
source(ixo^s) + wct(ixo^s, mphi_)**2 / wct(ixo^s, irho)
1513 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1516 where (wct(ixo^s, irho) > minrho)
1517 source(ixo^s) = -wct(ixo^s, mphi_) * wct(ixo^s, mr_) / wct(ixo^s, irho)
1518 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1522 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1528 call mpistop(
"Diffusion term not implemented yet with spherical geometries")
1532 call mpistop(
"Dust geom source terms not implemented yet with spherical geometries")
1536 h1x^l=ixo^l-
kr(1,^
d); {^nooned h2x^l=ixo^l-
kr(2,^
d);}
1539 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1540 *(
block%surfaceC(ixo^s, 1) -
block%surfaceC(h1x^s, 1)) &
1541 /
block%dvolume(ixo^s)
1547 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
1551 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1552 * (
block%surfaceC(ixo^s, 2) -
block%surfaceC(h2x^s, 2)) &
1553 /
block%dvolume(ixo^s)
1555 source(ixo^s) =
source(ixo^s) + (wct(ixo^s,
mom(3))**2 / wct(ixo^s,
rho_)) / tan(x(ixo^s, 2))
1558 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1562 source(ixo^s) = -(wct(ixo^s,
mom(3)) * wct(ixo^s, mr_)) / wct(ixo^s,
rho_)&
1563 - (wct(ixo^s,
mom(2)) * wct(ixo^s,
mom(3))) / wct(ixo^s,
rho_) / tan(x(ixo^s, 2))
1564 w(ixo^s,
mom(3)) = w(ixo^s,
mom(3)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1573 call mpistop(
"Rotating frame not implemented yet with dust")
1582 subroutine rhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1590 integer,
intent(in) :: ixI^L, ixO^L
1591 double precision,
intent(in) :: qdt,dtfactor
1592 double precision,
intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw),x(ixI^S, 1:ndim)
1593 double precision,
intent(inout) :: w(ixI^S, 1:nw)
1594 logical,
intent(in) :: qsourcesplit
1595 logical,
intent(inout) :: active
1597 double precision :: gravity_field(ixI^S, 1:ndim)
1598 integer :: idust, idim
1606 qsourcesplit,active,
rc_fl)
1621 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1625 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s,
dust_rho(idust))
1635 if(.not.qsourcesplit)
then
1650 integer,
intent(in) :: ixI^L, ixO^L
1651 double precision,
intent(in) :: qdt, x(ixI^S,1:ndim)
1652 double precision,
intent(in) :: wCT(ixI^S,1:nw)
1653 double precision,
intent(inout) :: w(ixI^S,1:nw)
1654 logical,
intent(in) :: qsourcesplit
1655 logical,
intent(inout) :: active
1656 double precision :: cmax(ixI^S)
1694 call mpistop(
'Radiation formalism unknown')
1712 integer,
intent(in) :: ixI^L, ixO^L
1713 double precision,
intent(in) :: dx^D, x(ixI^S, 1:^ND)
1714 double precision,
intent(in) :: w(ixI^S, 1:nw)
1715 double precision,
intent(inout) :: dtnew
1719 if (.not. dt_c)
then
1732 call mpistop(
'Radiation formalism unknown')
1756 integer,
intent(in) :: ixi^l, ixo^l
1757 double precision,
intent(in) :: w(ixi^s, nw)
1758 double precision :: ke(ixo^s)
1759 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1761 if (
present(inv_rho))
then
1762 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * inv_rho
1764 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) / w(ixo^s,
rho_)
1770 integer,
intent(in) :: ixi^l, ixo^l
1771 double precision,
intent(in) :: w(ixi^s, nw)
1772 double precision :: inv_rho(ixo^s)
1775 inv_rho = 1.0d0 / w(ixo^s,
rho_)
1785 logical,
intent(in) :: primitive
1786 integer,
intent(in) :: ixI^L,ixO^L
1787 double precision,
intent(inout) :: w(ixI^S,1:nw)
1788 double precision,
intent(in) :: x(ixI^S,1:ndim)
1789 character(len=*),
intent(in) :: subname
1792 logical :: flag(ixI^S,1:nw)
1794 call rhd_check_w(primitive, ixi^l, ixo^l, w, flag)
1802 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(idir)) = 0.0d0
1824 where(flag(ixo^s,
e_))
1849 -0.5d0*sum(w(ixi^s,
mom(:))**2, dim=ndim+1)/w(ixi^s,
rho_))
1852 +0.5d0*sum(w(ixi^s,
mom(:))**2, dim=ndim+1)/w(ixi^s,
rho_)
1866 if(.not.primitive)
then
1874 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir))/w(ixo^s,
rho_)
1886 integer,
intent(in) :: ixI^L, ixO^L
1887 double precision,
intent(in) :: w(ixI^S,1:nw)
1888 double precision,
intent(in) :: x(ixI^S,1:ndim)
1889 double precision,
intent(out):: Rfactor(ixI^S)
1891 double precision :: iz_H(ixO^S),iz_He(ixO^S)
1895 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
1901 integer,
intent(in) :: ixI^L, ixO^L
1902 double precision,
intent(in) :: w(ixI^S,1:nw)
1903 double precision,
intent(in) :: x(ixI^S,1:ndim)
1904 double precision,
intent(out):: Rfactor(ixI^S)
1914 integer,
intent(in) :: ixI^L, ixO^L
1915 double precision,
intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
1916 double precision,
intent(inout) :: w(ixI^S,1:nw)
1918 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 for including anisotropic flux limited diffusion (AFLD)-approximation in Radiation-hydrodynami...
subroutine, public afld_radforce_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine, public get_afld_energy_interact(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 This subroutine handles th...
subroutine, public afld_init(He_abundance, rhd_radiation_diffusion, afld_gamma)
Initialising FLD-module: Read opacities Initialise Multigrid adimensionalise kappa Add extra variable...
subroutine afld_get_diffcoef_central(w, wCT, x, ixIL, ixOL)
Calculates cell-centered diffusion coefficient to be used in multigrid.
subroutine, public get_afld_rad_force(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 This subroutine handles th...
subroutine, public afld_get_radpress(w, x, ixIL, ixOL, rad_pressure)
Calculate Radiation Pressure Returns Radiation Pressure as tensor.
Module with basic data types used in amrvac.
integer, parameter std_len
Default length for strings.
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.
double precision, parameter const_rad_a
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_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_init(g_rho, g_mom, g_energy)
subroutine, public dust_to_conserved(ixIL, ixOL, w, x)
Module for flux conservation near refinement boundaries.
Nicolas Moens Module for including flux limited diffusion (FLD)-approximation in Radiation-hydrodynam...
double precision, public fld_mu
mean particle mass
subroutine, public fld_radforce_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
character(len=8) fld_diff_scheme
Which method to solve diffusion part.
subroutine, public get_fld_rad_force(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 This subroutine handles th...
subroutine, public fld_get_radpress(w, x, ixIL, ixOL, rad_pressure)
Calculate Radiation Pressure Returns Radiation Pressure as tensor.
subroutine, public fld_init(He_abundance, rhd_radiation_diffusion, rhd_gamma)
Initialising FLD-module: Read opacities Initialise Multigrid adimensionalise kappa Add extra variable...
subroutine, public get_fld_energy_interact(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 This subroutine handles th...
subroutine fld_get_diffcoef_central(w, wCT, x, ixIL, ixOL)
Calculates cell-centered diffusion coefficient to be used in multigrid.
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.
integer, parameter bc_noinflow
double precision small_pressure
double precision unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
double precision unit_opacity
Physical scaling factor for Opacity.
integer, parameter unitpar
file handle for IO
integer, parameter bc_asymm
double precision global_time
The global simulation time.
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.
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
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.
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
double precision unit_velocity
Physical scaling factor for velocity.
double precision unit_temperature
Physical scaling factor for temperature.
double precision unit_radflux
Physical scaling factor for radiation flux.
integer, parameter bc_cont
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
integer nghostcells
Number of ghost cells surrounding a grid.
integer, parameter bc_symm
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.
logical use_multigrid
Use multigrid (only available in 2D and 3D)
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.
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 to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
type(mg_t) mg
Data structure containing the multigrid tree.
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)
Radiation-Hydrodynamics physics module Module aims at solving the Hydrodynamic equations toghether wi...
subroutine rhd_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
Add geometrical source terms to w.
logical, public, protected rhd_radiative_cooling
Whether radiative cooling is added.
integer, public, protected rhd_n_tracer
Number of tracer species.
subroutine, public rhd_to_primitive(ixIL, ixOL, w, x)
Transform conservative variables into primitive ones.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
type(tc_fluid), allocatable, public tc_fl
subroutine rhd_ei_to_e(ixIL, ixOL, w, x)
Transform internal energy to total energy.
logical, public, protected rhd_dust
Whether dust is added.
integer, public, protected rhd_trac_type
integer, public, protected te_
Indices of temperature.
subroutine rhd_ei_to_e1(ixIL, ixOL, w, x)
logical, public, protected rhd_rotating_frame
Whether rotating frame is activated.
double precision function, dimension(ixo^s), public rhd_kin_en(w, ixIL, ixOL, inv_rho)
subroutine rhd_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine, public rhd_get_pradiation(w, x, ixIL, ixOL, prad)
Calculate radiation pressure within ixO^L.
subroutine, public rhd_to_conserved(ixIL, ixOL, w, x)
Transform primitive variables into conservative ones.
logical, public, protected rhd_viscosity
Whether viscosity is added.
double precision, public kbmpmua4
kb/(m_p mu)* 1/a_rad**4,
subroutine rhd_get_temperature_from_eint(w, x, ixIL, ixOL, res)
Calculate temperature=p/rho when in e_ the internal energy is stored.
subroutine, public rhd_get_tgas(w, x, ixIL, ixOL, tgas)
Calculates gas temperature.
subroutine, public rhd_check_params
subroutine rfactor_from_temperature_ionization(w, x, ixIL, ixOL, Rfactor)
subroutine rhd_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
Calculate cmax_idim = csound + abs(v_idim) within ixO^L.
logical, public, protected rhd_trac
Whether TRAC method is used.
subroutine, public rhd_phys_init()
Initialize the module.
character(len=8), public rhd_radiation_formalism
Formalism to treat radiation.
integer, public, protected rho_
Index of the density (in the w array)
integer, public, protected r_e
Index of the radiation energy.
logical, public, protected rhd_radiation_diffusion
Treat radiation energy diffusion.
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
character(len=8), public rhd_pressure
In the case of no rhd_energy, how to compute pressure.
subroutine rhd_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
subroutine rhd_e_to_ei(ixIL, ixOL, w, x)
Transform total energy to internal energy.
subroutine rhd_add_radiation_source(qdt, ixIL, ixOL, wCT, w, x, qsourcesplit, active)
subroutine, public rhd_get_trad(w, x, ixIL, ixOL, trad)
Calculates radiation temperature.
subroutine rhd_get_tcutoff(ixIL, ixOL, w, x, tco_local, Tmax_local)
get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
subroutine rhd_write_info(fh)
Write this module's parameters to a snapsoht.
subroutine rfactor_from_constant_ionization(w, x, ixIL, ixOL, Rfactor)
logical, public, protected rhd_partial_ionization
Whether plasma is partially ionized.
subroutine rhos_to_e(ixIL, ixOL, w, x)
logical, public, protected rhd_force_diagonal
Allows overruling default corner filling (for debug mode, since otherwise corner primitives fail)
subroutine, public rhd_get_csound2(w, x, ixIL, ixOL, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. csound2=gamma*p/rho.
logical, public, protected rhd_radiation_force
Treat radiation fld_Rad_force.
logical, public, protected rhd_energy
Whether an energy equation is used.
subroutine, public rhd_get_ptot(w, x, ixIL, ixOL, ptot)
calculates the sum of the gas pressure and max Prad tensor element
double precision, public rhd_adiab
The adiabatic constant.
subroutine rhd_e_to_ei1(ixIL, ixOL, w, x)
Transform total energy to internal energy.
double precision, public, protected small_r_e
The smallest allowed radiation energy.
subroutine tc_params_read_rhd(fl)
subroutine rhd_get_rho(w, x, ixIL, ixOL, rho)
subroutine rhd_tc_handle_small_e(w, x, ixIL, ixOL, step)
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
subroutine rhd_handle_small_values(primitive, w, x, ixIL, ixOL, subname)
subroutine rhd_get_v(w, x, ixIL, ixOL, idim, v)
Calculate v_i = m_i / rho within ixO^L.
double precision function rhd_get_tc_dt_rhd(w, ixIL, ixOL, dxD, x)
subroutine rhd_get_temperature_from_etot(w, x, ixIL, ixOL, res)
Calculate temperature=p/rho when in e_ the total energy is stored.
subroutine rhd_get_flux(wC, w, x, ixIL, ixOL, idim, f)
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
subroutine rhd_get_cmax(w, x, ixIL, ixOL, idim, cmax)
Calculate cmax_idim = csound + abs(v_idim) within ixO^L.
subroutine rhd_update_temperature(ixIL, ixOL, wCT, w, x)
integer, public, protected e_
Index of the energy density (-1 if not present)
subroutine rhd_radio_acoustic_filter(x, ixIL, ixOL, prad_max)
Filter peaks in cmax due to radiation energy density, used for debugging.
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
subroutine, public rhd_set_mg_bounds
Set the boundaries for the diffusion of E.
subroutine e_to_rhos(ixIL, ixOL, w, x)
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
logical, public, protected rhd_thermal_conduction
Whether thermal conduction is added.
subroutine rhd_physical_units
double precision, public, protected rr
subroutine, public rhd_get_pthermal(w, x, ixIL, ixOL, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
subroutine rhd_sts_set_source_tc_rhd(ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux)
double precision function, dimension(ixo^s) rhd_inv_rho(w, ixIL, ixOL)
logical, public, protected rhd_radiation_advection
Treat radiation advection.
logical, public, protected rhd_particles
Whether particles module is added.
subroutine rhd_get_a2max(w, x, ixIL, ixOL, a2max)
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
subroutine rc_params_read(fl)
type(te_fluid), allocatable, public te_fl_rhd
logical, public, protected rhd_energy_interact
Treat radiation-gas energy interaction.
logical, public, protected eq_state_units
type(rc_fluid), allocatable, public rc_fl
logical, public, protected rhd_gravity
Whether gravity is added.
subroutine rhd_get_flux_cons(w, x, ixIL, ixOL, idim, f)
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
subroutine, public rhd_check_w(primitive, ixIL, ixOL, w, flag)
Returns logical argument flag where values are ok.
subroutine rhd_te_images()
double precision, public rhd_gamma
The adiabatic index.
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)
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(special_mg_bc), pointer usr_special_mg_bc
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)