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_
78 integer,
public,
protected ::
fip_ = -1
81 logical,
public,
protected ::
hd_fip = .false.
87 double precision,
public ::
hd_gamma = 5.d0/3.0d0
90 double precision :: gamma_1, inv_gamma_1
96 logical,
public,
protected ::
hd_trac = .false.
103 double precision,
public,
protected ::
h_ion_fr=1d0
113 double precision,
public,
protected ::
rr=1d0
153 subroutine hd_read_params(files)
155 character(len=*),
intent(in) :: files(:)
165 do n = 1,
size(files)
166 open(
unitpar, file=trim(files(n)), status=
"old")
167 read(
unitpar, hd_list,
end=111)
171 end subroutine hd_read_params
174 subroutine hd_write_info(fh)
176 integer,
intent(in) :: fh
177 integer,
parameter :: n_par = 1
178 double precision :: values(n_par)
179 character(len=name_len) :: names(n_par)
180 integer,
dimension(MPI_STATUS_SIZE) :: st
183 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
187 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
188 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
189 end subroutine hd_write_info
214 phys_internal_e = .false.
223 if(
mype==0)
write(*,*)
'WARNING: set hd_trac_type=1'
228 if(
mype==0)
write(*,*)
'WARNING: set hd_trac=F when ndim>=2'
236 if(
mype==0)
write(*,*)
'WARNING: set hd_thermal_conduction=F when hd_energy=F'
240 if(
mype==0)
write(*,*)
'WARNING: set hd_radiative_cooling=F when hd_energy=F'
246 if(
mype==0)
write(*,*)
'WARNING: set hd_partial_ionization=F when eq_state_units=F'
251 allocate(start_indices(number_species),stop_indices(number_species))
259 mom(:) = var_set_momentum(
ndir)
264 e_ = var_set_energy()
274 write(*,*)
'Warning: CAK force addition together with FLD radiation'
279 write(*,*)
'Warning: Optically thin cooling together with FLD radiation'
283 call mpistop(
'implicit dust addition not compatible with FLD radiation')
286 call mpistop(
'using FLD implies the use of an energy equation, set hd_energy=T')
289 r_e = var_set_radiation_energy()
299 phys_get_dt => hd_get_dt
300 phys_get_cmax => hd_get_cmax
301 phys_get_tcutoff => hd_get_tcutoff
302 phys_get_cbounds => hd_get_cbounds
303 phys_get_flux => hd_get_flux
304 phys_add_source_geom => hd_add_source_geom
305 phys_add_source => hd_add_source
311 phys_get_v => hd_get_v
312 phys_get_rho => hd_get_rho
313 phys_write_info => hd_write_info
314 phys_handle_small_values => hd_handle_small_values
317 call hd_physical_units()
324 fip_ = var_set_fluxvar(
'rho_fip',
'fip', need_bc=.false.)
333 tracer(itr) = var_set_fluxvar(
"trc",
"trp", itr, need_bc=.false.)
338 te_ = var_set_auxvar(
'Te',
'Te')
347 stop_indices(1)=nwflux
359 phys_update_temperature => hd_update_temperature
371 call mpistop(
"thermal conduction needs hd_energy=T")
382 tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
383 tc_fl%get_rho => hd_get_rho
390 call mpistop(
"radiative cooling needs hd_energy=T")
395 rc_fl%get_rho => hd_get_rho
406 phys_te_images => hd_te_images
422 if (.not.
allocated(flux_type))
then
423 allocate(flux_type(
ndir, nw))
424 flux_type = flux_default
425 else if (any(shape(flux_type) /= [
ndir, nw]))
then
426 call mpistop(
"phys_check error: flux_type has wrong shape")
430 allocate(iw_vector(nvector))
431 iw_vector(1) =
mom(1) - 1
438 subroutine hd_te_images
442 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
444 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
446 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
448 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
451 call mpistop(
"Error in synthesize emission: Unknown convert_type")
453 end subroutine hd_te_images
458 subroutine hd_sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
462 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
463 double precision,
intent(in) :: x(ixi^s,1:
ndim)
464 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
465 double precision,
intent(in) :: my_dt
466 logical,
intent(in) :: fix_conserve_at_step
468 end subroutine hd_sts_set_source_tc_hd
470 function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x)
result(dtnew)
476 integer,
intent(in) :: ixi^
l, ixo^
l
477 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
478 double precision,
intent(in) :: w(ixi^s,1:nw)
479 double precision :: dtnew
482 end function hd_get_tc_dt_hd
484 subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
489 integer,
intent(in) :: ixi^
l,ixo^
l
490 double precision,
intent(inout) :: w(ixi^s,1:nw)
491 double precision,
intent(in) :: x(ixi^s,1:
ndim)
492 integer,
intent(in) :: step
495 logical :: flag(ixi^s,1:nw)
496 character(len=140) :: error_msg
500 if(any(flag(ixo^s,
e_)))
then
510 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,
rho_)
512 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
516 end subroutine hd_tc_handle_small_e
519 subroutine tc_params_read_hd(fl)
521 type(tc_fluid),
intent(inout) :: fl
523 logical :: tc_saturate=.false.
524 double precision :: tc_k_para=0d0
526 namelist /tc_list/ tc_saturate, tc_k_para
530 read(
unitpar, tc_list,
end=111)
533 fl%tc_saturate = tc_saturate
534 fl%tc_k_para = tc_k_para
536 end subroutine tc_params_read_hd
538 subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
540 integer,
intent(in) :: ixi^
l, ixo^
l
541 double precision,
intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:
ndim)
542 double precision,
intent(out) :: rho(ixi^s)
544 rho(ixo^s) = w(ixo^s,
rho_)
546 end subroutine hd_get_rho
550 subroutine rc_params_read(fl)
554 type(rc_fluid),
intent(inout) :: fl
557 integer :: ncool = 4000
560 character(len=std_len) :: coolcurve=
'JCcorona'
563 logical :: tfix=.false.
569 logical :: rc_split=.false.
571 logical :: rad_damp=.false.
572 double precision :: rad_damp_height=0.5d0
573 double precision :: rad_damp_scale=0.15d0
574 logical :: rad_newton = .false.
575 double precision :: rad_newton_trad = 0.006d0
576 double precision :: rad_newton_rhosurf = 1.d4
577 double precision :: rad_newton_pthick = 25.d0
580 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split, &
581 rad_newton, rad_newton_trad, rad_newton_rhosurf, &
582 rad_newton_pthick, rad_damp, rad_damp_height, rad_damp_scale
586 read(
unitpar, rc_list,
end=111)
591 fl%coolcurve=coolcurve
595 fl%rad_damp = rad_damp
596 fl%rad_damp_height = rad_damp_height
597 fl%rad_damp_scale = rad_damp_scale
598 fl%rad_newton = rad_newton
599 fl%rad_newton_trad = rad_newton_trad
600 fl%rad_newton_rhosurf = rad_newton_rhosurf
601 fl%rad_newton_pthick = rad_newton_pthick
602 end subroutine rc_params_read
610 use mod_particles,
only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
613 double precision :: a,b,xfrac,yfrac
621 if (
hd_gamma <= 0.0d0)
call mpistop (
"Error: hd_gamma <= 0")
622 if (
hd_adiab < 0.0d0)
call mpistop (
"Error: hd_adiab < 0")
626 call mpistop (
"Error: hd_gamma <= 0 or hd_gamma == 1.0")
636 call mpistop(
'select IMEX scheme for implicit dust update')
644 call mpistop(
'select IMEX scheme for FLD radiation use')
647 call phys_set_mg_bounds()
649 if(.not.
fld_no_mg)
call mpistop(
'multigrid must have BCs for IMEX and FLD radiation use')
652 write(*,*)
'==FLD SETUP======================'
653 write(*,*)
'Using FLD with settings:'
658 write(*,*)
'Using FLD with settings: fld_kappa0=',
fld_kappa0
659 write(*,*)
'Using FLD with settings: fld_opal_table=',
fld_opal_table
661 write(*,*)
'Using FLD with settings: fld_bisect_tol=',
fld_bisect_tol
662 write(*,*)
'Using FLD with settings: fld_diff_tol=',
fld_diff_tol
666 print *,
'NORMALIZED arad_norm=',
arad_norm
667 print *,
'NORMALIZED c_norm=',
c_norm
674 print *,
'physical fld_kappa (in cgs or SI) =',
fld_kappa0
677 if(
fld_gamma/=
hd_gamma)
call mpistop(
"you must set fld_gamma and hd_gamma equal!")
678 write(*,*)
'===FLD SETUP====================='
682 write(*,*)
'====HD run with settings===================='
683 write(*,*)
'Using mod_hd_phys with settings:'
685 write(*,*)
'Dimensionality :',
ndim
686 write(*,*)
'vector components:',
ndir
688 write(*,*)
'number of variables nw=',nw
689 write(*,*)
' start index iwstart=',iwstart
690 write(*,*)
'number of vector variables=',nvector
691 write(*,*)
'number of stagger variables nws=',nws
692 write(*,*)
'number of variables with BCs=',nwgc
693 write(*,*)
'number of vars with fluxes=',nwflux
694 write(*,*)
'number of vars with flux + BC=',nwfluxbc
695 write(*,*)
'number of auxiliary variables=',nwaux
696 write(*,*)
'number of extra vars without flux=',nwextra
697 write(*,*)
'number of extra vars for wextra=',nw_extra
698 write(*,*)
'number of auxiliary I/O variables=',
nwauxio
712 write(*,*)
'*****Using particles: npayload,ngridvars :', npayload,ngridvars
713 write(*,*)
'*****Using particles: nusrpayload :', nusrpayload
714 write(*,*)
'*****Using particles: num_particles :', num_particles
715 write(*,*)
'*****Using particles: physics_type_particles=',physics_type_particles
718 write(*,*)
'number due to phys_wider_stencil=',phys_wider_stencil
719 write(*,*)
'==========================================='
720 print *,
'========EOS and UNITS==========='
726 print *,
'========EOS and UNITS==========='
742 print *,
' compare this to ',mp_si*(1.d0+4.d0*
he_abundance)
744 print *,
' compare this to ',mp_cgs*(1.d0+4.d0*
he_abundance)
748 print *,
' compare this to ',kb_si*(2.d0+3.d0*
he_abundance)
752 print *,
' compare this to ',kb_cgs*(2.d0+3.d0*
he_abundance)
760 print *,
'mass fraction hydrogen X is =',1/a,
' and this equals ', 1.d0/(1.d0+4.d0*
he_abundance)
761 print *,
'mass fraction helium Y is =',yfrac
762 print *,
' check that 1/mu', b/a,
' is equal to 2X+3Y/4=',2.d0*xfrac+3.d0*yfrac/4.d0
765 print *,
'========UNITS==========='
770 subroutine hd_physical_units
772 double precision :: mp,kb,c_lightspeed,xfrac,sigma_telectron
773 double precision :: a,b
780 sigma_telectron=sigma_te_si
786 sigma_telectron=sigma_te_cgs
887 end subroutine hd_physical_units
894 logical,
intent(in) :: primitive
895 integer,
intent(in) :: ixi^
l, ixo^
l
896 double precision,
intent(in) :: w(ixi^s, nw)
897 logical,
intent(inout) :: flag(ixi^s,1:nw)
898 double precision :: tmp(ixi^s)
907 half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_))
921 subroutine hd_bound_fip(primitive, ixI^L, ixO^L, w)
923 logical,
intent(in) :: primitive
924 integer,
intent(in) :: ixi^
l, ixo^
l
925 double precision,
intent(inout) :: w(ixi^s,1:nw)
927 double precision :: rho_safe(ixi^s), fip_prim(ixi^s)
935 fip_prim(ixo^s) = w(ixo^s,
fip_) / rho_safe(ixo^s)
936 fip_prim(ixo^s) = min(
maxfip, max(
minfip, fip_prim(ixo^s)))
937 w(ixo^s,
fip_) = rho_safe(ixo^s) * fip_prim(ixo^s)
939 end subroutine hd_bound_fip
945 integer,
intent(in) :: ixi^
l, ixo^
l
946 double precision,
intent(inout) :: w(ixi^s, nw)
947 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
951 if (
hd_fip)
call hd_bound_fip(.true., ixi^
l, ixo^
l, w)
953 {
do ix^db=ixomin^db,ixomax^db\}
956 w(ix^
d,
e_)=w(ix^
d,
e_)*inv_gamma_1+&
965 call dust_to_conserved(ixi^l, ixo^l, w, x)
974 integer,
intent(in) :: ixi^
l, ixo^
l
975 double precision,
intent(inout) :: w(ixi^s, nw)
976 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
978 double precision :: inv_rho
982 call hd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'hd_to_primitive')
985 {
do ix^db=ixomin^db,ixomax^db\}
986 inv_rho = 1.d0/w(ix^
d,
rho_)
998 if (
hd_fip)
call hd_bound_fip(.true., ixi^l, ixo^l, w)
1002 call dust_to_primitive(ixi^l, ixo^l, w, x)
1010 integer,
intent(in) :: ixi^
l, ixo^
l
1011 double precision,
intent(inout) :: w(ixi^s, nw)
1012 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1015 w(ixo^s,
e_)=w(ixo^s,
e_)+half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
1022 integer,
intent(in) :: ixi^
l, ixo^
l
1023 double precision,
intent(inout) :: w(ixi^s, nw)
1024 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1027 w(ixo^s,
e_)=w(ixo^s,
e_)-half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
1032 subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
1034 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1035 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
1036 double precision,
intent(out) :: v(ixi^s)
1038 v(ixo^s) = w(ixo^s,
mom(idim)) / w(ixo^s,
rho_)
1039 end subroutine hd_get_v_idim
1042 subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
1045 integer,
intent(in) :: ixi^
l, ixo^
l
1046 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:^nd)
1047 double precision,
intent(out) :: v(ixi^s,1:
ndir)
1052 v(ixo^s,idir) = w(ixo^s,
mom(idir)) / w(ixo^s,
rho_)
1055 end subroutine hd_get_v
1058 subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
1063 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1065 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
1066 double precision,
intent(inout) :: cmax(ixi^s)
1076 cmax(ixo^s)=dabs(w(ixo^s,
mom(idim)))+dsqrt(
hd_gamma*cmax(ixo^s)/w(ixo^s,
rho_))
1082 end subroutine hd_get_cmax
1085 subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
1087 integer,
intent(in) :: ixi^
l,ixo^
l
1088 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1090 double precision,
intent(inout) :: w(ixi^s,1:nw)
1091 double precision,
intent(out) :: tco_local, tmax_local
1093 double precision,
parameter :: trac_delta=0.25d0
1094 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
1095 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
1096 integer :: jxo^
l,hxo^
l
1097 integer :: jxp^
l,hxp^
l,ixp^
l
1098 logical :: lrlt(ixi^s)
1102 te(ixi^s)=w(ixi^s,
p_)/(te(ixi^s)*w(ixi^s,
rho_))
1105 tmax_local=maxval(te(ixo^s))
1112 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
1114 where(lts(ixo^s) > trac_delta)
1117 if(any(lrlt(ixo^s)))
then
1118 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
1129 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
1130 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
1131 w(ixo^s,
tcoff_)=te(ixo^s)*&
1132 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
1134 call mpistop(
"hd_trac_type not allowed for 1D simulation")
1137 end subroutine hd_get_tcutoff
1140 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
1145 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1147 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
1149 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
1150 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1152 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
1155 double precision :: wmean(ixi^s,
nw)
1156 double precision,
dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1164 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1165 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1166 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1167 umean(ixo^s)=(wlp(ixo^s,
mom(idim))*tmp1(ixo^s)+wrp(ixo^s,
mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1179 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1180 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1181 (wrp(ixo^s,
mom(idim))-wlp(ixo^s,
mom(idim)))**2
1183 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1184 if(
present(cmin))
then
1185 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1186 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1188 {
do ix^db=ixomin^db,ixomax^db\}
1189 cmin(ix^
d,1)=sign(one,cmin(ix^
d,1))*max(abs(cmin(ix^
d,1)),hspeed(ix^
d,1))
1190 cmax(ix^
d,1)=sign(one,cmax(ix^
d,1))*max(abs(cmax(ix^
d,1)),hspeed(ix^
d,1))
1194 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1198 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1199 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1210 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1211 tmp1(ixo^s)=wmean(ixo^s,
mom(idim))/wmean(ixo^s,
rho_)
1214 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1216 if(
present(cmin))
then
1217 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1218 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1219 if(h_correction)
then
1220 {
do ix^db=ixomin^db,ixomax^db\}
1221 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1222 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1226 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1230 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1243 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1244 if(
present(cmin))
then
1245 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))-csoundl(ixo^s)
1246 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1247 if(h_correction)
then
1248 {
do ix^db=ixomin^db,ixomax^db\}
1249 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1250 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1254 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1257 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1258 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1262 end subroutine hd_get_cbounds
1268 integer,
intent(in) :: ixi^
l, ixo^
l
1269 double precision,
intent(in) :: w(ixi^s,nw)
1270 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1271 double precision,
intent(out) :: csound2(ixi^s)
1284 integer,
intent(in) :: ixi^
l, ixo^
l
1285 double precision,
intent(in) :: w(ixi^s, 1:nw)
1286 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1287 double precision,
intent(out):: pth(ixi^s)
1291 pth(ixo^s) = (
hd_gamma - 1.0d0) * (w(ixo^s,
e_) - &
1302 {
do ix^db= ixo^lim^db\}
1308 {
do ix^db= ixo^lim^db\}
1310 write(*,*)
"Error: small value of gas pressure",pth(ix^
d),&
1311 " encountered when call hd_get_pthermal"
1313 write(*,*)
"Location: ", x(ix^
d,:)
1314 write(*,*)
"Cell number: ", ix^
d
1316 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^
d,iw)
1320 write(*,*)
"Saving status at the previous time step"
1333 integer,
intent(in) :: ixi^
l, ixo^
l
1334 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
1335 double precision,
intent(out):: csound(ixi^s)
1337 double precision :: wprim(ixi^s, nw)
1339 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1351 integer,
intent(in) :: ixi^
l, ixo^
l
1352 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
1353 double precision,
intent(out):: csound(ixi^s)
1356 double precision :: inv_rho
1357 double precision :: prad_tensor(ixi^s, 1:
ndim, 1:
ndim)
1358 double precision :: prad_max(ixi^s)
1362 {
do ix^db=ixomin^db,ixomax^db \}
1363 inv_rho=1.d0/w(ix^
d,
rho_)
1364 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1368 if(minval(csound(ixo^s))<smalldouble)
then
1369 print *,
'issue with squared speed and rad pressure'
1370 print *,minval(csound(ixo^s))
1371 print *,minval(prad_max(ixo^s))
1372 call mpistop(
"negative squared speed in get_csrad2 for dt")
1383 integer,
intent(in) :: ixi^
l, ixo^
l
1384 double precision,
intent(in) :: w(ixi^s, 1:nw)
1385 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1386 double precision,
intent(out):: prad(ixi^s, 1:
ndim, 1:
ndim)
1396 integer,
intent(in) :: ixi^
l, ixo^
l
1397 double precision,
intent(in) :: w(ixi^s, 1:nw)
1398 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1399 double precision,
intent(out):: pth_plus_prad(ixi^s)
1401 double precision :: wprim(ixi^s, 1:nw)
1402 double precision :: prad_tensor(ixi^s, 1:
ndim, 1:
ndim)
1403 double precision :: prad_max(ixi^s)
1406 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1409 {
do ix^
d = ixomin^
d,ixomax^
d\}
1410 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1412 pth_plus_prad(ixo^s) = wprim(ixo^s,
p_) + prad_max(ixo^s)
1420 integer,
intent(in) :: ixi^
l, ixo^
l
1421 double precision,
intent(in) :: w(ixi^s, 1:nw)
1422 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1423 double precision,
intent(out):: trad(ixi^s)
1432 integer,
intent(in) :: ixi^
l, ixo^
l
1433 double precision,
intent(in) :: w(ixi^s, 1:nw)
1434 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1435 double precision,
intent(out):: res(ixi^s)
1437 double precision :: r(ixi^s)
1441 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1445 subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1447 integer,
intent(in) :: ixi^
l, ixo^
l
1448 double precision,
intent(in) :: w(ixi^s, 1:nw)
1449 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1450 double precision,
intent(out):: res(ixi^s)
1452 double precision :: r(ixi^s)
1455 res(ixo^s) = (
hd_gamma - 1.0d0) * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1456 end subroutine hd_get_temperature_from_eint
1461 integer,
intent(in) :: ixi^
l, ixo^
l
1462 double precision,
intent(in) :: w(ixi^s, 1:nw)
1463 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1464 double precision,
intent(out):: res(ixi^s)
1466 double precision :: r(ixi^s)
1469 res(ixo^s)=w(ixo^s,
p_)/(r(ixo^s)*w(ixo^s,
rho_))
1474 subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1478 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1480 double precision,
intent(in) :: wc(ixi^s, 1:nw)
1482 double precision,
intent(in) :: w(ixi^s, 1:nw)
1483 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1484 double precision,
intent(out) :: f(ixi^s, nwflux)
1486 double precision :: pth(ixi^s)
1490 {
do ix^db=ixomin^db,ixomax^db\}
1500 {
do ix^db=ixomin^db,ixomax^db\}
1503 ^
c&f(ix^d,
m^
c_)=w(ix^d,
mom(idim))*wc(ix^d,
m^
c_)\
1504 f(ix^d,
mom(idim))=f(ix^d,
mom(idim))+pth(ix^d)
1509 {
do ix^db=ixomin^db,ixomax^db\}
1511 f(ix^d,
r_e)=w(ix^d,
mom(idim))*wc(ix^d,
r_e)
1520 f(ixo^s,
fip_) = w(ixo^s,
mom(idim)) * wc(ixo^s,
fip_)
1525 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1528 end subroutine hd_get_flux
1535 subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1541 integer,
intent(in) :: ixi^
l, ixo^
l
1542 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:
ndim)
1543 double precision,
intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
1544 double precision :: pth(ixi^s),
source(ixi^s), minrho
1545 integer :: iw,idir, h1x^
l{^nooned, h2x^
l}
1546 integer :: mr_,mphi_
1547 integer :: irho, ifluid, n_fluids
1548 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1570 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1574 do ifluid = 0, n_fluids-1
1576 if (ifluid == 0)
then
1600 where (wct(ixo^s, irho) > minrho)
1601 source(ixo^s) =
source(ixo^s) + wct(ixo^s,mphi_)*wprim(ixo^s,mphi_)
1602 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt*
source(ixo^s)/x(ixo^s,
r_)
1605 where (wct(ixo^s, irho) > minrho)
1606 source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s, mr_)
1607 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1611 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1616 call mpistop(
"Dust geom source terms not implemented yet with spherical geometries")
1620 h1x^
l=ixo^
l-
kr(1,^
d); {^nooned h2x^
l=ixo^
l-
kr(2,^
d);}
1622 pth(ixo^s)=wprim(ixo^s,
p_)
1631 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1632 *(
block%surfaceC(ixo^s, 1) -
block%surfaceC(h1x^s, 1)) &
1633 /
block%dvolume(ixo^s)
1637 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
1641 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1642 * (
block%surfaceC(ixo^s, 2) -
block%surfaceC(h2x^s, 2)) &
1643 /
block%dvolume(ixo^s)
1645 source(ixo^s) =
source(ixo^s) + (wprim(ixo^s,
mom(3))**2 * wprim(ixo^s,
rho_)) / tan(x(ixo^s, 2))
1647 source(ixo^s) =
source(ixo^s) - (wprim(ixo^s,
mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)
1648 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1652 source(ixo^s) = -(wprim(ixo^s,
mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)&
1653 - (wprim(ixo^s,
mom(2)) * wprim(ixo^s,
mom(3))) * wprim(ixo^s,
rho_) / tan(x(ixo^s, 2))
1654 w(ixo^s,
mom(3)) = w(ixo^s,
mom(3)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1661 call mpistop(
"Rotating frame not implemented yet with dust")
1667 end subroutine hd_add_source_geom
1670 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1679 integer,
intent(in) :: ixi^
l, ixo^
l
1680 double precision,
intent(in) :: qdt, dtfactor
1681 double precision,
intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:
ndim)
1682 double precision,
intent(inout) :: w(ixi^s, 1:nw)
1683 logical,
intent(in) :: qsourcesplit
1684 logical,
intent(inout) :: active
1686 double precision :: gravity_field(ixi^s, 1:
ndim)
1687 integer :: idust, idim
1695 qsourcesplit,active,
rc_fl)
1714 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s,
dust_rho(idust))
1726 call hd_add_radiation_source(qdt,ixi^
l,ixo^
l,wct,wctprim,w,x,qsourcesplit,active)
1730 if(.not.qsourcesplit)
then
1732 call hd_update_temperature(ixi^
l,ixo^
l,wct,w,x)
1736 end subroutine hd_add_source
1738 subroutine hd_add_radiation_source(qdt,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1744 integer,
intent(in) :: ixi^
l, ixo^
l
1745 double precision,
intent(in) :: qdt, x(ixi^s,1:
ndim)
1746 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
1747 double precision,
intent(inout) :: w(ixi^s,1:nw)
1748 logical,
intent(in) :: qsourcesplit
1749 logical,
intent(inout) :: active
1755 end subroutine hd_add_radiation_source
1757 subroutine hd_get_dt(wprim, ixI^L, ixO^L, dtnew, dx^D, x)
1765 integer,
intent(in) :: ixi^
l, ixo^
l
1766 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
1767 double precision,
intent(in) :: wprim(ixi^s, 1:nw)
1768 double precision,
intent(inout) :: dtnew
1792 end subroutine hd_get_dt
1796 integer,
intent(in) :: ixi^
l, ixo^
l
1797 double precision,
intent(in) :: w(ixi^s, nw)
1798 double precision :: ke(ixo^s)
1799 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1801 if (
present(inv_rho))
then
1802 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * inv_rho
1804 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) / w(ixo^s,
rho_)
1808 function hd_inv_rho(w, ixI^L, ixO^L)
result(inv_rho)
1810 integer,
intent(in) :: ixi^
l, ixo^
l
1811 double precision,
intent(in) :: w(ixi^s, nw)
1812 double precision :: inv_rho(ixo^s)
1815 inv_rho = 1.0d0 / w(ixo^s,
rho_)
1816 end function hd_inv_rho
1818 subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1825 logical,
intent(in) :: primitive
1826 integer,
intent(in) :: ixi^
l,ixo^
l
1827 double precision,
intent(inout) :: w(ixi^s,1:nw)
1828 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1829 character(len=*),
intent(in) :: subname
1832 logical :: flag(ixi^s,1:nw)
1842 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(idir)) = 0.0d0
1864 where(flag(ixo^s,
e_))
1889 -0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_))
1892 +0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_)
1908 if(.not.primitive)
then
1916 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir))/w(ixo^s,
rho_)
1923 if (
hd_fip)
call hd_bound_fip(primitive, ixi^
l, ixo^
l, w)
1924 end subroutine hd_handle_small_values
1926 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1929 integer,
intent(in) :: ixi^
l, ixo^
l
1930 double precision,
intent(in) :: w(ixi^s,1:nw)
1931 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1932 double precision,
intent(out):: rfactor(ixi^s)
1934 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1938 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
1940 end subroutine rfactor_from_temperature_ionization
1942 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1944 integer,
intent(in) :: ixi^
l, ixo^
l
1945 double precision,
intent(in) :: w(ixi^s,1:nw)
1946 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1947 double precision,
intent(out):: rfactor(ixi^s)
1951 end subroutine rfactor_from_constant_ionization
1953 subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1957 integer,
intent(in) :: ixi^
l, ixo^
l
1958 double precision,
intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:
ndim)
1959 double precision,
intent(inout) :: w(ixi^s,1:nw)
1961 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1970 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_fip
Whether FIP passive scalar is enabled.
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 fip_
Index of the FIP passive scalar rho*fip in conserved form, fip in primitive form.
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.