24 logical,
public,
protected ::
hd_dust = .false.
54 integer,
public,
protected ::
rho_
57 integer,
allocatable,
public,
protected ::
mom(:)
60 integer,
public,
protected :: ^
c&m^C_
63 integer,
allocatable,
public,
protected ::
tracer(:)
66 integer,
public,
protected ::
e_
69 integer,
public,
protected ::
p_
72 integer,
public,
protected ::
r_e
75 integer,
public,
protected ::
te_
81 double precision,
public ::
hd_gamma = 5.d0/3.0d0
84 double precision :: gamma_1, inv_gamma_1
90 logical,
public,
protected ::
hd_trac = .false.
97 double precision,
public,
protected ::
h_ion_fr=1d0
107 double precision,
public,
protected ::
rr=1d0
147 subroutine hd_read_params(files)
149 character(len=*),
intent(in) :: files(:)
159 do n = 1,
size(files)
160 open(
unitpar, file=trim(files(n)), status=
"old")
161 read(
unitpar, hd_list,
end=111)
165 end subroutine hd_read_params
168 subroutine hd_write_info(fh)
170 integer,
intent(in) :: fh
171 integer,
parameter :: n_par = 1
172 double precision :: values(n_par)
173 character(len=name_len) :: names(n_par)
174 integer,
dimension(MPI_STATUS_SIZE) :: st
177 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
181 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
182 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
183 end subroutine hd_write_info
208 phys_internal_e = .false.
217 if(
mype==0)
write(*,*)
'WARNING: set hd_trac_type=1'
222 if(
mype==0)
write(*,*)
'WARNING: set hd_trac=F when ndim>=2'
230 if(
mype==0)
write(*,*)
'WARNING: set hd_thermal_conduction=F when hd_energy=F'
234 if(
mype==0)
write(*,*)
'WARNING: set hd_radiative_cooling=F when hd_energy=F'
240 if(
mype==0)
write(*,*)
'WARNING: set hd_partial_ionization=F when eq_state_units=F'
245 allocate(start_indices(number_species),stop_indices(number_species))
253 mom(:) = var_set_momentum(
ndir)
258 e_ = var_set_energy()
268 write(*,*)
'Warning: CAK force addition together with FLD radiation'
273 write(*,*)
'Warning: Optically thin cooling together with FLD radiation'
277 call mpistop(
'implicit dust addition not compatible with FLD radiation')
280 call mpistop(
'using FLD implies the use of an energy equation, set hd_energy=T')
283 r_e = var_set_radiation_energy()
293 phys_get_dt => hd_get_dt
294 phys_get_cmax => hd_get_cmax
295 phys_get_tcutoff => hd_get_tcutoff
296 phys_get_cbounds => hd_get_cbounds
297 phys_get_flux => hd_get_flux
298 phys_add_source_geom => hd_add_source_geom
299 phys_add_source => hd_add_source
305 phys_get_v => hd_get_v
306 phys_get_rho => hd_get_rho
307 phys_write_info => hd_write_info
308 phys_handle_small_values => hd_handle_small_values
311 call hd_physical_units()
321 tracer(itr) = var_set_fluxvar(
"trc",
"trp", itr, need_bc=.false.)
326 te_ = var_set_auxvar(
'Te',
'Te')
335 stop_indices(1)=nwflux
347 phys_update_temperature => hd_update_temperature
359 call mpistop(
"thermal conduction needs hd_energy=T")
370 tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
371 tc_fl%get_rho => hd_get_rho
378 call mpistop(
"radiative cooling needs hd_energy=T")
382 rc_fl%get_rho => hd_get_rho
393 phys_te_images => hd_te_images
409 if (.not.
allocated(flux_type))
then
410 allocate(flux_type(
ndir, nw))
411 flux_type = flux_default
412 else if (any(shape(flux_type) /= [
ndir, nw]))
then
413 call mpistop(
"phys_check error: flux_type has wrong shape")
417 allocate(iw_vector(nvector))
418 iw_vector(1) =
mom(1) - 1
425 subroutine hd_te_images
429 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
431 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
433 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
435 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
438 call mpistop(
"Error in synthesize emission: Unknown convert_type")
440 end subroutine hd_te_images
445 subroutine hd_sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
449 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
450 double precision,
intent(in) :: x(ixi^s,1:
ndim)
451 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
452 double precision,
intent(in) :: my_dt
453 logical,
intent(in) :: fix_conserve_at_step
455 end subroutine hd_sts_set_source_tc_hd
457 function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x)
result(dtnew)
463 integer,
intent(in) :: ixi^
l, ixo^
l
464 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
465 double precision,
intent(in) :: w(ixi^s,1:nw)
466 double precision :: dtnew
469 end function hd_get_tc_dt_hd
471 subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
476 integer,
intent(in) :: ixi^
l,ixo^
l
477 double precision,
intent(inout) :: w(ixi^s,1:nw)
478 double precision,
intent(in) :: x(ixi^s,1:
ndim)
479 integer,
intent(in) :: step
482 logical :: flag(ixi^s,1:nw)
483 character(len=140) :: error_msg
487 if(any(flag(ixo^s,
e_)))
then
497 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,
rho_)
499 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
503 end subroutine hd_tc_handle_small_e
506 subroutine tc_params_read_hd(fl)
508 type(tc_fluid),
intent(inout) :: fl
510 logical :: tc_saturate=.false.
511 double precision :: tc_k_para=0d0
513 namelist /tc_list/ tc_saturate, tc_k_para
517 read(
unitpar, tc_list,
end=111)
520 fl%tc_saturate = tc_saturate
521 fl%tc_k_para = tc_k_para
523 end subroutine tc_params_read_hd
525 subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
527 integer,
intent(in) :: ixi^
l, ixo^
l
528 double precision,
intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:
ndim)
529 double precision,
intent(out) :: rho(ixi^s)
531 rho(ixo^s) = w(ixo^s,
rho_)
533 end subroutine hd_get_rho
537 subroutine rc_params_read(fl)
541 type(rc_fluid),
intent(inout) :: fl
544 integer :: ncool = 4000
547 character(len=std_len) :: coolcurve=
'JCcorona'
550 logical :: tfix=.false.
556 logical :: rc_split=.false.
559 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split
563 read(
unitpar, rc_list,
end=111)
568 fl%coolcurve=coolcurve
572 end subroutine rc_params_read
580 use mod_particles,
only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
583 double precision :: a,b,xfrac,yfrac
591 if (
hd_gamma <= 0.0d0)
call mpistop (
"Error: hd_gamma <= 0")
592 if (
hd_adiab < 0.0d0)
call mpistop (
"Error: hd_adiab < 0")
596 call mpistop (
"Error: hd_gamma <= 0 or hd_gamma == 1.0")
606 call mpistop(
'select IMEX scheme for implicit dust update')
614 call mpistop(
'select IMEX scheme for FLD radiation use')
617 call phys_set_mg_bounds()
619 if(.not.
fld_no_mg)
call mpistop(
'multigrid must have BCs for IMEX and FLD radiation use')
622 write(*,*)
'==FLD SETUP======================'
623 write(*,*)
'Using FLD with settings:'
628 write(*,*)
'Using FLD with settings: fld_kappa0=',
fld_kappa0
629 write(*,*)
'Using FLD with settings: fld_opal_table=',
fld_opal_table
631 write(*,*)
'Using FLD with settings: fld_bisect_tol=',
fld_bisect_tol
632 write(*,*)
'Using FLD with settings: fld_diff_tol=',
fld_diff_tol
636 print *,
'NORMALIZED arad_norm=',
arad_norm
637 print *,
'NORMALIZED c_norm=',
c_norm
644 print *,
'physical fld_kappa (in cgs or SI) =',
fld_kappa0
647 if(
fld_gamma/=
hd_gamma)
call mpistop(
"you must set fld_gamma and hd_gamma equal!")
648 write(*,*)
'===FLD SETUP====================='
652 write(*,*)
'====HD run with settings===================='
653 write(*,*)
'Using mod_hd_phys with settings:'
655 write(*,*)
'Dimensionality :',
ndim
656 write(*,*)
'vector components:',
ndir
658 write(*,*)
'number of variables nw=',nw
659 write(*,*)
' start index iwstart=',iwstart
660 write(*,*)
'number of vector variables=',nvector
661 write(*,*)
'number of stagger variables nws=',nws
662 write(*,*)
'number of variables with BCs=',nwgc
663 write(*,*)
'number of vars with fluxes=',nwflux
664 write(*,*)
'number of vars with flux + BC=',nwfluxbc
665 write(*,*)
'number of auxiliary variables=',nwaux
666 write(*,*)
'number of extra vars without flux=',nwextra
667 write(*,*)
'number of extra vars for wextra=',nw_extra
668 write(*,*)
'number of auxiliary I/O variables=',
nwauxio
682 write(*,*)
'*****Using particles: npayload,ngridvars :', npayload,ngridvars
683 write(*,*)
'*****Using particles: nusrpayload :', nusrpayload
684 write(*,*)
'*****Using particles: num_particles :', num_particles
685 write(*,*)
'*****Using particles: physics_type_particles=',physics_type_particles
688 write(*,*)
'number due to phys_wider_stencil=',phys_wider_stencil
689 write(*,*)
'==========================================='
690 print *,
'========EOS and UNITS==========='
696 print *,
'========EOS and UNITS==========='
712 print *,
' compare this to ',mp_si*(1.d0+4.d0*
he_abundance)
714 print *,
' compare this to ',mp_cgs*(1.d0+4.d0*
he_abundance)
718 print *,
' compare this to ',kb_si*(2.d0+3.d0*
he_abundance)
722 print *,
' compare this to ',kb_cgs*(2.d0+3.d0*
he_abundance)
730 print *,
'mass fraction hydrogen X is =',1/a,
' and this equals ', 1.d0/(1.d0+4.d0*
he_abundance)
731 print *,
'mass fraction helium Y is =',yfrac
732 print *,
' check that 1/mu', b/a,
' is equal to 2X+3Y/4=',2.d0*xfrac+3.d0*yfrac/4.d0
735 print *,
'========UNITS==========='
740 subroutine hd_physical_units
742 double precision :: mp,kb,c_lightspeed,xfrac,sigma_telectron
743 double precision :: a,b
750 sigma_telectron=sigma_te_si
756 sigma_telectron=sigma_te_cgs
857 end subroutine hd_physical_units
864 logical,
intent(in) :: primitive
865 integer,
intent(in) :: ixi^
l, ixo^
l
866 double precision,
intent(in) :: w(ixi^s, nw)
867 logical,
intent(inout) :: flag(ixi^s,1:nw)
868 double precision :: tmp(ixi^s)
877 half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_))
895 integer,
intent(in) :: ixi^
l, ixo^
l
896 double precision,
intent(inout) :: w(ixi^s, nw)
897 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
901 {
do ix^db=ixomin^db,ixomax^db\}
904 w(ix^
d,
e_)=w(ix^
d,
e_)*inv_gamma_1+&
912 call dust_to_conserved(ixi^l, ixo^l, w, x)
921 integer,
intent(in) :: ixi^
l, ixo^
l
922 double precision,
intent(inout) :: w(ixi^s, nw)
923 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
925 double precision :: inv_rho
929 call hd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'hd_to_primitive')
932 {
do ix^db=ixomin^db,ixomax^db\}
933 inv_rho = 1.d0/w(ix^
d,
rho_)
946 call dust_to_primitive(ixi^l, ixo^l, w, x)
954 integer,
intent(in) :: ixi^
l, ixo^
l
955 double precision,
intent(inout) :: w(ixi^s, nw)
956 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
959 w(ixo^s,
e_)=w(ixo^s,
e_)+half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
966 integer,
intent(in) :: ixi^
l, ixo^
l
967 double precision,
intent(inout) :: w(ixi^s, nw)
968 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
971 w(ixo^s,
e_)=w(ixo^s,
e_)-half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
976 subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
978 integer,
intent(in) :: ixi^
l, ixo^
l, idim
979 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
980 double precision,
intent(out) :: v(ixi^s)
982 v(ixo^s) = w(ixo^s,
mom(idim)) / w(ixo^s,
rho_)
983 end subroutine hd_get_v_idim
986 subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
989 integer,
intent(in) :: ixi^
l, ixo^
l
990 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:^nd)
991 double precision,
intent(out) :: v(ixi^s,1:
ndir)
996 v(ixo^s,idir) = w(ixo^s,
mom(idir)) / w(ixo^s,
rho_)
999 end subroutine hd_get_v
1002 subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
1007 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1009 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
1010 double precision,
intent(inout) :: cmax(ixi^s)
1020 cmax(ixo^s)=dabs(w(ixo^s,
mom(idim)))+dsqrt(
hd_gamma*cmax(ixo^s)/w(ixo^s,
rho_))
1026 end subroutine hd_get_cmax
1029 subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
1031 integer,
intent(in) :: ixi^
l,ixo^
l
1032 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1034 double precision,
intent(inout) :: w(ixi^s,1:nw)
1035 double precision,
intent(out) :: tco_local, tmax_local
1037 double precision,
parameter :: trac_delta=0.25d0
1038 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
1039 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
1040 integer :: jxo^
l,hxo^
l
1041 integer :: jxp^
l,hxp^
l,ixp^
l
1042 logical :: lrlt(ixi^s)
1046 te(ixi^s)=w(ixi^s,
p_)/(te(ixi^s)*w(ixi^s,
rho_))
1049 tmax_local=maxval(te(ixo^s))
1056 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
1058 where(lts(ixo^s) > trac_delta)
1061 if(any(lrlt(ixo^s)))
then
1062 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
1073 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
1074 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
1075 w(ixo^s,
tcoff_)=te(ixo^s)*&
1076 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
1078 call mpistop(
"hd_trac_type not allowed for 1D simulation")
1081 end subroutine hd_get_tcutoff
1084 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
1089 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1091 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
1093 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
1094 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1096 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
1099 double precision :: wmean(ixi^s,
nw)
1100 double precision,
dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1108 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1109 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1110 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1111 umean(ixo^s)=(wlp(ixo^s,
mom(idim))*tmp1(ixo^s)+wrp(ixo^s,
mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1123 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1124 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1125 (wrp(ixo^s,
mom(idim))-wlp(ixo^s,
mom(idim)))**2
1127 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1128 if(
present(cmin))
then
1129 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1130 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1132 {
do ix^db=ixomin^db,ixomax^db\}
1133 cmin(ix^
d,1)=sign(one,cmin(ix^
d,1))*max(abs(cmin(ix^
d,1)),hspeed(ix^
d,1))
1134 cmax(ix^
d,1)=sign(one,cmax(ix^
d,1))*max(abs(cmax(ix^
d,1)),hspeed(ix^
d,1))
1138 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1142 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1143 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1154 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1155 tmp1(ixo^s)=wmean(ixo^s,
mom(idim))/wmean(ixo^s,
rho_)
1158 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1160 if(
present(cmin))
then
1161 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1162 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1163 if(h_correction)
then
1164 {
do ix^db=ixomin^db,ixomax^db\}
1165 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1166 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1170 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1174 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1187 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1188 if(
present(cmin))
then
1189 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))-csoundl(ixo^s)
1190 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1191 if(h_correction)
then
1192 {
do ix^db=ixomin^db,ixomax^db\}
1193 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1194 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1198 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1201 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1202 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1206 end subroutine hd_get_cbounds
1212 integer,
intent(in) :: ixi^
l, ixo^
l
1213 double precision,
intent(in) :: w(ixi^s,nw)
1214 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1215 double precision,
intent(out) :: csound2(ixi^s)
1228 integer,
intent(in) :: ixi^
l, ixo^
l
1229 double precision,
intent(in) :: w(ixi^s, 1:nw)
1230 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1231 double precision,
intent(out):: pth(ixi^s)
1235 pth(ixo^s) = (
hd_gamma - 1.0d0) * (w(ixo^s,
e_) - &
1246 {
do ix^db= ixo^lim^db\}
1252 {
do ix^db= ixo^lim^db\}
1254 write(*,*)
"Error: small value of gas pressure",pth(ix^
d),&
1255 " encountered when call hd_get_pthermal"
1257 write(*,*)
"Location: ", x(ix^
d,:)
1258 write(*,*)
"Cell number: ", ix^
d
1260 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^
d,iw)
1264 write(*,*)
"Saving status at the previous time step"
1277 integer,
intent(in) :: ixi^
l, ixo^
l
1278 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
1279 double precision,
intent(out):: csound(ixi^s)
1281 double precision :: wprim(ixi^s, nw)
1283 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1295 integer,
intent(in) :: ixi^
l, ixo^
l
1296 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
1297 double precision,
intent(out):: csound(ixi^s)
1300 double precision :: inv_rho
1301 double precision :: prad_tensor(ixi^s, 1:
ndim, 1:
ndim)
1302 double precision :: prad_max(ixi^s)
1306 {
do ix^db=ixomin^db,ixomax^db \}
1307 inv_rho=1.d0/w(ix^
d,
rho_)
1308 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1312 if(minval(csound(ixo^s))<smalldouble)
then
1313 print *,
'issue with squared speed and rad pressure'
1314 print *,minval(csound(ixo^s))
1315 print *,minval(prad_max(ixo^s))
1316 call mpistop(
"negative squared speed in get_csrad2 for dt")
1327 integer,
intent(in) :: ixi^
l, ixo^
l
1328 double precision,
intent(in) :: w(ixi^s, 1:nw)
1329 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1330 double precision,
intent(out):: prad(ixi^s, 1:
ndim, 1:
ndim)
1340 integer,
intent(in) :: ixi^
l, ixo^
l
1341 double precision,
intent(in) :: w(ixi^s, 1:nw)
1342 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1343 double precision,
intent(out):: pth_plus_prad(ixi^s)
1345 double precision :: wprim(ixi^s, 1:nw)
1346 double precision :: prad_tensor(ixi^s, 1:
ndim, 1:
ndim)
1347 double precision :: prad_max(ixi^s)
1350 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1353 {
do ix^
d = ixomin^
d,ixomax^
d\}
1354 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1356 pth_plus_prad(ixo^s) = wprim(ixo^s,
p_) + prad_max(ixo^s)
1364 integer,
intent(in) :: ixi^
l, ixo^
l
1365 double precision,
intent(in) :: w(ixi^s, 1:nw)
1366 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1367 double precision,
intent(out):: trad(ixi^s)
1376 integer,
intent(in) :: ixi^
l, ixo^
l
1377 double precision,
intent(in) :: w(ixi^s, 1:nw)
1378 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1379 double precision,
intent(out):: res(ixi^s)
1381 double precision :: r(ixi^s)
1385 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1389 subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1391 integer,
intent(in) :: ixi^
l, ixo^
l
1392 double precision,
intent(in) :: w(ixi^s, 1:nw)
1393 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1394 double precision,
intent(out):: res(ixi^s)
1396 double precision :: r(ixi^s)
1399 res(ixo^s) = (
hd_gamma - 1.0d0) * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1400 end subroutine hd_get_temperature_from_eint
1405 integer,
intent(in) :: ixi^
l, ixo^
l
1406 double precision,
intent(in) :: w(ixi^s, 1:nw)
1407 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1408 double precision,
intent(out):: res(ixi^s)
1410 double precision :: r(ixi^s)
1413 res(ixo^s)=w(ixo^s,
p_)/(r(ixo^s)*w(ixo^s,
rho_))
1418 subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1422 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1424 double precision,
intent(in) :: wc(ixi^s, 1:nw)
1426 double precision,
intent(in) :: w(ixi^s, 1:nw)
1427 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1428 double precision,
intent(out) :: f(ixi^s, nwflux)
1430 double precision :: pth(ixi^s)
1434 {
do ix^db=ixomin^db,ixomax^db\}
1444 {
do ix^db=ixomin^db,ixomax^db\}
1447 ^
c&f(ix^d,
m^
c_)=w(ix^d,
mom(idim))*wc(ix^d,
m^
c_)\
1448 f(ix^d,
mom(idim))=f(ix^d,
mom(idim))+pth(ix^d)
1453 {
do ix^db=ixomin^db,ixomax^db\}
1455 f(ix^d,
r_e)=w(ix^d,
mom(idim))*wc(ix^d,
r_e)
1465 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1468 end subroutine hd_get_flux
1475 subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1481 integer,
intent(in) :: ixi^
l, ixo^
l
1482 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:
ndim)
1483 double precision,
intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
1484 double precision :: pth(ixi^s),
source(ixi^s), minrho
1485 integer :: iw,idir, h1x^
l{^nooned, h2x^
l}
1486 integer :: mr_,mphi_
1487 integer :: irho, ifluid, n_fluids
1488 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1510 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1514 do ifluid = 0, n_fluids-1
1516 if (ifluid == 0)
then
1540 where (wct(ixo^s, irho) > minrho)
1541 source(ixo^s) =
source(ixo^s) + wct(ixo^s,mphi_)*wprim(ixo^s,mphi_)
1542 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt*
source(ixo^s)/x(ixo^s,
r_)
1545 where (wct(ixo^s, irho) > minrho)
1546 source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s, mr_)
1547 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1551 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1556 call mpistop(
"Dust geom source terms not implemented yet with spherical geometries")
1560 h1x^
l=ixo^
l-
kr(1,^
d); {^nooned h2x^
l=ixo^
l-
kr(2,^
d);}
1562 pth(ixo^s)=wprim(ixo^s,
p_)
1571 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1572 *(
block%surfaceC(ixo^s, 1) -
block%surfaceC(h1x^s, 1)) &
1573 /
block%dvolume(ixo^s)
1577 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
1581 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1582 * (
block%surfaceC(ixo^s, 2) -
block%surfaceC(h2x^s, 2)) &
1583 /
block%dvolume(ixo^s)
1585 source(ixo^s) =
source(ixo^s) + (wprim(ixo^s,
mom(3))**2 * wprim(ixo^s,
rho_)) / tan(x(ixo^s, 2))
1587 source(ixo^s) =
source(ixo^s) - (wprim(ixo^s,
mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)
1588 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1592 source(ixo^s) = -(wprim(ixo^s,
mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)&
1593 - (wprim(ixo^s,
mom(2)) * wprim(ixo^s,
mom(3))) * wprim(ixo^s,
rho_) / tan(x(ixo^s, 2))
1594 w(ixo^s,
mom(3)) = w(ixo^s,
mom(3)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1601 call mpistop(
"Rotating frame not implemented yet with dust")
1607 end subroutine hd_add_source_geom
1610 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1619 integer,
intent(in) :: ixi^
l, ixo^
l
1620 double precision,
intent(in) :: qdt, dtfactor
1621 double precision,
intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:
ndim)
1622 double precision,
intent(inout) :: w(ixi^s, 1:nw)
1623 logical,
intent(in) :: qsourcesplit
1624 logical,
intent(inout) :: active
1626 double precision :: gravity_field(ixi^s, 1:
ndim)
1627 integer :: idust, idim
1635 qsourcesplit,active,
rc_fl)
1654 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s,
dust_rho(idust))
1666 call hd_add_radiation_source(qdt,ixi^
l,ixo^
l,wct,wctprim,w,x,qsourcesplit,active)
1670 if(.not.qsourcesplit)
then
1672 call hd_update_temperature(ixi^
l,ixo^
l,wct,w,x)
1676 end subroutine hd_add_source
1678 subroutine hd_add_radiation_source(qdt,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1684 integer,
intent(in) :: ixi^
l, ixo^
l
1685 double precision,
intent(in) :: qdt, x(ixi^s,1:
ndim)
1686 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
1687 double precision,
intent(inout) :: w(ixi^s,1:nw)
1688 logical,
intent(in) :: qsourcesplit
1689 logical,
intent(inout) :: active
1695 end subroutine hd_add_radiation_source
1697 subroutine hd_get_dt(wprim, ixI^L, ixO^L, dtnew, dx^D, x)
1705 integer,
intent(in) :: ixi^
l, ixo^
l
1706 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
1707 double precision,
intent(in) :: wprim(ixi^s, 1:nw)
1708 double precision,
intent(inout) :: dtnew
1732 end subroutine hd_get_dt
1736 integer,
intent(in) :: ixi^
l, ixo^
l
1737 double precision,
intent(in) :: w(ixi^s, nw)
1738 double precision :: ke(ixo^s)
1739 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1741 if (
present(inv_rho))
then
1742 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * inv_rho
1744 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) / w(ixo^s,
rho_)
1748 function hd_inv_rho(w, ixI^L, ixO^L)
result(inv_rho)
1750 integer,
intent(in) :: ixi^
l, ixo^
l
1751 double precision,
intent(in) :: w(ixi^s, nw)
1752 double precision :: inv_rho(ixo^s)
1755 inv_rho = 1.0d0 / w(ixo^s,
rho_)
1756 end function hd_inv_rho
1758 subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1765 logical,
intent(in) :: primitive
1766 integer,
intent(in) :: ixi^
l,ixo^
l
1767 double precision,
intent(inout) :: w(ixi^s,1:nw)
1768 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1769 character(len=*),
intent(in) :: subname
1772 logical :: flag(ixi^s,1:nw)
1782 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(idir)) = 0.0d0
1804 where(flag(ixo^s,
e_))
1829 -0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_))
1832 +0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_)
1848 if(.not.primitive)
then
1856 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir))/w(ixo^s,
rho_)
1863 end subroutine hd_handle_small_values
1865 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1868 integer,
intent(in) :: ixi^
l, ixo^
l
1869 double precision,
intent(in) :: w(ixi^s,1:nw)
1870 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1871 double precision,
intent(out):: rfactor(ixi^s)
1873 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1877 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
1879 end subroutine rfactor_from_temperature_ionization
1881 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1883 integer,
intent(in) :: ixi^
l, ixo^
l
1884 double precision,
intent(in) :: w(ixi^s,1:nw)
1885 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1886 double precision,
intent(out):: rfactor(ixi^s)
1890 end subroutine rfactor_from_constant_ionization
1892 subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1896 integer,
intent(in) :: ixi^
l, ixo^
l
1897 double precision,
intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:
ndim)
1898 double precision,
intent(inout) :: w(ixi^s,1:nw)
1900 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1909 end subroutine hd_update_temperature
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_init(phys_gamma)
Initialize the module.
subroutine cak_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
Check time step for total radiation contribution.
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_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_evaluate_implicit(qtc, psa)
inplace update of psa==>F_im(psa)
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_to_conserved(ixil, ixol, w, x)
integer, public, protected dust_n_species
The number of dust species.
subroutine, public dust_get_flux_prim(w, x, ixil, ixol, idim, f)
subroutine, public dust_check_w(ixil, ixol, w, flag)
integer, dimension(:), allocatable, public, protected dust_rho
Indices of the dust densities.
subroutine, public dust_get_cmax(w, x, ixil, ixol, idim, cmax, cmin)
subroutine, public dust_check_params()
subroutine, public dust_get_cmax_prim(w, x, ixil, ixol, idim, cmax, cmin)
subroutine, public dust_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
Get dt related to dust and gas stopping time (Laibe 2011)
subroutine, public dust_init(g_rho, g_mom, g_energy)
subroutine, public dust_implicit_update(dtfactor, qdt, qtc, psb, psa)
Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
Module for flux conservation near refinement boundaries.
Module for flux limited diffusion (FLD)-approximation in Radiation-(Magneto)hydrodynamics simulations...
subroutine, public fld_get_radpress(w, x, ixil, ixol, rad_pressure)
Returns Radiation Pressure as tensor NOTE: w is primitive on entry.
double precision, public fld_bisect_tol
Tolerance for bisection method for Energy sourceterms This is a percentage of the minimum of gas- and...
double precision, public fld_diff_tol
Tolerance for radiative Energy diffusion.
double precision, public fld_gamma
A copy of (m)hd_gamma.
character(len=40) fld_fluxlimiter
flux limiter choice
character(len=40) fld_opal_table
double precision, public fld_kappa0
Opacity value when using constant opacity.
character(len=40) fld_opacity_law
switches for opacity
character(len=40) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
subroutine, public add_fld_rad_force(qdt, ixil, ixol, wct, wctprim, w, x, 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...
logical fld_radforce_split
source split for energy interact and radforce:
subroutine, public fld_init(r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
subroutine, public fld_radforce_get_dt(w, ixil, ixol, dtnew, dxd, x)
get dt limit for radiation force and FLD explicit source additions NOTE: w is primitive on entry
integer nth_for_diff_mg
diffusion coefficient stencil control
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 const_kappae
double precision arad_norm
Normalised radiation constant.
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
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.
double precision const_rad_a
Physical factors useful for radiation fld.
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.
logical slab
Cartesian geometry or not.
double precision const_sigmasb
integer nwauxio
Number of auxiliary variables that are only included in output.
double precision unit_velocity
Physical scaling factor for velocity.
double precision small_r_e
double precision c_norm
Normalised speed of light.
double precision unit_temperature
Physical scaling factor for temperature.
double precision unit_radflux
Physical scaling factor for radiation flux.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
integer nghostcells
Number of ghost cells surrounding a grid.
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.
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 unit_erad
Physical scaling factor for radiation energy density.
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(wprim, ixil, ixol, dtnew, dxd, x)
subroutine gravity_init()
Initialize the module.
subroutine gravity_add_source(qdt, ixil, ixol, wct, wctprim, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
Hydrodynamics physics module.
integer, public, protected m
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, public hd_get_pthermal_plus_pradiation(w, x, ixil, ixol, pth_plus_prad)
calculates the sum of the gas pressure and max Prad tensor element NOTE: only for diagnostic purposes...
subroutine, public hd_ei_to_e(ixil, ixol, w, x)
Transform internal energy to total energy.
subroutine, public hd_get_temperature_from_prim(w, x, ixil, ixol, res)
Calculate temperature=p/rho when in we work in primitives (hence e_ is p_)
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
double precision, public hd_gamma
The adiabatic index.
subroutine, public hd_get_temperature_from_etot(w, x, ixil, ixol, res)
Calculate temperature=p/rho when in e_ the total energy is stored.
integer, public, protected hd_trac_type
logical, public, protected hd_particles
Whether particles module is added.
logical, public, protected hd_radiation_fld
Whether radiation-gas interaction is handled using flux limited diffusion.
type(tc_fluid), allocatable, public tc_fl
subroutine, public hd_check_w(primitive, ixil, ixol, w, flag)
Returns logical argument flag where values are ok.
logical, public, protected hd_viscosity
Whether viscosity is added.
integer, public, protected r_e
Index of the radiation energy (when fld active)
procedure(sub_get_pthermal), pointer, public hd_get_rfactor
integer, public, protected c
Indices of the momentum density for the form of better vectorization.
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.
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
integer, public, protected te_
Indices of temperature.
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
subroutine, public hd_get_pradiation_from_prim(w, x, ixil, ixol, prad)
Calculate radiation pressure within ixO^L NOTE: w is primitive on entry here! NOTE: used in FLD modul...
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
double precision function, dimension(ixo^s), public hd_kin_en(w, ixil, ixol, inv_rho)
subroutine, public hd_to_conserved(ixil, ixol, w, x)
Transform primitive variables into conservative ones.
logical, public, protected hd_cak_force
Whether CAK radiation line force is activated.
subroutine, public hd_phys_init()
Initialize the module.
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
subroutine, public hd_get_csrad2(w, x, ixil, ixol, csound)
Calculate modified squared sound speed for FLD NOTE: only for diagnostic purposes,...
logical, public, protected hd_thermal_conduction
Whether thermal conduction is added.
logical, public, protected eq_state_units
double precision, public hd_adiab
The adiabatic constant.
subroutine, public hd_get_trad(w, x, ixil, ixol, trad)
Calculates radiation temperature.
subroutine, public hd_to_primitive(ixil, ixol, w, x)
Transform conservative variables into primitive ones.
subroutine, public hd_get_csrad2_prim(w, x, ixil, ixol, csound)
Calculate modified squared sound speed for FLD NOTE: w is primitive on entry here!...
integer, public, protected rho_
Index of the density (in the w array)
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.
integer, public, protected c_
type(rc_fluid), allocatable, public rc_fl
logical, public, protected hd_dust_implicit
Whether dust is added using and implicit update in IMEX.
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
logical, public, protected hd_rotating_frame
Whether rotating frame is activated.
subroutine, public hd_get_pthermal(w, x, ixil, ixol, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
subroutine, public hd_e_to_ei(ixil, ixol, w, x)
Transform total energy to internal energy.
module ionization degree - get ionization degree for given temperature
subroutine ionization_degree_from_temperature(ixil, ixol, te, iz_h, iz_he)
subroutine ionization_degree_init()
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
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine radiative_cooling_init(fl, read_params)
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.
subroutine, public small_values_average(ixil, ixol, w, x, w_flag, windex)
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
subroutine, public small_values_error(wprim, x, ixil, ixol, w_flag, subname)
logical, dimension(:), allocatable, public small_values_fix_iw
Whether to apply small value fixes to certain variables.
character(len=20), public small_values_method
How to handle small values.
Generic supertimestepping method which can be used for multiple source terms in the governing equatio...
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 tc_get_hd_params(fl, read_hd_params)
Init TC coefficients: HD case.
double precision function, public get_tc_dt_hd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (hd implementation) Note: also used in 1D MHD (or for neutrals i...
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_hd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, 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
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
The module add viscous source terms and check time step.
subroutine, public viscosity_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
procedure(sub_add_source), pointer, public viscosity_add_source
subroutine, public viscosity_init(phys_wider_stencil)
Initialize the module.