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
146 subroutine hd_read_params(files)
148 character(len=*),
intent(in) :: files(:)
158 do n = 1,
size(files)
159 open(
unitpar, file=trim(files(n)), status=
"old")
160 read(
unitpar, hd_list,
end=111)
164 end subroutine hd_read_params
167 subroutine hd_write_info(fh)
169 integer,
intent(in) :: fh
170 integer,
parameter :: n_par = 1
171 double precision :: values(n_par)
172 character(len=name_len) :: names(n_par)
173 integer,
dimension(MPI_STATUS_SIZE) :: st
176 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
180 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
181 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
182 end subroutine hd_write_info
207 phys_internal_e = .false.
216 if(
mype==0)
write(*,*)
'WARNING: set hd_trac_type=1'
221 if(
mype==0)
write(*,*)
'WARNING: set hd_trac=F when ndim>=2'
229 if(
mype==0)
write(*,*)
'WARNING: set hd_thermal_conduction=F when hd_energy=F'
233 if(
mype==0)
write(*,*)
'WARNING: set hd_radiative_cooling=F when hd_energy=F'
239 if(
mype==0)
write(*,*)
'WARNING: set hd_partial_ionization=F when eq_state_units=F'
244 allocate(start_indices(number_species),stop_indices(number_species))
252 mom(:) = var_set_momentum(
ndir)
257 e_ = var_set_energy()
267 write(*,*)
'Warning: CAK force addition together with FLD radiation'
272 write(*,*)
'Warning: Optically thin cooling together with FLD radiation'
276 call mpistop(
'implicit dust addition not compatible with FLD radiation')
279 call mpistop(
'using FLD implies the use of cgs units')
282 call mpistop(
'using FLD implies the use of an energy equation, set hd_energy=T')
285 r_e = var_set_radiation_energy()
295 phys_get_dt => hd_get_dt
296 phys_get_cmax => hd_get_cmax
297 phys_get_tcutoff => hd_get_tcutoff
298 phys_get_cbounds => hd_get_cbounds
299 phys_get_flux => hd_get_flux
300 phys_add_source_geom => hd_add_source_geom
301 phys_add_source => hd_add_source
307 phys_get_v => hd_get_v
308 phys_get_rho => hd_get_rho
309 phys_write_info => hd_write_info
310 phys_handle_small_values => hd_handle_small_values
313 call hd_physical_units()
323 tracer(itr) = var_set_fluxvar(
"trc",
"trp", itr, need_bc=.false.)
328 te_ = var_set_auxvar(
'Te',
'Te')
337 stop_indices(1)=nwflux
349 phys_update_temperature => hd_update_temperature
361 call mpistop(
"thermal conduction needs hd_energy=T")
371 tc_fl%get_temperature_from_conserved => hd_get_temperature_from_etot
372 tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
373 tc_fl%get_rho => hd_get_rho
380 call mpistop(
"radiative cooling needs hd_energy=T")
384 rc_fl%get_rho => hd_get_rho
395 phys_te_images => hd_te_images
411 if (.not.
allocated(flux_type))
then
412 allocate(flux_type(
ndir, nw))
413 flux_type = flux_default
414 else if (any(shape(flux_type) /= [
ndir, nw]))
then
415 call mpistop(
"phys_check error: flux_type has wrong shape")
419 allocate(iw_vector(nvector))
420 iw_vector(1) =
mom(1) - 1
427 subroutine hd_te_images
431 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
433 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
435 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
437 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
440 call mpistop(
"Error in synthesize emission: Unknown convert_type")
442 end subroutine hd_te_images
447 subroutine hd_sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
451 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
452 double precision,
intent(in) :: x(ixi^s,1:
ndim)
453 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
454 double precision,
intent(in) :: my_dt
455 logical,
intent(in) :: fix_conserve_at_step
457 end subroutine hd_sts_set_source_tc_hd
459 function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x)
result(dtnew)
465 integer,
intent(in) :: ixi^
l, ixo^
l
466 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
467 double precision,
intent(in) :: w(ixi^s,1:nw)
468 double precision :: dtnew
471 end function hd_get_tc_dt_hd
473 subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
478 integer,
intent(in) :: ixi^
l,ixo^
l
479 double precision,
intent(inout) :: w(ixi^s,1:nw)
480 double precision,
intent(in) :: x(ixi^s,1:
ndim)
481 integer,
intent(in) :: step
484 logical :: flag(ixi^s,1:nw)
485 character(len=140) :: error_msg
489 if(any(flag(ixo^s,
e_)))
then
499 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,
rho_)
501 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
505 end subroutine hd_tc_handle_small_e
508 subroutine tc_params_read_hd(fl)
510 type(tc_fluid),
intent(inout) :: fl
512 logical :: tc_saturate=.false.
513 double precision :: tc_k_para=0d0
515 namelist /tc_list/ tc_saturate, tc_k_para
519 read(
unitpar, tc_list,
end=111)
522 fl%tc_saturate = tc_saturate
523 fl%tc_k_para = tc_k_para
525 end subroutine tc_params_read_hd
527 subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
529 integer,
intent(in) :: ixi^
l, ixo^
l
530 double precision,
intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:
ndim)
531 double precision,
intent(out) :: rho(ixi^s)
533 rho(ixo^s) = w(ixo^s,
rho_)
535 end subroutine hd_get_rho
539 subroutine rc_params_read(fl)
543 type(rc_fluid),
intent(inout) :: fl
546 integer :: ncool = 4000
549 character(len=std_len) :: coolcurve=
'JCcorona'
552 logical :: tfix=.false.
558 logical :: rc_split=.false.
561 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split
565 read(
unitpar, rc_list,
end=111)
570 fl%coolcurve=coolcurve
574 end subroutine rc_params_read
582 use mod_particles,
only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
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(*,*)
'========================'
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
635 write(*,*)
'========================'
639 write(*,*)
'====HD run with settings===================='
640 write(*,*)
'Using mod_hd_phys with settings:'
642 write(*,*)
'Dimensionality :',
ndim
643 write(*,*)
'vector components:',
ndir
645 write(*,*)
'number of variables nw=',nw
646 write(*,*)
' start index iwstart=',iwstart
647 write(*,*)
'number of vector variables=',nvector
648 write(*,*)
'number of stagger variables nws=',nws
649 write(*,*)
'number of variables with BCs=',nwgc
650 write(*,*)
'number of vars with fluxes=',nwflux
651 write(*,*)
'number of vars with flux + BC=',nwfluxbc
652 write(*,*)
'number of auxiliary variables=',nwaux
653 write(*,*)
'number of extra vars without flux=',nwextra
654 write(*,*)
'number of extra vars for wextra=',nw_extra
655 write(*,*)
'number of auxiliary I/O variables=',
nwauxio
669 write(*,*)
'*****Using particles: npayload,ngridvars :', npayload,ngridvars
670 write(*,*)
'*****Using particles: nusrpayload :', nusrpayload
671 write(*,*)
'*****Using particles: num_particles :', num_particles
672 write(*,*)
'*****Using particles: physics_type_particles=',physics_type_particles
675 write(*,*)
'number due to phys_wider_stencil=',phys_wider_stencil
676 write(*,*)
'==========================================='
681 subroutine hd_physical_units
683 double precision :: mp,kb
684 double precision :: a,b
776 end subroutine hd_physical_units
783 logical,
intent(in) :: primitive
784 integer,
intent(in) :: ixi^
l, ixo^
l
785 double precision,
intent(in) :: w(ixi^s, nw)
786 logical,
intent(inout) :: flag(ixi^s,1:nw)
787 double precision :: tmp(ixi^s)
796 half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_))
814 integer,
intent(in) :: ixi^
l, ixo^
l
815 double precision,
intent(inout) :: w(ixi^s, nw)
816 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
820 {
do ix^db=ixomin^db,ixomax^db\}
823 w(ix^
d,
e_)=w(ix^
d,
e_)*inv_gamma_1+&
831 call dust_to_conserved(ixi^l, ixo^l, w, x)
840 integer,
intent(in) :: ixi^
l, ixo^
l
841 double precision,
intent(inout) :: w(ixi^s, nw)
842 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
844 double precision :: inv_rho
848 call hd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'hd_to_primitive')
851 {
do ix^db=ixomin^db,ixomax^db\}
852 inv_rho = 1.d0/w(ix^
d,
rho_)
865 call dust_to_primitive(ixi^l, ixo^l, w, x)
873 integer,
intent(in) :: ixi^
l, ixo^
l
874 double precision,
intent(inout) :: w(ixi^s, nw)
875 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
878 w(ixo^s,
e_)=w(ixo^s,
e_)+half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
885 integer,
intent(in) :: ixi^
l, ixo^
l
886 double precision,
intent(inout) :: w(ixi^s, nw)
887 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
890 w(ixo^s,
e_)=w(ixo^s,
e_)-half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
895 subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
897 integer,
intent(in) :: ixi^
l, ixo^
l, idim
898 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
899 double precision,
intent(out) :: v(ixi^s)
901 v(ixo^s) = w(ixo^s,
mom(idim)) / w(ixo^s,
rho_)
902 end subroutine hd_get_v_idim
905 subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
908 integer,
intent(in) :: ixi^
l, ixo^
l
909 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:^nd)
910 double precision,
intent(out) :: v(ixi^s,1:
ndir)
915 v(ixo^s,idir) = w(ixo^s,
mom(idir)) / w(ixo^s,
rho_)
918 end subroutine hd_get_v
921 subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
926 integer,
intent(in) :: ixi^
l, ixo^
l, idim
928 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
929 double precision,
intent(inout) :: cmax(ixi^s)
939 cmax(ixo^s)=dabs(w(ixo^s,
mom(idim)))+dsqrt(
hd_gamma*cmax(ixo^s)/w(ixo^s,
rho_))
945 end subroutine hd_get_cmax
948 subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
950 integer,
intent(in) :: ixi^
l,ixo^
l
951 double precision,
intent(in) :: x(ixi^s,1:
ndim)
953 double precision,
intent(inout) :: w(ixi^s,1:nw)
954 double precision,
intent(out) :: tco_local, tmax_local
956 double precision,
parameter :: trac_delta=0.25d0
957 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
958 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
959 integer :: jxo^
l,hxo^
l
960 integer :: jxp^
l,hxp^
l,ixp^
l
961 logical :: lrlt(ixi^s)
965 te(ixi^s)=w(ixi^s,
p_)/(te(ixi^s)*w(ixi^s,
rho_))
968 tmax_local=maxval(te(ixo^s))
975 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
977 where(lts(ixo^s) > trac_delta)
980 if(any(lrlt(ixo^s)))
then
981 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
992 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
993 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
994 w(ixo^s,
tcoff_)=te(ixo^s)*&
995 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
997 call mpistop(
"hd_trac_type not allowed for 1D simulation")
1000 end subroutine hd_get_tcutoff
1003 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
1008 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1010 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
1012 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
1013 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1015 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
1018 double precision :: wmean(ixi^s,
nw)
1019 double precision,
dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1027 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1028 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1029 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1030 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)
1073 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1074 tmp1(ixo^s)=wmean(ixo^s,
mom(idim))/wmean(ixo^s,
rho_)
1077 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1079 if(
present(cmin))
then
1080 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1081 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1082 if(h_correction)
then
1083 {
do ix^db=ixomin^db,ixomax^db\}
1084 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1085 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1089 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1093 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1106 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1107 if(
present(cmin))
then
1108 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))-csoundl(ixo^s)
1109 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1110 if(h_correction)
then
1111 {
do ix^db=ixomin^db,ixomax^db\}
1112 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1113 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1117 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
1120 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1121 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1125 end subroutine hd_get_cbounds
1131 integer,
intent(in) :: ixi^
l, ixo^
l
1132 double precision,
intent(in) :: w(ixi^s,nw)
1133 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1134 double precision,
intent(out) :: csound2(ixi^s)
1147 integer,
intent(in) :: ixi^
l, ixo^
l
1148 double precision,
intent(in) :: w(ixi^s, 1:nw)
1149 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1150 double precision,
intent(out):: pth(ixi^s)
1154 pth(ixo^s) = (
hd_gamma - 1.0d0) * (w(ixo^s,
e_) - &
1165 {
do ix^db= ixo^lim^db\}
1171 {
do ix^db= ixo^lim^db\}
1173 write(*,*)
"Error: small value of gas pressure",pth(ix^
d),&
1174 " encountered when call hd_get_pthermal"
1176 write(*,*)
"Location: ", x(ix^
d,:)
1177 write(*,*)
"Cell number: ", ix^
d
1179 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^
d,iw)
1183 write(*,*)
"Saving status at the previous time step"
1196 integer,
intent(in) :: ixi^
l, ixo^
l
1197 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
1198 double precision,
intent(out):: csound(ixi^s)
1200 double precision :: wprim(ixi^s, nw)
1202 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1214 integer,
intent(in) :: ixi^
l, ixo^
l
1215 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
1216 double precision,
intent(out):: csound(ixi^s)
1219 double precision :: inv_rho
1220 double precision :: prad_tensor(ixi^s, 1:
ndim, 1:
ndim)
1221 double precision :: prad_max(ixi^s)
1225 {
do ix^db=ixomin^db,ixomax^db \}
1226 inv_rho=1.d0/w(ix^
d,
rho_)
1227 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1231 if(minval(csound(ixo^s))<smalldouble)
then
1232 print *,
'issue with squared speed and rad pressure'
1233 print *,minval(csound(ixo^s))
1234 print *,minval(prad_max(ixo^s))
1235 call mpistop(
"negative squared speed in get_csrad2 for dt")
1246 integer,
intent(in) :: ixi^
l, ixo^
l
1247 double precision,
intent(in) :: w(ixi^s, 1:nw)
1248 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1249 double precision,
intent(out):: prad(ixi^s, 1:
ndim, 1:
ndim)
1259 integer,
intent(in) :: ixi^
l, ixo^
l
1260 double precision,
intent(in) :: w(ixi^s, 1:nw)
1261 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1262 double precision,
intent(out):: pth_plus_prad(ixi^s)
1264 double precision :: wprim(ixi^s, 1:nw)
1265 double precision :: prad_tensor(ixi^s, 1:
ndim, 1:
ndim)
1266 double precision :: prad_max(ixi^s)
1269 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1272 {
do ix^
d = ixomin^
d,ixomax^
d\}
1273 prad_max(ix^
d) = maxval(prad_tensor(ix^
d,:,:))
1275 pth_plus_prad(ixo^s) = wprim(ixo^s,
p_) + prad_max(ixo^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):: trad(ixi^s)
1295 subroutine hd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1297 integer,
intent(in) :: ixi^
l, ixo^
l
1298 double precision,
intent(in) :: w(ixi^s, 1:nw)
1299 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1300 double precision,
intent(out):: res(ixi^s)
1302 double precision :: r(ixi^s)
1306 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1307 end subroutine hd_get_temperature_from_etot
1310 subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1312 integer,
intent(in) :: ixi^
l, ixo^
l
1313 double precision,
intent(in) :: w(ixi^s, 1:nw)
1314 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1315 double precision,
intent(out):: res(ixi^s)
1317 double precision :: r(ixi^s)
1320 res(ixo^s) = (
hd_gamma - 1.0d0) * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1321 end subroutine hd_get_temperature_from_eint
1326 integer,
intent(in) :: ixi^
l, ixo^
l
1327 double precision,
intent(in) :: w(ixi^s, 1:nw)
1328 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1329 double precision,
intent(out):: res(ixi^s)
1331 double precision :: r(ixi^s)
1334 res(ixo^s)=w(ixo^s,
p_)/(r(ixo^s)*w(ixo^s,
rho_))
1339 subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1343 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1345 double precision,
intent(in) :: wc(ixi^s, 1:nw)
1347 double precision,
intent(in) :: w(ixi^s, 1:nw)
1348 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1349 double precision,
intent(out) :: f(ixi^s, nwflux)
1351 double precision :: pth(ixi^s)
1355 {
do ix^db=ixomin^db,ixomax^db\}
1365 {
do ix^db=ixomin^db,ixomax^db\}
1368 ^
c&f(ix^d,
m^
c_)=w(ix^d,
mom(idim))*wc(ix^d,
m^
c_)\
1369 f(ix^d,
mom(idim))=f(ix^d,
mom(idim))+pth(ix^d)
1374 {
do ix^db=ixomin^db,ixomax^db\}
1376 f(ix^d,
r_e)=w(ix^d,
mom(idim))*wc(ix^d,
r_e)
1386 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1389 end subroutine hd_get_flux
1396 subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1402 integer,
intent(in) :: ixi^
l, ixo^
l
1403 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:
ndim)
1404 double precision,
intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
1405 double precision :: pth(ixi^s),
source(ixi^s), minrho
1406 integer :: iw,idir, h1x^
l{^nooned, h2x^
l}
1407 integer :: mr_,mphi_
1408 integer :: irho, ifluid, n_fluids
1409 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1431 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1435 do ifluid = 0, n_fluids-1
1437 if (ifluid == 0)
then
1461 where (wct(ixo^s, irho) > minrho)
1462 source(ixo^s) =
source(ixo^s) + wct(ixo^s,mphi_)*wprim(ixo^s,mphi_)
1463 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt*
source(ixo^s)/x(ixo^s,
r_)
1466 where (wct(ixo^s, irho) > minrho)
1467 source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s, mr_)
1468 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1472 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1477 call mpistop(
"Dust geom source terms not implemented yet with spherical geometries")
1481 h1x^
l=ixo^
l-
kr(1,^
d); {^nooned h2x^
l=ixo^
l-
kr(2,^
d);}
1483 pth(ixo^s)=wprim(ixo^s,
p_)
1492 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1493 *(
block%surfaceC(ixo^s, 1) -
block%surfaceC(h1x^s, 1)) &
1494 /
block%dvolume(ixo^s)
1498 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
1502 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1503 * (
block%surfaceC(ixo^s, 2) -
block%surfaceC(h2x^s, 2)) &
1504 /
block%dvolume(ixo^s)
1506 source(ixo^s) =
source(ixo^s) + (wprim(ixo^s,
mom(3))**2 * wprim(ixo^s,
rho_)) / tan(x(ixo^s, 2))
1508 source(ixo^s) =
source(ixo^s) - (wprim(ixo^s,
mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)
1509 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1513 source(ixo^s) = -(wprim(ixo^s,
mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)&
1514 - (wprim(ixo^s,
mom(2)) * wprim(ixo^s,
mom(3))) * wprim(ixo^s,
rho_) / tan(x(ixo^s, 2))
1515 w(ixo^s,
mom(3)) = w(ixo^s,
mom(3)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1522 call mpistop(
"Rotating frame not implemented yet with dust")
1528 end subroutine hd_add_source_geom
1531 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1540 integer,
intent(in) :: ixi^
l, ixo^
l
1541 double precision,
intent(in) :: qdt, dtfactor
1542 double precision,
intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:
ndim)
1543 double precision,
intent(inout) :: w(ixi^s, 1:nw)
1544 logical,
intent(in) :: qsourcesplit
1545 logical,
intent(inout) :: active
1547 double precision :: gravity_field(ixi^s, 1:
ndim)
1548 integer :: idust, idim
1556 qsourcesplit,active,
rc_fl)
1575 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s,
dust_rho(idust))
1587 call hd_add_radiation_source(qdt,ixi^
l,ixo^
l,wct,wctprim,w,x,qsourcesplit,active)
1591 if(.not.qsourcesplit)
then
1593 call hd_update_temperature(ixi^
l,ixo^
l,wct,w,x)
1597 end subroutine hd_add_source
1599 subroutine hd_add_radiation_source(qdt,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1605 integer,
intent(in) :: ixi^
l, ixo^
l
1606 double precision,
intent(in) :: qdt, x(ixi^s,1:
ndim)
1607 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
1608 double precision,
intent(inout) :: w(ixi^s,1:nw)
1609 logical,
intent(in) :: qsourcesplit
1610 logical,
intent(inout) :: active
1615 if(
fix_small_values)
call hd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'fld_add_radiation')
1617 end subroutine hd_add_radiation_source
1619 subroutine hd_get_dt(wprim, ixI^L, ixO^L, dtnew, dx^D, x)
1627 integer,
intent(in) :: ixi^
l, ixo^
l
1628 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
1629 double precision,
intent(in) :: wprim(ixi^s, 1:nw)
1630 double precision,
intent(inout) :: dtnew
1654 end subroutine hd_get_dt
1658 integer,
intent(in) :: ixi^
l, ixo^
l
1659 double precision,
intent(in) :: w(ixi^s, nw)
1660 double precision :: ke(ixo^s)
1661 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1663 if (
present(inv_rho))
then
1664 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * inv_rho
1666 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) / w(ixo^s,
rho_)
1670 function hd_inv_rho(w, ixI^L, ixO^L)
result(inv_rho)
1672 integer,
intent(in) :: ixi^
l, ixo^
l
1673 double precision,
intent(in) :: w(ixi^s, nw)
1674 double precision :: inv_rho(ixo^s)
1677 inv_rho = 1.0d0 / w(ixo^s,
rho_)
1678 end function hd_inv_rho
1680 subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1687 logical,
intent(in) :: primitive
1688 integer,
intent(in) :: ixi^
l,ixo^
l
1689 double precision,
intent(inout) :: w(ixi^s,1:nw)
1690 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1691 character(len=*),
intent(in) :: subname
1694 logical :: flag(ixi^s,1:nw)
1704 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(idir)) = 0.0d0
1726 where(flag(ixo^s,
e_))
1751 -0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_))
1754 +0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_)
1770 if(.not.primitive)
then
1778 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir))/w(ixo^s,
rho_)
1785 end subroutine hd_handle_small_values
1787 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1790 integer,
intent(in) :: ixi^
l, ixo^
l
1791 double precision,
intent(in) :: w(ixi^s,1:nw)
1792 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1793 double precision,
intent(out):: rfactor(ixi^s)
1795 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1799 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
1801 end subroutine rfactor_from_temperature_ionization
1803 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1805 integer,
intent(in) :: ixi^
l, ixo^
l
1806 double precision,
intent(in) :: w(ixi^s,1:nw)
1807 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1808 double precision,
intent(out):: rfactor(ixi^s)
1812 end subroutine rfactor_from_constant_ionization
1814 subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1818 integer,
intent(in) :: ixi^
l, ixo^
l
1819 double precision,
intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:
ndim)
1820 double precision,
intent(inout) :: w(ixi^s,1:nw)
1822 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1831 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.
double precision, parameter const_rad_a
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.
character(len=8) fld_opacity_law
switches for opacity
character(len=16) fld_fluxlimiter
flux limiter choice
character(len=50) fld_opal_table
double precision, public fld_kappa0
Opacity per unit of unit_density.
subroutine, public fld_init(he_abundance, r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
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...
character(len=8) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
logical fld_radforce_split
source split for energy interact and radforce:
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 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.
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.
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 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
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.
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 Note use of cgs units here in factor const_rad_a.
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 viscosity_init(phys_wider_stencil)
Initialize the module.
subroutine viscosity_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
procedure(sub_add_source), pointer, public viscosity_add_source