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
374 call mpistop(
"thermal conduction needs rhd_energy=T")
375 phys_req_diagonal = .true.
394 call mpistop(
"radiative cooling needs rhd_energy=T")
400 rc_fl%get_var_Rfactor => rhd_get_rfactor
407 te_fl_rhd%get_var_Rfactor => rhd_get_rfactor
423 phys_req_diagonal = .true.
427 if (.not.
allocated(flux_type))
then
428 allocate(flux_type(
ndir, nw))
429 flux_type = flux_default
430 else if (any(shape(flux_type) /= [
ndir, nw]))
then
431 call mpistop(
"phys_check error: flux_type has wrong shape")
435 allocate(iw_vector(nvector))
436 iw_vector(1) =
mom(1) - 1
450 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
452 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
454 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
457 call mpistop(
"Error in synthesize emission: Unknown convert_type")
468 integer,
intent(in) :: ixI^L, ixO^L, igrid, nflux
469 double precision,
intent(in) :: x(ixI^S,1:ndim)
470 double precision,
intent(inout) :: wres(ixI^S,1:nw), w(ixI^S,1:nw)
471 double precision,
intent(in) :: my_dt
472 logical,
intent(in) :: fix_conserve_at_step
484 integer,
intent(in) :: ixi^
l, ixo^
l
485 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
486 double precision,
intent(in) :: w(ixi^s,1:nw)
487 double precision :: dtnew
498 integer,
intent(in) :: ixI^L,ixO^L
499 double precision,
intent(inout) :: w(ixI^S,1:nw)
500 double precision,
intent(in) :: x(ixI^S,1:ndim)
501 integer,
intent(in) :: step
504 logical :: flag(ixI^S,1:nw)
505 character(len=140) :: error_msg
508 where(w(ixo^s,
e_)<small_e) flag(ixo^s,
e_)=.true.
509 if(any(flag(ixo^s,
e_)))
then
512 where(flag(ixo^s,
e_)) w(ixo^s,
e_)=small_e
519 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,
rho_)
521 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
531 type(tc_fluid),
intent(inout) :: fl
533 logical :: tc_saturate=.false.
534 double precision :: tc_k_para=0d0
536 namelist /tc_list/ tc_saturate, tc_k_para
538 do n = 1,
size(par_files)
539 open(
unitpar, file=trim(par_files(n)), status=
"old")
540 read(
unitpar, tc_list,
end=111)
543 fl%tc_saturate = tc_saturate
544 fl%tc_k_para = tc_k_para
550 integer,
intent(in) :: ixI^L, ixO^L
551 double precision,
intent(in) :: w(ixI^S,1:nw),x(ixI^S,1:ndim)
552 double precision,
intent(out) :: rho(ixI^S)
554 rho(ixo^s) = w(ixo^s,
rho_)
564 type(rc_fluid),
intent(inout) :: fl
567 integer :: ncool = 4000
568 double precision :: cfrac=0.1d0
571 character(len=std_len) :: coolcurve=
'JCcorona'
574 character(len=std_len) :: coolmethod=
'exact'
577 logical :: Tfix=.false.
583 logical :: rc_split=.false.
586 namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
590 read(
unitpar, rc_list,
end=111)
595 fl%coolcurve=coolcurve
596 fl%coolmethod=coolmethod
609 if (
rhd_gamma <= 0.0d0)
call mpistop (
"Error: rhd_gamma <= 0")
610 if (
rhd_adiab < 0.0d0)
call mpistop (
"Error: rhd_adiab < 0")
616 call mpistop (
"Error: rhd_gamma <= 0 or rhd_gamma == 1.0")
626 call mpistop(
"Use an IMEX scheme when doing FLD")
645 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
646 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
649 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
650 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
654 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
655 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
663 call mpistop(
"divE_multigrid warning: unknown b.c. ")
670 double precision :: mp,kB
671 double precision :: a,b
728 logical,
intent(in) :: primitive
729 integer,
intent(in) :: ixi^
l, ixo^
l
730 double precision,
intent(in) :: w(ixi^s, nw)
731 logical,
intent(inout) :: flag(ixi^s,1:nw)
732 double precision :: tmp(ixi^s)
758 integer,
intent(in) :: ixi^
l, ixo^
l
759 double precision,
intent(inout) :: w(ixi^s, nw)
760 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
761 double precision :: invgam
771 w(ixo^s,
e_) = w(ixo^s,
e_) * invgam + &
772 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * w(ixo^s,
rho_)
777 w(ixo^s,
mom(idir)) = w(ixo^s,
rho_) * w(ixo^s,
mom(idir))
790 integer,
intent(in) :: ixi^
l, ixo^
l
791 double precision,
intent(inout) :: w(ixi^s, nw)
792 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
794 double precision :: inv_rho(ixo^s)
800 inv_rho = 1.0d0 / w(ixo^s,
rho_)
810 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir)) * inv_rho
823 integer,
intent(in) :: ixI^L, ixO^L
824 double precision,
intent(inout) :: w(ixI^S, nw)
825 double precision,
intent(in) :: x(ixI^S, 1:ndim)
828 w(ixo^s,
e_)=w(ixo^s,
e_)&
836 integer,
intent(in) :: ixI^L, ixO^L
837 double precision,
intent(inout) :: w(ixI^S, nw)
838 double precision,
intent(in) :: x(ixI^S, 1:ndim)
841 w(ixo^s,
e_)=w(ixo^s,
e_)&
849 integer,
intent(in) :: ixI^L, ixO^L
850 double precision :: w(ixI^S, nw)
851 double precision,
intent(in) :: x(ixI^S, 1:ndim)
857 call mpistop(
"energy from entropy can not be used with -eos = iso !")
864 integer,
intent(in) :: ixI^L, ixO^L
865 double precision :: w(ixI^S, nw)
866 double precision,
intent(in) :: x(ixI^S, 1:ndim)
872 call mpistop(
"entropy from energy can not be used with -eos = iso !")
879 integer,
intent(in) :: ixI^L, ixO^L, idim
880 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
881 double precision,
intent(out) :: v(ixI^S)
883 v(ixo^s) = w(ixo^s,
mom(idim)) / w(ixo^s,
rho_)
892 integer,
intent(in) :: ixI^L, ixO^L, idim
894 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
895 double precision,
intent(inout) :: cmax(ixI^S)
910 call mpistop(
'rhd_pressure unknown, use Trad or adiabatic')
915 cmax(ixo^s)=dabs(w(ixo^s,
mom(idim)))+dsqrt(
rhd_gamma*cmax(ixo^s)/w(ixo^s,
rho_))
926 integer,
intent(in) :: ixI^L, ixO^L
927 double precision,
intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
928 double precision,
intent(inout) :: a2max(ndim)
929 double precision :: a2(ixI^S,ndim,nw)
930 integer :: gxO^L,hxO^L,jxO^L,kxO^L,i,j
935 hxo^l=ixo^l-
kr(i,^
d);
936 gxo^l=hxo^l-
kr(i,^
d);
937 jxo^l=ixo^l+
kr(i,^
d);
938 kxo^l=jxo^l+
kr(i,^
d);
939 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
940 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
941 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/
dxlevel(i)**2
948 integer,
intent(in) :: ixI^L,ixO^L
949 double precision,
intent(in) :: x(ixI^S,1:ndim)
950 double precision,
intent(inout) :: w(ixI^S,1:nw)
951 double precision,
intent(out) :: tco_local, Tmax_local
953 double precision,
parameter :: trac_delta=0.25d0
954 double precision :: tmp1(ixI^S),Te(ixI^S),lts(ixI^S)
955 double precision :: ltr(ixI^S),ltrc,ltrp,tcoff(ixI^S)
956 integer :: jxO^L,hxO^L
957 integer :: jxP^L,hxP^L,ixP^L
958 logical :: lrlt(ixI^S)
964 tmax_local=maxval(te(ixo^s))
971 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
973 where(lts(ixo^s) > trac_delta)
976 if(any(lrlt(ixo^s)))
then
977 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
988 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
989 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
990 w(ixo^s,
tcoff_)=te(ixo^s)*&
991 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
993 call mpistop(
"mrhd_trac_type not allowed for 1D simulation")
999 subroutine rhd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
1004 integer,
intent(in) :: ixI^L, ixO^L, idim
1006 double precision,
intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
1008 double precision,
intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
1009 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1010 double precision,
intent(inout) :: cmax(ixI^S,1:number_species)
1011 double precision,
intent(inout),
optional :: cmin(ixI^S,1:number_species)
1012 double precision,
intent(in) :: Hspeed(ixI^S,1:number_species)
1014 double precision :: wmean(ixI^S,nw)
1015 double precision,
dimension(ixI^S) :: umean, dmean, csoundL, csoundR, tmp1,tmp2,tmp3
1023 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1024 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1025 tmp3(ixo^s)=1.d0/(dsqrt(wlp(ixo^s,
rho_))+dsqrt(wrp(ixo^s,
rho_)))
1026 umean(ixo^s)=(wlp(ixo^s,
mom(idim))*tmp1(ixo^s)+wrp(ixo^s,
mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1042 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1043 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1044 (wrp(ixo^s,
mom(idim))-wlp(ixo^s,
mom(idim)))**2
1046 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1047 if(
present(cmin))
then
1048 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1049 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1051 {
do ix^db=ixomin^db,ixomax^db\}
1052 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1053 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1057 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1061 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1062 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1066 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1067 tmp1(ixo^s)=wmean(ixo^s,
mom(idim))/wmean(ixo^s,
rho_)
1069 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1071 if(
present(cmin))
then
1072 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1073 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1074 if(h_correction)
then
1075 {
do ix^db=ixomin^db,ixomax^db\}
1076 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1077 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1081 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1085 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1102 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1103 if(
present(cmin))
then
1104 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))-csoundl(ixo^s)
1105 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1106 if(h_correction)
then
1107 {
do ix^db=ixomin^db,ixomax^db\}
1108 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1109 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1113 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1116 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1117 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1127 integer,
intent(in) :: ixi^
l, ixo^
l
1128 double precision,
intent(in) :: w(ixi^s,nw)
1129 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1130 double precision,
intent(out) :: csound2(ixi^s)
1133 csound2(ixo^s)=max(
rhd_gamma,4.d0/3.d0)*csound2(ixo^s)/w(ixo^s,
rho_)
1144 integer,
intent(in) :: ixi^
l, ixo^
l
1145 double precision,
intent(in) :: w(ixi^s, 1:nw)
1146 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1147 double precision,
intent(out):: pth(ixi^s)
1151 pth(ixi^s) = (
rhd_gamma - 1.d0) * (w(ixi^s,
e_) - &
1152 0.5d0 * sum(w(ixi^s,
mom(:))**2, dim=
ndim+1) / w(ixi^s,
rho_))
1164 call mpistop(
'rhd_pressure unknown, use Trad or adiabatic')
1172 {
do ix^db= ixo^lim^db\}
1178 {
do ix^db= ixo^lim^db\}
1180 write(*,*)
"Error: small value of gas pressure",pth(ix^
d),&
1181 " encountered when call rhd_get_pthermal"
1183 write(*,*)
"Location: ", x(ix^
d,:)
1184 write(*,*)
"Cell number: ", ix^
d
1186 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^
d,iw)
1190 write(*,*)
"Saving status at the previous time step"
1204 integer,
intent(in) :: ixi^
l, ixo^
l
1205 double precision,
intent(in) :: w(ixi^s, 1:nw)
1206 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1207 double precision,
intent(out):: prad(ixo^s, 1:
ndim, 1:
ndim)
1215 call mpistop(
'Radiation formalism unknown')
1223 integer,
intent(in) :: ixi^
l, ixo^
l
1224 double precision,
intent(in) :: w(ixi^s, 1:nw)
1225 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1226 double precision :: pth(ixi^s)
1227 double precision :: prad_tensor(ixo^s, 1:
ndim, 1:
ndim)
1228 double precision :: prad_max(ixo^s)
1229 double precision,
intent(out):: ptot(ixi^s)
1235 {
do ix^
d = ixomin^
d,ixomax^
d\}
1236 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1240 if (radio_acoustic_filter)
then
1244 ptot(ixo^s) = pth(ixo^s) + prad_max(ixo^s)
1252 integer,
intent(in) :: ixI^L, ixO^L
1253 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1254 double precision,
intent(inout) :: prad_max(ixO^S)
1256 double precision :: tmp_prad(ixI^S)
1257 integer :: ix^D, filter, idim
1259 if (size_ra_filter .lt. 1)
call mpistop(
"ra filter of size < 1 makes no sense")
1260 if (size_ra_filter .gt.
nghostcells)
call mpistop(
"ra filter of size < nghostcells makes no sense")
1262 tmp_prad(ixi^s) = zero
1263 tmp_prad(ixo^s) = prad_max(ixo^s)
1265 do filter = 1,size_ra_filter
1268 {
do ix^d = ixomin^d,ixomax^d\}
1269 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d+filter*
kr(idim,^d)))
1270 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d-filter*
kr(idim,^d)))
1279 integer,
intent(in) :: ixI^L, ixO^L
1280 double precision,
intent(in) :: w(ixI^S, 1:nw)
1281 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1282 double precision,
intent(out):: res(ixI^S)
1284 double precision :: R(ixI^S)
1286 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1288 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1295 integer,
intent(in) :: ixI^L, ixO^L
1296 double precision,
intent(in) :: w(ixI^S, 1:nw)
1297 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1298 double precision,
intent(out):: res(ixI^S)
1299 double precision :: R(ixI^S)
1301 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1302 res(ixo^s) = (
rhd_gamma - 1.0d0) * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1309 integer,
intent(in) :: ixi^
l, ixo^
l
1310 double precision,
intent(in) :: w(ixi^s, 1:nw)
1311 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1312 double precision :: pth(ixi^s)
1313 double precision,
intent(out):: tgas(ixi^s)
1316 tgas(ixi^s) = pth(ixi^s)/w(ixi^s,
rho_)
1325 integer,
intent(in) :: ixi^
l, ixo^
l
1326 double precision,
intent(in) :: w(ixi^s, 1:nw)
1327 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1328 double precision,
intent(out):: trad(ixi^s)
1340 integer,
intent(in) :: ixI^L, ixO^L
1341 double precision,
intent(inout) :: w(ixI^S, nw)
1342 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1345 w(ixo^s,
e_)=w(ixo^s,
e_)&
1354 integer,
intent(in) :: ixI^L, ixO^L
1355 double precision,
intent(inout) :: w(ixI^S, nw)
1356 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1359 w(ixo^s,
e_)=w(ixo^s,
e_)&
1369 integer,
intent(in) :: ixI^L, ixO^L, idim
1370 double precision,
intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
1371 double precision,
intent(out) :: f(ixI^S, nwflux)
1372 double precision :: pth(ixI^S), v(ixI^S),frame_vel(ixI^S)
1373 integer :: idir, itr
1376 call rhd_get_v(w, x, ixi^l, ixo^l, idim, v)
1378 f(ixo^s,
rho_) = v(ixo^s) * w(ixo^s,
rho_)
1382 f(ixo^s,
mom(idir)) = v(ixo^s) * w(ixo^s,
mom(idir))
1385 f(ixo^s,
mom(idim)) = f(ixo^s,
mom(idim)) + pth(ixo^s)
1389 f(ixo^s,
e_) = v(ixo^s) * (w(ixo^s,
e_) + pth(ixo^s))
1393 f(ixo^s,
r_e) = v(ixo^s) * w(ixo^s,
r_e)
1395 f(ixo^s,
r_e) = zero
1399 f(ixo^s,
tracer(itr)) = v(ixo^s) * w(ixo^s,
tracer(itr))
1415 integer,
intent(in) :: ixI^L, ixO^L, idim
1417 double precision,
intent(in) :: wC(ixI^S, 1:nw)
1419 double precision,
intent(in) :: w(ixI^S, 1:nw)
1420 double precision,
intent(in) :: x(ixI^S, 1:ndim)
1421 double precision,
intent(out) :: f(ixI^S, nwflux)
1422 double precision :: pth(ixI^S),frame_vel(ixI^S)
1423 integer :: idir, itr
1426 pth(ixo^s) = w(ixo^s,
p_)
1431 f(ixo^s,
rho_) = w(ixo^s,
mom(idim)) * w(ixo^s,
rho_)
1435 f(ixo^s,
mom(idir)) = w(ixo^s,
mom(idim)) * wc(ixo^s,
mom(idir))
1438 f(ixo^s,
mom(idim)) = f(ixo^s,
mom(idim)) + pth(ixo^s)
1442 f(ixo^s,
e_) = w(ixo^s,
mom(idim)) * (wc(ixo^s,
e_) + w(ixo^s,
p_))
1446 f(ixo^s,
r_e) = w(ixo^s,
mom(idim)) * wc(ixo^s,
r_e)
1448 f(ixo^s,
r_e) = zero
1479 integer,
intent(in) :: ixI^L, ixO^L
1480 double precision,
intent(in) :: qdt, dtfactor, x(ixI^S, 1:ndim)
1481 double precision,
intent(inout) :: wCT(ixI^S, 1:nw), wprim(ixI^S,1:nw), w(ixI^S, 1:nw)
1485 double precision :: pth(ixI^S), source(ixI^S), minrho
1486 integer :: iw,idir, h1x^L{^NOONED, h2x^L}
1487 integer :: mr_,mphi_
1488 integer :: irho, ifluid, n_fluids
1489 double precision :: exp_factor(ixI^S), del_exp_factor(ixI^S), exp_factor_primitive(ixI^S)
1511 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1516 call mpistop(
"Diffusion term not implemented yet with cylkindrical geometries")
1519 do ifluid = 0, n_fluids-1
1521 if (ifluid == 0)
then
1545 source(ixo^s) =
source(ixo^s) + wprim(ixo^s, mphi_)**2 * wprim(ixo^s, irho)
1546 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1548 source(ixo^s) = -wprim(ixo^s, mphi_) * wprim(ixo^s, mr_) * wprim(ixo^s, irho)
1549 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1552 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1558 call mpistop(
"Diffusion term not implemented yet with spherical geometries")
1562 call mpistop(
"Dust geom source terms not implemented yet with spherical geometries")
1566 h1x^l=ixo^l-
kr(1,^
d); {^nooned h2x^l=ixo^l-
kr(2,^
d);}
1569 pth(ixo^s) =wprim(ixo^s,
p_)
1577 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1578 *(
block%surfaceC(ixo^s, 1) -
block%surfaceC(h1x^s, 1)) &
1579 /
block%dvolume(ixo^s)
1583 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
1587 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1588 * (
block%surfaceC(ixo^s, 2) -
block%surfaceC(h2x^s, 2)) &
1589 /
block%dvolume(ixo^s)
1591 source(ixo^s) =
source(ixo^s) + (wprim(ixo^s,
mom(3))**2 * wprim(ixo^s,
rho_)) / tan(x(ixo^s, 2))
1593 source(ixo^s) =
source(ixo^s) - (wprim(ixo^s,
mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)
1594 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1598 source(ixo^s) = -(wprim(ixo^s,
mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)&
1599 - (wprim(ixo^s,
mom(2)) * wprim(ixo^s,
mom(3))) * wprim(ixo^s,
rho_) / tan(x(ixo^s, 2))
1600 w(ixo^s,
mom(3)) = w(ixo^s,
mom(3)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1609 call mpistop(
"Rotating frame not implemented yet with dust")
1618 subroutine rhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1626 integer,
intent(in) :: ixI^L, ixO^L
1627 double precision,
intent(in) :: qdt,dtfactor
1628 double precision,
intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw),x(ixI^S, 1:ndim)
1629 double precision,
intent(inout) :: w(ixI^S, 1:nw)
1630 logical,
intent(in) :: qsourcesplit
1631 logical,
intent(inout) :: active
1633 double precision :: gravity_field(ixI^S, 1:ndim)
1634 integer :: idust, idim
1642 qsourcesplit,active,
rc_fl)
1657 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1661 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s,
dust_rho(idust))
1671 if(.not.qsourcesplit)
then
1686 integer,
intent(in) :: ixI^L, ixO^L
1687 double precision,
intent(in) :: qdt, x(ixI^S,1:ndim)
1688 double precision,
intent(in) :: wCT(ixI^S,1:nw)
1689 double precision,
intent(inout) :: w(ixI^S,1:nw)
1690 logical,
intent(in) :: qsourcesplit
1691 logical,
intent(inout) :: active
1692 double precision :: cmax(ixI^S)
1730 call mpistop(
'Radiation formalism unknown')
1748 integer,
intent(in) :: ixI^L, ixO^L
1749 double precision,
intent(in) :: dx^D, x(ixI^S, 1:^ND)
1750 double precision,
intent(in) :: w(ixI^S, 1:nw)
1751 double precision,
intent(inout) :: dtnew
1755 if (.not. dt_c)
then
1768 call mpistop(
'Radiation formalism unknown')
1792 integer,
intent(in) :: ixi^l, ixo^l
1793 double precision,
intent(in) :: w(ixi^s, nw)
1794 double precision :: ke(ixo^s)
1795 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1797 if (
present(inv_rho))
then
1798 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * inv_rho
1800 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) / w(ixo^s,
rho_)
1806 integer,
intent(in) :: ixi^l, ixo^l
1807 double precision,
intent(in) :: w(ixi^s, nw)
1808 double precision :: inv_rho(ixo^s)
1811 inv_rho = 1.0d0 / w(ixo^s,
rho_)
1821 logical,
intent(in) :: primitive
1822 integer,
intent(in) :: ixI^L,ixO^L
1823 double precision,
intent(inout) :: w(ixI^S,1:nw)
1824 double precision,
intent(in) :: x(ixI^S,1:ndim)
1825 character(len=*),
intent(in) :: subname
1828 logical :: flag(ixI^S,1:nw)
1830 call rhd_check_w(primitive, ixi^l, ixo^l, w, flag)
1838 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(idir)) = 0.0d0
1860 where(flag(ixo^s,
e_))
1885 -0.5d0*sum(w(ixi^s,
mom(:))**2, dim=ndim+1)/w(ixi^s,
rho_))
1888 +0.5d0*sum(w(ixi^s,
mom(:))**2, dim=ndim+1)/w(ixi^s,
rho_)
1902 if(.not.primitive)
then
1910 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir))/w(ixo^s,
rho_)
1922 integer,
intent(in) :: ixI^L, ixO^L
1923 double precision,
intent(in) :: w(ixI^S,1:nw)
1924 double precision,
intent(in) :: x(ixI^S,1:ndim)
1925 double precision,
intent(out):: Rfactor(ixI^S)
1927 double precision :: iz_H(ixO^S),iz_He(ixO^S)
1931 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
1937 integer,
intent(in) :: ixI^L, ixO^L
1938 double precision,
intent(in) :: w(ixI^S,1:nw)
1939 double precision,
intent(in) :: x(ixI^S,1:ndim)
1940 double precision,
intent(out):: Rfactor(ixI^S)
1950 integer,
intent(in) :: ixI^L, ixO^L
1951 double precision,
intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
1952 double precision,
intent(inout) :: w(ixI^S,1:nw)
1954 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_get_cmax_prim(w, x, ixIL, ixOL, idim, cmax, cmin)
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.
double precision, 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 for MHD or 1D HD.
logical fix_small_values
fix small values with average or replace methods
logical crash
Save a snapshot before crash a run met unphysical values.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
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
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...
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+)
subroutine rhd_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, wprim, w, x)
Add geometrical source terms to w.
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 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 tc_init_params(phys_gamma)
double precision function, public get_tc_dt_hd(w, ixIL, ixOL, dxD, x, fl)
Get the explicit timestep for the TC (hd implementation)
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)