24 logical,
public,
protected ::
hd_dust = .false.
48 integer,
public,
protected ::
rho_
51 integer,
allocatable,
public,
protected ::
mom(:)
54 integer,
public,
protected :: ^
c&m^C_
57 integer,
allocatable,
public,
protected ::
tracer(:)
60 integer,
public,
protected ::
e_
63 integer,
public,
protected ::
p_
66 integer,
public,
protected ::
te_
72 double precision,
public ::
hd_gamma = 5.d0/3.0d0
75 double precision :: gamma_1, inv_gamma_1
81 double precision,
protected :: small_e
84 logical,
public,
protected ::
hd_trac = .false.
91 double precision,
public,
protected ::
h_ion_fr=1d0
94 double precision,
public,
protected ::
he_ion_fr=1d0
101 double precision,
public,
protected ::
rr=1d0
123 subroutine hd_read_params(files)
125 character(len=*),
intent(in) :: files(:)
134 do n = 1,
size(files)
135 open(
unitpar, file=trim(files(n)), status=
"old")
136 read(
unitpar, hd_list,
end=111)
140 end subroutine hd_read_params
143 subroutine hd_write_info(fh)
145 integer,
intent(in) :: fh
146 integer,
parameter :: n_par = 1
147 double precision :: values(n_par)
148 character(len=name_len) :: names(n_par)
149 integer,
dimension(MPI_STATUS_SIZE) :: st
152 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
156 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
157 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
158 end subroutine hd_write_info
183 phys_internal_e = .false.
192 if(
mype==0)
write(*,*)
'WARNING: set hd_trac_type=1'
197 if(
mype==0)
write(*,*)
'WARNING: set hd_trac=F when ndim>=2'
205 if(
mype==0)
write(*,*)
'WARNING: set hd_thermal_conduction=F when hd_energy=F'
209 if(
mype==0)
write(*,*)
'WARNING: set hd_radiative_cooling=F when hd_energy=F'
215 if(
mype==0)
write(*,*)
'WARNING: set hd_partial_ionization=F when eq_state_units=F'
220 allocate(start_indices(number_species),stop_indices(number_species))
228 mom(:) = var_set_momentum(
ndir)
233 e_ = var_set_energy()
240 phys_get_dt => hd_get_dt
241 phys_get_cmax => hd_get_cmax
242 phys_get_a2max => hd_get_a2max
243 phys_get_tcutoff => hd_get_tcutoff
244 phys_get_cbounds => hd_get_cbounds
245 phys_get_flux => hd_get_flux
246 phys_add_source_geom => hd_add_source_geom
247 phys_add_source => hd_add_source
253 phys_get_v => hd_get_v
254 phys_get_rho => hd_get_rho
255 phys_write_info => hd_write_info
256 phys_handle_small_values => hd_handle_small_values
259 call hd_physical_units()
269 tracer(itr) = var_set_fluxvar(
"trc",
"trp", itr, need_bc=.false.)
276 stop_indices(1)=nwflux
280 te_ = var_set_auxvar(
'Te',
'Te')
294 hd_get_rfactor=>rfactor_from_temperature_ionization
295 phys_update_temperature => hd_update_temperature
299 hd_get_rfactor=>rfactor_from_constant_ionization
305 call mpistop(
"thermal conduction needs hd_energy=T")
315 tc_fl%get_temperature_from_conserved => hd_get_temperature_from_etot
316 tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
317 tc_fl%get_rho => hd_get_rho
324 call mpistop(
"radiative cooling needs hd_energy=T")
328 rc_fl%get_rho => hd_get_rho
330 rc_fl%get_var_Rfactor => hd_get_rfactor
337 te_fl_hd%get_var_Rfactor => hd_get_rfactor
339 phys_te_images => hd_te_images
359 if (.not.
allocated(flux_type))
then
360 allocate(flux_type(
ndir, nw))
361 flux_type = flux_default
362 else if (any(shape(flux_type) /= [
ndir, nw]))
then
363 call mpistop(
"phys_check error: flux_type has wrong shape")
367 allocate(iw_vector(nvector))
368 iw_vector(1) =
mom(1) - 1
375 subroutine hd_te_images
379 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
381 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
383 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
385 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
388 call mpistop(
"Error in synthesize emission: Unknown convert_type")
390 end subroutine hd_te_images
395 subroutine hd_sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
399 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
400 double precision,
intent(in) :: x(ixi^s,1:
ndim)
401 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
402 double precision,
intent(in) :: my_dt
403 logical,
intent(in) :: fix_conserve_at_step
405 end subroutine hd_sts_set_source_tc_hd
407 function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x)
result(dtnew)
413 integer,
intent(in) :: ixi^
l, ixo^
l
414 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
415 double precision,
intent(in) :: w(ixi^s,1:nw)
416 double precision :: dtnew
419 end function hd_get_tc_dt_hd
421 subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
426 integer,
intent(in) :: ixi^
l,ixo^
l
427 double precision,
intent(inout) :: w(ixi^s,1:nw)
428 double precision,
intent(in) :: x(ixi^s,1:
ndim)
429 integer,
intent(in) :: step
432 logical :: flag(ixi^s,1:nw)
433 character(len=140) :: error_msg
436 where(w(ixo^s,
e_)<small_e) flag(ixo^s,
e_)=.true.
437 if(any(flag(ixo^s,
e_)))
then
440 where(flag(ixo^s,
e_)) w(ixo^s,
e_)=small_e
447 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,
rho_)
449 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
453 end subroutine hd_tc_handle_small_e
456 subroutine tc_params_read_hd(fl)
458 type(tc_fluid),
intent(inout) :: fl
460 logical :: tc_saturate=.false.
461 double precision :: tc_k_para=0d0
463 namelist /tc_list/ tc_saturate, tc_k_para
467 read(
unitpar, tc_list,
end=111)
470 fl%tc_saturate = tc_saturate
471 fl%tc_k_para = tc_k_para
473 end subroutine tc_params_read_hd
475 subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
477 integer,
intent(in) :: ixi^
l, ixo^
l
478 double precision,
intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:
ndim)
479 double precision,
intent(out) :: rho(ixi^s)
481 rho(ixo^s) = w(ixo^s,
rho_)
483 end subroutine hd_get_rho
487 subroutine rc_params_read(fl)
491 type(rc_fluid),
intent(inout) :: fl
494 integer :: ncool = 4000
497 character(len=std_len) :: coolcurve=
'JCcorona'
500 logical :: tfix=.false.
506 logical :: rc_split=.false.
509 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split
513 read(
unitpar, rc_list,
end=111)
518 fl%coolcurve=coolcurve
522 end subroutine rc_params_read
530 if (
hd_gamma <= 0.0d0)
call mpistop (
"Error: hd_gamma <= 0")
531 if (
hd_adiab < 0.0d0)
call mpistop (
"Error: hd_adiab < 0")
535 call mpistop (
"Error: hd_gamma <= 0 or hd_gamma == 1.0")
549 subroutine hd_physical_units
551 double precision :: mp,kb
552 double precision :: a,b
640 end subroutine hd_physical_units
647 logical,
intent(in) :: primitive
648 integer,
intent(in) :: ixi^
l, ixo^
l
649 double precision,
intent(in) :: w(ixi^s, nw)
650 logical,
intent(inout) :: flag(ixi^s,1:nw)
651 double precision :: tmp(ixi^s)
660 half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_))
675 integer,
intent(in) :: ixi^
l, ixo^
l
676 double precision,
intent(inout) :: w(ixi^s, nw)
677 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
681 {
do ix^db=ixomin^db,ixomax^db\}
684 w(ix^
d,
e_)=w(ix^
d,
e_)*inv_gamma_1+&
692 call dust_to_conserved(ixi^l, ixo^l, w, x)
701 integer,
intent(in) :: ixi^
l, ixo^
l
702 double precision,
intent(inout) :: w(ixi^s, nw)
703 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
705 double precision :: inv_rho
709 call hd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'hd_to_primitive')
712 {
do ix^db=ixomin^db,ixomax^db\}
713 inv_rho = 1.d0/w(ix^
d,
rho_)
726 call dust_to_primitive(ixi^l, ixo^l, w, x)
734 integer,
intent(in) :: ixi^
l, ixo^
l
735 double precision,
intent(inout) :: w(ixi^s, nw)
736 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
739 w(ixo^s,
e_)=w(ixo^s,
e_)+half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
746 integer,
intent(in) :: ixi^
l, ixo^
l
747 double precision,
intent(inout) :: w(ixi^s, nw)
748 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
751 w(ixo^s,
e_)=w(ixo^s,
e_)-half*(^
c&w(ixo^s,
m^
c_)**2+)/w(ixo^s,
rho_)
756 subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
758 integer,
intent(in) :: ixi^
l, ixo^
l, idim
759 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
760 double precision,
intent(out) :: v(ixi^s)
762 v(ixo^s) = w(ixo^s,
mom(idim)) / w(ixo^s,
rho_)
763 end subroutine hd_get_v_idim
766 subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
769 integer,
intent(in) :: ixi^
l, ixo^
l
770 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:^nd)
771 double precision,
intent(out) :: v(ixi^s,1:
ndir)
776 v(ixo^s,idir) = w(ixo^s,
mom(idir)) / w(ixo^s,
rho_)
779 end subroutine hd_get_v
782 subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
787 integer,
intent(in) :: ixi^
l, ixo^
l, idim
789 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s, 1:
ndim)
790 double precision,
intent(inout) :: cmax(ixi^s)
800 cmax(ixo^s)=dabs(w(ixo^s,
mom(idim)))+dsqrt(
hd_gamma*cmax(ixo^s)/w(ixo^s,
rho_))
806 end subroutine hd_get_cmax
808 subroutine hd_get_a2max(w,x,ixI^L,ixO^L,a2max)
811 integer,
intent(in) :: ixi^
l, ixo^
l
812 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
813 double precision,
intent(inout) :: a2max(
ndim)
814 double precision :: a2(ixi^s,
ndim,nw)
815 integer :: gxo^
l,hxo^
l,jxo^
l,kxo^
l,i,j
820 hxo^
l=ixo^
l-
kr(i,^
d);
821 gxo^
l=hxo^
l-
kr(i,^
d);
822 jxo^
l=ixo^
l+
kr(i,^
d);
823 kxo^
l=jxo^
l+
kr(i,^
d);
824 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
825 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
826 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/
dxlevel(i)**2
828 end subroutine hd_get_a2max
831 subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
833 integer,
intent(in) :: ixi^
l,ixo^
l
834 double precision,
intent(in) :: x(ixi^s,1:
ndim)
836 double precision,
intent(inout) :: w(ixi^s,1:nw)
837 double precision,
intent(out) :: tco_local, tmax_local
839 double precision,
parameter :: trac_delta=0.25d0
840 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
841 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
842 integer :: jxo^
l,hxo^
l
843 integer :: jxp^
l,hxp^
l,ixp^
l
844 logical :: lrlt(ixi^s)
847 call hd_get_rfactor(w,x,ixi^
l,ixi^
l,te)
848 te(ixi^s)=w(ixi^s,
p_)/(te(ixi^s)*w(ixi^s,
rho_))
851 tmax_local=maxval(te(ixo^s))
858 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
860 where(lts(ixo^s) > trac_delta)
863 if(any(lrlt(ixo^s)))
then
864 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
875 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
876 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
877 w(ixo^s,
tcoff_)=te(ixo^s)*&
878 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
880 call mpistop(
"hd_trac_type not allowed for 1D simulation")
883 end subroutine hd_get_tcutoff
886 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
891 integer,
intent(in) :: ixi^
l, ixo^
l, idim
893 double precision,
intent(in) :: wlc(ixi^s,
nw), wrc(ixi^s,
nw)
895 double precision,
intent(in) :: wlp(ixi^s,
nw), wrp(ixi^s,
nw)
896 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
898 double precision,
intent(inout),
optional :: cmin(ixi^s,1:
number_species)
901 double precision :: wmean(ixi^s,
nw)
902 double precision,
dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
910 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
911 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
912 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
913 umean(ixo^s)=(wlp(ixo^s,
mom(idim))*tmp1(ixo^s)+wrp(ixo^s,
mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
923 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
924 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
925 (wrp(ixo^s,
mom(idim))-wlp(ixo^s,
mom(idim)))**2
927 dmean(ixo^s)=dsqrt(dmean(ixo^s))
928 if(
present(cmin))
then
929 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
930 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
932 {
do ix^db=ixomin^db,ixomax^db\}
933 cmin(ix^
d,1)=sign(one,cmin(ix^
d,1))*max(abs(cmin(ix^
d,1)),hspeed(ix^
d,1))
934 cmax(ix^
d,1)=sign(one,cmax(ix^
d,1))*max(abs(cmax(ix^
d,1)),hspeed(ix^
d,1))
938 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
942 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
943 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
947 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
948 tmp1(ixo^s)=wmean(ixo^s,
mom(idim))/wmean(ixo^s,
rho_)
950 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
952 if(
present(cmin))
then
953 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
954 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
955 if(h_correction)
then
956 {
do ix^db=ixomin^db,ixomax^db\}
957 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
958 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
962 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
966 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
977 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
978 if(
present(cmin))
then
979 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))-csoundl(ixo^s)
980 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
981 if(h_correction)
then
982 {
do ix^db=ixomin^db,ixomax^db\}
983 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
984 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
988 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(idim)),wrp(ixo^s,
mom(idim)))+csoundl(ixo^s)
991 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
992 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
996 end subroutine hd_get_cbounds
1002 integer,
intent(in) :: ixi^
l, ixo^
l
1003 double precision,
intent(in) :: w(ixi^s,nw)
1004 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1005 double precision,
intent(out) :: csound2(ixi^s)
1018 integer,
intent(in) :: ixi^
l, ixo^
l
1019 double precision,
intent(in) :: w(ixi^s, 1:nw)
1020 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1021 double precision,
intent(out):: pth(ixi^s)
1025 pth(ixo^s) = (
hd_gamma - 1.0d0) * (w(ixo^s,
e_) - &
1036 {
do ix^db= ixo^lim^db\}
1042 {
do ix^db= ixo^lim^db\}
1044 write(*,*)
"Error: small value of gas pressure",pth(ix^
d),&
1045 " encountered when call hd_get_pthermal"
1047 write(*,*)
"Location: ", x(ix^
d,:)
1048 write(*,*)
"Cell number: ", ix^
d
1050 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^
d,iw)
1054 write(*,*)
"Saving status at the previous time step"
1063 subroutine hd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1065 integer,
intent(in) :: ixi^
l, ixo^
l
1066 double precision,
intent(in) :: w(ixi^s, 1:nw)
1067 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1068 double precision,
intent(out):: res(ixi^s)
1070 double precision :: r(ixi^s)
1072 call hd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1074 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1075 end subroutine hd_get_temperature_from_etot
1078 subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1080 integer,
intent(in) :: ixi^
l, ixo^
l
1081 double precision,
intent(in) :: w(ixi^s, 1:nw)
1082 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1083 double precision,
intent(out):: res(ixi^s)
1085 double precision :: r(ixi^s)
1087 call hd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1088 res(ixo^s) = (
hd_gamma - 1.0d0) * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1089 end subroutine hd_get_temperature_from_eint
1092 subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1096 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1098 double precision,
intent(in) :: wc(ixi^s, 1:nw)
1100 double precision,
intent(in) :: w(ixi^s, 1:nw)
1101 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1102 double precision,
intent(out) :: f(ixi^s, nwflux)
1104 double precision :: pth(ixi^s)
1108 {
do ix^db=ixomin^db,ixomax^db\}
1118 {
do ix^db=ixomin^db,ixomax^db\}
1121 ^
c&f(ix^d,
m^
c_)=w(ix^d,
mom(idim))*wc(ix^d,
m^
c_)\
1122 f(ix^d,
mom(idim))=f(ix^d,
mom(idim))+pth(ix^d)
1132 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1135 end subroutine hd_get_flux
1144 subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1150 integer,
intent(in) :: ixi^
l, ixo^
l
1151 double precision,
intent(in) :: qdt, dtfactor, x(ixi^s, 1:
ndim)
1152 double precision,
intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
1156 double precision :: pth(ixi^s),
source(ixi^s), minrho
1157 integer :: iw,idir, h1x^
l{^nooned, h2x^
l}
1158 integer :: mr_,mphi_
1159 integer :: irho, ifluid, n_fluids
1160 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1182 source(ixo^s) =
source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1186 do ifluid = 0, n_fluids-1
1188 if (ifluid == 0)
then
1212 where (wct(ixo^s, irho) > minrho)
1213 source(ixo^s) =
source(ixo^s) + wct(ixo^s,mphi_)*wprim(ixo^s,mphi_)
1214 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt*
source(ixo^s)/x(ixo^s,
r_)
1217 where (wct(ixo^s, irho) > minrho)
1218 source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s, mr_)
1219 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1223 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s,
r_)
1228 call mpistop(
"Dust geom source terms not implemented yet with spherical geometries")
1232 h1x^
l=ixo^
l-
kr(1,^
d); {^nooned h2x^
l=ixo^
l-
kr(2,^
d);}
1234 pth(ixo^s)=wprim(ixo^s,
p_)
1243 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1244 *(
block%surfaceC(ixo^s, 1) -
block%surfaceC(h1x^s, 1)) &
1245 /
block%dvolume(ixo^s)
1249 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt *
source(ixo^s) / x(ixo^s, 1)
1253 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1254 * (
block%surfaceC(ixo^s, 2) -
block%surfaceC(h2x^s, 2)) &
1255 /
block%dvolume(ixo^s)
1257 source(ixo^s) =
source(ixo^s) + (wprim(ixo^s,
mom(3))**2 * wprim(ixo^s,
rho_)) / tan(x(ixo^s, 2))
1259 source(ixo^s) =
source(ixo^s) - (wprim(ixo^s,
mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)
1260 w(ixo^s,
mom(2)) = w(ixo^s,
mom(2)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1264 source(ixo^s) = -(wprim(ixo^s,
mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s,
rho_)&
1265 - (wprim(ixo^s,
mom(2)) * wprim(ixo^s,
mom(3))) * wprim(ixo^s,
rho_) / tan(x(ixo^s, 2))
1266 w(ixo^s,
mom(3)) = w(ixo^s,
mom(3)) + qdt *
source(ixo^s) / x(ixo^s, 1)
1273 call mpistop(
"Rotating frame not implemented yet with dust")
1279 end subroutine hd_add_source_geom
1282 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1291 integer,
intent(in) :: ixi^
l, ixo^
l
1292 double precision,
intent(in) :: qdt, dtfactor
1293 double precision,
intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:
ndim)
1294 double precision,
intent(inout) :: w(ixi^s, 1:nw)
1295 logical,
intent(in) :: qsourcesplit
1296 logical,
intent(inout) :: active
1298 double precision :: gravity_field(ixi^s, 1:
ndim)
1299 integer :: idust, idim
1307 qsourcesplit,active,
rc_fl)
1326 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s,
dust_rho(idust))
1337 if(.not.qsourcesplit)
then
1339 call hd_update_temperature(ixi^
l,ixo^
l,wct,w,x)
1343 end subroutine hd_add_source
1345 subroutine hd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
1352 integer,
intent(in) :: ixi^
l, ixo^
l
1353 double precision,
intent(in) ::
dx^
d, x(ixi^s, 1:^nd)
1354 double precision,
intent(in) :: w(ixi^s, 1:nw)
1355 double precision,
intent(inout) :: dtnew
1375 end subroutine hd_get_dt
1379 integer,
intent(in) :: ixi^
l, ixo^
l
1380 double precision,
intent(in) :: w(ixi^s, nw)
1381 double precision :: ke(ixo^s)
1382 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1384 if (
present(inv_rho))
then
1385 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) * inv_rho
1387 ke = 0.5d0 * sum(w(ixo^s,
mom(:))**2, dim=
ndim+1) / w(ixo^s,
rho_)
1391 function hd_inv_rho(w, ixI^L, ixO^L)
result(inv_rho)
1393 integer,
intent(in) :: ixi^
l, ixo^
l
1394 double precision,
intent(in) :: w(ixi^s, nw)
1395 double precision :: inv_rho(ixo^s)
1398 inv_rho = 1.0d0 / w(ixo^s,
rho_)
1399 end function hd_inv_rho
1401 subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1408 logical,
intent(in) :: primitive
1409 integer,
intent(in) :: ixi^
l,ixo^
l
1410 double precision,
intent(inout) :: w(ixi^s,1:nw)
1411 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1412 character(len=*),
intent(in) :: subname
1415 logical :: flag(ixi^s,1:nw)
1425 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(idir)) = 0.0d0
1442 where(flag(ixo^s,
e_))
1467 -0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_))
1470 +0.5d0*sum(w(ixi^s,
mom(:))**2, dim=
ndim+1)/w(ixi^s,
rho_)
1482 if(.not.primitive)
then
1490 w(ixo^s,
mom(idir)) = w(ixo^s,
mom(idir))/w(ixo^s,
rho_)
1497 end subroutine hd_handle_small_values
1499 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1502 integer,
intent(in) :: ixi^
l, ixo^
l
1503 double precision,
intent(in) :: w(ixi^s,1:nw)
1504 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1505 double precision,
intent(out):: rfactor(ixi^s)
1507 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1511 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
1513 end subroutine rfactor_from_temperature_ionization
1515 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1517 integer,
intent(in) :: ixi^
l, ixo^
l
1518 double precision,
intent(in) :: w(ixi^s,1:nw)
1519 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1520 double precision,
intent(out):: rfactor(ixi^s)
1524 end subroutine rfactor_from_constant_ionization
1526 subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1530 integer,
intent(in) :: ixi^
l, ixo^
l
1531 double precision,
intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:
ndim)
1532 double precision,
intent(inout) :: w(ixi^s,1:nw)
1534 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1543 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_get_dt(w, ixil, ixol, dtnew, dxd, x)
Check time step for total radiation contribution.
subroutine cak_init(phys_gamma)
Initialize the module.
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)
subroutine, public dust_get_dt(w, ixil, ixol, dtnew, dxd, x)
Get dt related to dust and gas stopping time (Laibe 2011)
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_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 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.
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.
double precision unit_velocity
Physical scaling factor for velocity.
double precision unit_temperature
Physical scaling factor for temperature.
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
logical phys_trac
Use TRAC for MHD or 1D HD.
logical fix_small_values
fix small values with average or replace methods
logical crash
Save a snapshot before crash a run met unphysical values.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
double precision small_density
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module for including gravity in (magneto)hydrodynamics simulations.
logical grav_split
source split or not
subroutine gravity_get_dt(w, ixil, ixol, dtnew, dxd, x)
subroutine gravity_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_ei_to_e(ixil, ixol, w, x)
Transform internal energy to total energy.
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.
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 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.
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.
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_to_primitive(ixil, ixol, w, x)
Transform conservative variables into primitive ones.
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_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(w, ixil, ixol, dtnew, dxd, x)
procedure(sub_add_source), pointer, public viscosity_add_source