54 integer,
public,
protected ::
rho_
57 integer,
allocatable,
public,
protected ::
mom(:)
60 integer,
public,
protected ::
e_
63 integer,
public,
protected ::
p_
66 integer,
public,
protected ::
te_
71 integer,
public,
protected ::
q_
80 double precision,
protected :: small_e
89 double precision,
public,
protected ::
h_ion_fr=1d0
92 double precision,
public,
protected ::
he_ion_fr=1d0
99 double precision,
public,
protected ::
rr=1d0
106 double precision :: gamma_1, inv_gamma_1
111 function fun_kin_en(w, ixI^L, ixO^L, inv_rho)
result(ke)
113 integer,
intent(in) :: ixi^
l, ixo^
l
114 double precision,
intent(in) :: w(ixi^s, nw)
115 double precision :: ke(ixo^s)
116 double precision,
intent(in),
optional :: inv_rho(ixo^s)
117 end function fun_kin_en
145 subroutine ffhd_read_params(files)
148 character(len=*),
intent(in) :: files(:)
156 do n = 1,
size(files)
157 open(
unitpar, file=trim(files(n)), status=
"old")
158 read(
unitpar, ffhd_list,
end=111)
161 end subroutine ffhd_read_params
164 subroutine ffhd_write_info(fh)
166 integer,
intent(in) :: fh
167 integer,
parameter :: n_par = 1
168 double precision :: values(n_par)
169 character(len=name_len) :: names(n_par)
170 integer,
dimension(MPI_STATUS_SIZE) :: st
173 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
177 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
178 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
179 end subroutine ffhd_write_info
197 if(
mype==0)
write(*,*)
'WARNING: set ffhd_thermal_conduction=F when ffhd_energy=F'
201 if(
mype==0)
write(*,*)
'WARNING: set ffhd_hyperbolic_thermal_conduction=F when ffhd_energy=F'
205 if(
mype==0)
write(*,*)
'WARNING: set ffhd_radiative_cooling=F when ffhd_energy=F'
209 if(
mype==0)
write(*,*)
'WARNING: set ffhd_trac=F when ffhd_energy=F'
213 if(
mype==0)
write(*,*)
'WARNING: set ffhd_partial_ionization=F when ffhd_energy=F'
219 if(
mype==0)
write(*,*)
'WARNING: set ffhd_partial_ionization=F when eq_state_units=F'
225 if(
mype==0)
write(*,*)
'WARNING: turn off parabolic TC when using hyperbolic TC'
228 physics_type =
"ffhd"
230 phys_internal_e=.false.
242 if(
mype==0)
write(*,*)
'WARNING: reset ffhd_trac_type=1 for 1D simulation'
247 if(
mype==0)
write(*,*)
'WARNING: set ffhd_trac_mask==bigdouble for global TRAC method'
251 allocate(start_indices(number_species),stop_indices(number_species))
257 mom(:) = var_set_momentum(1)
261 e_ = var_set_energy()
279 stop_indices(1)=nwflux
283 te_ = var_set_auxvar(
'Te',
'Te')
304 if(.not.
allocated(flux_type))
then
305 allocate(flux_type(
ndir, nwflux))
306 flux_type = flux_default
307 else if(any(shape(flux_type) /= [
ndir, nwflux]))
then
308 call mpistop(
"phys_check error: flux_type has wrong shape")
311 phys_get_dt => ffhd_get_dt
312 phys_get_cmax => ffhd_get_cmax_origin
313 phys_get_tcutoff => ffhd_get_tcutoff
314 phys_get_cbounds => ffhd_get_cbounds
315 phys_to_primitive => ffhd_to_primitive_origin
317 phys_to_conserved => ffhd_to_conserved_origin
319 phys_get_flux => ffhd_get_flux
320 phys_get_v => ffhd_get_v_origin
324 phys_add_source_geom => ffhd_add_source_geom
325 phys_add_source => ffhd_add_source
326 phys_check_params => ffhd_check_params
327 phys_write_info => ffhd_write_info
328 phys_handle_small_values => ffhd_handle_small_values_origin
329 ffhd_handle_small_values => ffhd_handle_small_values_origin
330 phys_check_w => ffhd_check_w_origin
333 phys_get_pthermal => ffhd_get_pthermal_iso
336 phys_get_pthermal => ffhd_get_pthermal_origin
342 ffhd_get_rfactor=>rfactor_from_temperature_ionization
343 phys_update_temperature => ffhd_update_temperature
347 ffhd_get_rfactor=>rfactor_from_constant_ionization
357 call ffhd_physical_units()
369 call mpistop(
"thermal conduction needs ffhd_energy=T")
372 call mpistop(
"hyperbolic thermal conduction needs ffhd_energy=T")
375 call mpistop(
"radiative cooling needs ffhd_energy=T")
385 call add_sts_method(ffhd_get_tc_dt_ffhd,ffhd_sts_set_source_tc_ffhd,
e_,1,
e_,1,.false.)
386 tc_fl%get_temperature_from_conserved => ffhd_get_temperature_from_etot
387 tc_fl%get_temperature_from_eint => ffhd_get_temperature_from_eint
402 rc_fl%get_var_Rfactor => ffhd_get_rfactor
405 rc_fl%has_equi = .false.
406 rc_fl%subtract_equi = .false.
413 te_fl_ffhd%get_var_Rfactor => ffhd_get_rfactor
414 phys_te_images => ffhd_te_images
427 subroutine ffhd_te_images
432 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
434 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
436 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
438 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
441 call mpistop(
"Error in synthesize emission: Unknown convert_type")
443 end subroutine ffhd_te_images
446 subroutine ffhd_sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
450 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
451 double precision,
intent(in) :: x(ixi^s,1:
ndim)
452 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
453 double precision,
intent(in) :: my_dt
454 logical,
intent(in) :: fix_conserve_at_step
456 end subroutine ffhd_sts_set_source_tc_ffhd
458 function ffhd_get_tc_dt_ffhd(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 ffhd_get_tc_dt_ffhd
473 subroutine ffhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
476 integer,
intent(in) :: ixi^
l,ixo^
l
477 double precision,
intent(inout) :: w(ixi^s,1:nw)
478 double precision,
intent(in) :: x(ixi^s,1:
ndim)
479 integer,
intent(in) :: step
480 character(len=140) :: error_msg
482 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
483 call ffhd_handle_small_ei(w,x,ixi^
l,ixo^
l,
e_,error_msg)
484 end subroutine ffhd_tc_handle_small_e
486 subroutine tc_params_read_ffhd(fl)
488 type(tc_fluid),
intent(inout) :: fl
491 logical :: tc_saturate=.false.
492 double precision :: tc_k_para=0d0
493 character(len=std_len) :: tc_slope_limiter=
"MC"
495 namelist /tc_list/ tc_saturate, tc_slope_limiter, tc_k_para
499 read(
unitpar, tc_list,
end=111)
503 fl%tc_saturate = tc_saturate
504 fl%tc_k_para = tc_k_para
505 select case(tc_slope_limiter)
507 fl%tc_slope_limiter = 0
510 fl%tc_slope_limiter = 1
513 fl%tc_slope_limiter = 2
516 fl%tc_slope_limiter = 3
519 fl%tc_slope_limiter = 4
521 call mpistop(
"Unknown tc_slope_limiter, choose MC, minmod")
523 end subroutine tc_params_read_ffhd
525 subroutine rc_params_read(fl)
528 type(rc_fluid),
intent(inout) :: fl
530 integer :: ncool = 4000
533 character(len=std_len) :: coolcurve=
'JCcorona'
536 logical :: tfix=.false.
542 logical :: rc_split=.false.
543 logical :: rad_cut=.false.
544 double precision :: rad_cut_hgt=0.5d0
545 double precision :: rad_cut_dey=0.15d0
547 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split, rad_cut, rad_cut_hgt, rad_cut_dey
551 read(
unitpar, rc_list,
end=111)
556 fl%coolcurve=coolcurve
561 fl%rad_cut_hgt=rad_cut_hgt
562 fl%rad_cut_dey=rad_cut_dey
563 end subroutine rc_params_read
565 subroutine ffhd_check_params
572 if (
ffhd_gamma <= 0.0d0)
call mpistop (
"Error: ffhd_gamma <= 0")
573 if (
ffhd_adiab < 0.0d0)
call mpistop (
"Error: ffhd_adiab < 0")
577 call mpistop (
"Error: ffhd_gamma <= 0 or ffhd_gamma == 1")
578 inv_gamma_1=1.d0/gamma_1
583 call mpistop(
"usr_set_equi_vars has to be implemented in the user file")
585 end subroutine ffhd_check_params
587 subroutine ffhd_physical_units()
589 double precision :: mp,kb
590 double precision :: a,b
677 end subroutine ffhd_physical_units
679 subroutine ffhd_check_w_origin(primitive,ixI^L,ixO^L,w,flag)
681 logical,
intent(in) :: primitive
682 integer,
intent(in) :: ixi^
l, ixo^
l
683 double precision,
intent(in) :: w(ixi^s,nw)
684 double precision :: tmp(ixi^s)
685 logical,
intent(inout) :: flag(ixi^s,1:nw)
695 where(tmp(ixo^s) < small_e) flag(ixo^s,
e_) = .true.
698 end subroutine ffhd_check_w_origin
700 subroutine ffhd_to_conserved_origin(ixI^L,ixO^L,w,x)
702 integer,
intent(in) :: ixi^
l, ixo^
l
703 double precision,
intent(inout) :: w(ixi^s, nw)
704 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
707 w(ixo^s,
e_)=w(ixo^s,
p_)*inv_gamma_1+half*w(ixo^s,
mom(1))**2*w(ixo^s,
rho_)
710 end subroutine ffhd_to_conserved_origin
712 subroutine ffhd_to_primitive_origin(ixI^L,ixO^L,w,x)
714 integer,
intent(in) :: ixi^
l, ixo^
l
715 double precision,
intent(inout) :: w(ixi^s, nw)
716 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
719 call ffhd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'ffhd_to_primitive_origin')
722 w(ixo^s,
mom(1)) = w(ixo^s,
mom(1))/w(ixo^s,
rho_)
724 w(ixo^s,
p_)=gamma_1*(w(ixo^s,
e_)-half*w(ixo^s,
rho_)*w(ixo^s,
mom(1))**2)
726 end subroutine ffhd_to_primitive_origin
730 integer,
intent(in) :: ixi^
l, ixo^
l
731 double precision,
intent(inout) :: w(ixi^s, nw)
732 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
739 integer,
intent(in) :: ixi^
l, ixo^
l
740 double precision,
intent(inout) :: w(ixi^s, nw)
741 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
745 call ffhd_handle_small_ei(w,x,ixi^
l,ixi^
l,
e_,
'ffhd_e_to_ei')
749 subroutine ffhd_handle_small_values_origin(primitive, w, x, ixI^L, ixO^L, subname)
752 logical,
intent(in) :: primitive
753 integer,
intent(in) :: ixi^
l,ixo^
l
754 double precision,
intent(inout) :: w(ixi^s,1:nw)
755 double precision,
intent(in) :: x(ixi^s,1:
ndim)
756 character(len=*),
intent(in) :: subname
758 logical :: flag(ixi^s,1:nw)
759 double precision :: tmp2(ixi^s)
761 call phys_check_w(primitive, ixi^
l, ixi^
l, w, flag)
768 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(1)) = 0.0d0
774 where(flag(ixo^s,
e_))
791 if(.not.primitive)
then
800 end subroutine ffhd_handle_small_values_origin
802 subroutine ffhd_get_v_origin(w,x,ixI^L,ixO^L,v)
804 integer,
intent(in) :: ixi^
l, ixo^
l
805 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
806 double precision,
intent(out) :: v(ixi^s,
ndir)
807 double precision :: rho(ixi^s)
811 rho(ixo^s)=1.d0/rho(ixo^s)
813 v(ixo^s,
ndir) = w(ixo^s,
mom(1))*
block%B0(ixo^s,idir,0)*rho(ixo^s)
815 end subroutine ffhd_get_v_origin
819 integer,
intent(in) :: ixi^
l, ixo^
l, idim
820 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
821 double precision,
intent(out) :: v(ixi^s)
822 double precision :: rho(ixi^s)
825 v(ixo^s) = (w(ixo^s,
mom(1))*
block%B0(ixo^s,idim,0)) / rho(ixo^s)
828 subroutine ffhd_get_cmax_origin(wprim,x,ixI^L,ixO^L,idim,cmax)
830 integer,
intent(in) :: ixi^
l, ixo^
l, idim
832 double precision,
intent(in) :: wprim(ixi^s, nw), x(ixi^s,1:
ndim)
833 double precision,
intent(inout) :: cmax(ixi^s)
840 cmax(ixo^s)=dabs(wprim(ixo^s,
mom(1))*
block%B0(ixo^s,idim,0))+cmax(ixo^s)
842 end subroutine ffhd_get_cmax_origin
844 subroutine ffhd_get_tcutoff(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
847 integer,
intent(in) :: ixi^
l,ixo^
l
848 double precision,
intent(in) :: x(ixi^s,1:
ndim)
849 double precision,
intent(inout) :: w(ixi^s,1:nw)
850 double precision,
intent(out) :: tco_local,tmax_local
851 double precision,
parameter :: trac_delta=0.25d0
852 double precision :: te(ixi^s),lts(ixi^s)
853 double precision,
dimension(ixI^S,1:ndim) :: gradt
854 double precision :: bdir(
ndim)
855 double precision :: ltrc,ltrp,altr
856 integer :: idims,ix^
d,jxo^
l,hxo^
l,ixa^
d,ixb^
d
857 integer :: jxp^
l,hxp^
l,ixp^
l,ixq^
l
858 logical :: lrlt(ixi^s)
862 tmax_local=maxval(te(ixo^s))
872 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
874 where(lts(ixo^s) > trac_delta)
877 if(any(lrlt(ixo^s)))
then
878 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
889 lts(ixp^s)=0.5d0*dabs(te(jxp^s)-te(hxp^s))/te(ixp^s)
890 lts(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
891 lts(ixo^s)=0.25d0*(lts(jxo^s)+two*lts(ixo^s)+lts(hxo^s))
892 block%wextra(ixo^s,
tcoff_)=te(ixo^s)*lts(ixo^s)**0.4d0
894 call mpistop(
"ffhd_trac_type not allowed for 1D simulation")
909 call gradient(te,ixi^
l,ixo^
l,idims,gradt(ixi^s,idims))
915 ixb^
d=(ixomin^
d+ixomax^
d-1)/2+ixa^
d;
918 if(sum(bdir(:)**2) .gt. zero)
then
919 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
921 block%special_values(3:ndim+2)=bdir(1:ndim)
923 {
do ix^db=ixomin^db,ixomax^db\}
926 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2))*&
927 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
931 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2),block%ds(ix^d,3))*&
932 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
934 if(lts(ix^d)>trac_delta)
then
935 block%special_values(1)=max(block%special_values(1),te(ix^d))
938 block%special_values(2)=tmax_local
956 call gradient(te,ixi^l,ixq^l,idims,gradt(ixi^s,idims))
957 call gradientf(te,x,ixi^l,hxp^l,idims,gradt(ixi^s,idims),nghostcells,.true.)
958 call gradientf(te,x,ixi^l,jxp^l,idims,gradt(ixi^s,idims),nghostcells,.false.)
960 {
do ix^db=ixpmin^db,ixpmax^db\}
962 lts(ix^d)=abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
964 lts(ix^d)=min(^d&block%ds({ix^d},^d))*lts(ix^d)
965 lts(ix^d)=max(one,(exp(lts(ix^d))/ltrc)**ltrp)
970 {
do ix^db=ixpmin^db,ixpmax^db\}
972 altr=0.25d0*((lts(ix1-1,ix2)+two*lts(ix^d)+lts(ix1+1,ix2))*block%B0(ix^d,1,0)**2+&
973 (lts(ix1,ix2-1)+two*lts(ix^d)+lts(ix1,ix2+1))*block%B0(ix^d,2,0)**2)
974 block%wextra(ix^d,
tcoff_)=te(ix^d)*altr**0.4d0
977 altr=0.25d0*((lts(ix1-1,ix2,ix3)+two*lts(ix^d)+lts(ix1+1,ix2,ix3))*block%B0(ix^d,1,0)**2+&
978 (lts(ix1,ix2-1,ix3)+two*lts(ix^d)+lts(ix1,ix2+1,ix3))*block%B0(ix^d,2,0)**2+&
979 (lts(ix1,ix2,ix3-1)+two*lts(ix^d)+lts(ix1,ix2,ix3+1))*block%B0(ix^d,3,0)**2)
980 block%wextra(ix^d,
tcoff_)=te(ix^d)*altr**0.4d0
986 call mpistop(
"unknown ffhd_trac_type")
989 end subroutine ffhd_get_tcutoff
991 subroutine ffhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
993 integer,
intent(in) :: ixi^
l, ixo^
l, idim
994 double precision,
intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
995 double precision,
intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
996 double precision,
intent(in) :: x(ixi^s,1:
ndim)
997 double precision,
intent(inout) :: cmax(ixi^s,1:number_species)
998 double precision,
intent(inout),
optional :: cmin(ixi^s,1:number_species)
999 double precision,
intent(in) :: hspeed(ixi^s,1:number_species)
1000 double precision :: wmean(ixi^s,nw)
1001 double precision,
dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1007 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1008 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1009 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1010 umean(ixo^s)=(wlp(ixo^s,
mom(1))*tmp1(ixo^s)&
1011 +wrp(ixo^s,
mom(1))*tmp2(ixo^s))*tmp3(ixo^s)
1012 umean(ixo^s)=umean(ixo^s)*
block%B0(ixo^s,idim,idim)
1015 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1016 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1017 ((wrp(ixo^s,
mom(1))-wlp(ixo^s,
mom(1)))*
block%B0(ixo^s,idim,idim))**2
1018 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1019 if(
present(cmin))
then
1020 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1021 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1023 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1026 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1027 tmp1(ixo^s)=wmean(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim)/wmean(ixo^s,
rho_)
1029 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1030 if(
present(cmin))
then
1031 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1032 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1034 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1040 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1041 if(
present(cmin))
then
1042 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1043 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))-csoundl(ixo^s)
1044 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1045 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1047 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1048 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1051 end subroutine ffhd_get_cbounds
1053 subroutine ffhd_get_pthermal_iso(w,x,ixI^L,ixO^L,pth)
1056 integer,
intent(in) :: ixi^
l, ixo^
l
1057 double precision,
intent(in) :: w(ixi^s,nw)
1058 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1059 double precision,
intent(out):: pth(ixi^s)
1063 end subroutine ffhd_get_pthermal_iso
1065 subroutine ffhd_get_pthermal_origin(w,x,ixI^L,ixO^L,pth)
1068 integer,
intent(in) :: ixi^
l, ixo^
l
1069 double precision,
intent(in) :: w(ixi^s,nw)
1070 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1071 double precision,
intent(out):: pth(ixi^s)
1076 {
do ix^db= ixo^lim^db\}
1081 elseif(check_small_values)
then
1082 {
do ix^db= ixo^lim^db\}
1083 if(pth(ix^d)<small_pressure)
then
1084 write(*,*)
"Error: small value of gas pressure",pth(ix^d),&
1085 " encountered when call ffhd_get_pthermal"
1086 write(*,*)
"Iteration: ", it,
" Time: ", global_time
1087 write(*,*)
"Location: ", x(ix^d,:)
1088 write(*,*)
"Cell number: ", ix^d
1090 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^d,iw)
1093 if(trace_small_values)
write(*,*) dsqrt(pth(ix^d)-bigdouble)
1094 write(*,*)
"Saving status at the previous time step"
1099 end subroutine ffhd_get_pthermal_origin
1101 subroutine ffhd_get_temperature_from_te(w, x, ixI^L, ixO^L, res)
1103 integer,
intent(in) :: ixi^
l, ixo^
l
1104 double precision,
intent(in) :: w(ixi^s, 1:nw)
1105 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1106 double precision,
intent(out):: res(ixi^s)
1108 res(ixo^s) = w(ixo^s,
te_)
1109 end subroutine ffhd_get_temperature_from_te
1111 subroutine ffhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1113 integer,
intent(in) :: ixi^
l, ixo^
l
1114 double precision,
intent(in) :: w(ixi^s, 1:nw)
1115 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1116 double precision,
intent(out):: res(ixi^s)
1117 double precision :: r(ixi^s)
1119 call ffhd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1120 res(ixo^s) = gamma_1 * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1121 end subroutine ffhd_get_temperature_from_eint
1123 subroutine ffhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1125 integer,
intent(in) :: ixi^
l, ixo^
l
1126 double precision,
intent(in) :: w(ixi^s, 1:nw)
1127 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1128 double precision,
intent(out):: res(ixi^s)
1130 double precision :: r(ixi^s)
1132 call ffhd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1134 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1135 end subroutine ffhd_get_temperature_from_etot
1139 integer,
intent(in) :: ixi^
l, ixo^
l
1140 double precision,
intent(in) :: w(ixi^s,nw)
1141 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1142 double precision,
intent(out) :: csound2(ixi^s)
1143 double precision :: rho(ixi^s)
1148 csound2(ixo^s)=
ffhd_gamma*csound2(ixo^s)/rho(ixo^s)
1154 subroutine ffhd_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
1157 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1159 double precision,
intent(in) :: wc(ixi^s,nw)
1161 double precision,
intent(in) :: w(ixi^s,nw)
1162 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1163 double precision,
intent(out) :: f(ixi^s,nwflux)
1164 double precision :: ptotal(ixo^s)
1169 ptotal(ixo^s)=w(ixo^s,
p_)
1175 f(ixo^s,
mom(1))=(wc(ixo^s,
mom(1))*w(ixo^s,
mom(1))+ptotal(ixo^s))*
block%B0(ixo^s,idim,idim)
1179 f(ixo^s,
e_)=w(ixo^s,
mom(1))*(wc(ixo^s,
e_)+ptotal(ixo^s))*
block%B0(ixo^s,idim,idim)
1181 f(ixo^s,
e_)=f(ixo^s,
e_)+w(ixo^s,
q_)*
block%B0(ixo^s,idim,idim)
1185 end subroutine ffhd_get_flux
1187 subroutine ffhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1191 integer,
intent(in) :: ixi^
l, ixo^
l
1192 double precision,
intent(in) :: qdt,dtfactor
1193 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:
ndim)
1194 double precision,
intent(inout) :: w(ixi^s,1:nw)
1195 logical,
intent(in) :: qsourcesplit
1196 logical,
intent(inout) :: active
1198 if (.not. qsourcesplit)
then
1200 call add_punitb(qdt,ixi^
l,ixo^
l,wct,w,x,wctprim)
1202 call add_hypertc_source(qdt,ixi^
l,ixo^
l,wct,w,x,wctprim)
1208 w,x,qsourcesplit,active,
rc_fl)
1218 if(.not.qsourcesplit)
then
1220 call ffhd_update_temperature(ixi^
l,ixo^
l,wct,w,x)
1223 end subroutine ffhd_add_source
1225 subroutine add_punitb(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1228 integer,
intent(in) :: ixi^
l,ixo^
l
1229 double precision,
intent(in) :: qdt
1230 double precision,
intent(in) :: wct(ixi^s,1:nw),x(ixi^s,1:
ndim)
1231 double precision,
intent(in) :: wctprim(ixi^s,1:nw)
1232 double precision,
intent(inout) :: w(ixi^s,1:nw)
1234 integer :: idims,hxo^
l
1235 double precision :: divb(ixi^s)
1236 double precision :: rhovpar(ixi^s),gradrhov(ixi^s)
1241 hxo^
l=ixo^
l-
kr(idims,^
d);
1242 divb(ixo^s)=divb(ixo^s)+(
block%B0(ixo^s,idims,idims)-
block%B0(hxo^s,idims,idims))/
dxlevel(idims)
1247 w(ixo^s,
mom(1))=w(ixo^s,
mom(1))+qdt*wctprim(ixo^s,
p_)*divb(ixo^s)
1248 end subroutine add_punitb
1252 integer,
intent(in) :: ixi^
l, ixo^
l
1253 double precision,
intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:
ndim)
1254 double precision,
intent(out) :: rho(ixi^s)
1256 rho(ixo^s) = w(ixo^s,
rho_)
1259 subroutine ffhd_handle_small_ei(w, x, ixI^L, ixO^L, ie, subname)
1262 integer,
intent(in) :: ixi^
l,ixo^
l, ie
1263 double precision,
intent(inout) :: w(ixi^s,1:nw)
1264 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1265 character(len=*),
intent(in) :: subname
1267 logical :: flag(ixi^s,1:nw)
1268 double precision :: rho(ixi^s)
1271 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
1272 if(any(flag(ixo^s,ie)))
then
1275 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
1279 w(ixo^s,
e_)=w(ixo^s,
e_)*gamma_1
1281 w(ixo^s,
mom(1)) = w(ixo^s,
mom(1))/rho(ixo^s)
1285 end subroutine ffhd_handle_small_ei
1287 subroutine ffhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1290 integer,
intent(in) :: ixi^
l, ixo^
l
1291 double precision,
intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:
ndim)
1292 double precision,
intent(inout) :: w(ixi^s,1:nw)
1293 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1299 end subroutine ffhd_update_temperature
1301 subroutine ffhd_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
1305 integer,
intent(in) :: ixi^
l, ixo^
l
1306 double precision,
intent(inout) :: dtnew
1307 double precision,
intent(in) ::
dx^
d
1308 double precision,
intent(in) :: w(ixi^s,1:nw)
1309 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1316 end subroutine ffhd_get_dt
1318 subroutine ffhd_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
1320 integer,
intent(in) :: ixi^
l, ixo^
l
1321 double precision,
intent(in) :: qdt, dtfactor,x(ixi^s,1:
ndim)
1322 double precision,
intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw),w(ixi^s,1:nw)
1325 end subroutine ffhd_add_source_geom
1327 function ffhd_kin_en_origin(w, ixI^L, ixO^L, inv_rho)
result(ke)
1329 integer,
intent(in) :: ixi^
l, ixo^
l
1330 double precision,
intent(in) :: w(ixi^s, nw)
1331 double precision :: ke(ixo^s)
1332 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1334 if(
present(inv_rho))
then
1335 ke(ixo^s)=0.5d0*w(ixo^s,
mom(1))**2*inv_rho(ixo^s)
1337 ke(ixo^s)=0.5d0*w(ixo^s,
mom(1))**2/w(ixo^s,
rho_)
1339 end function ffhd_kin_en_origin
1341 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1344 integer,
intent(in) :: ixi^
l, ixo^
l
1345 double precision,
intent(in) :: w(ixi^s,1:nw)
1346 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1347 double precision,
intent(out):: rfactor(ixi^s)
1348 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1351 rfactor(ixo^s)=(1.d0+iz_h(ixo^s)+0.1d0*(1.d0+iz_he(ixo^s)*(1.d0+iz_he(ixo^s))))/(2.d0+3.d0*
he_abundance)
1352 end subroutine rfactor_from_temperature_ionization
1354 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1356 integer,
intent(in) :: ixi^
l, ixo^
l
1357 double precision,
intent(in) :: w(ixi^s,1:nw)
1358 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1359 double precision,
intent(out):: rfactor(ixi^s)
1362 end subroutine rfactor_from_constant_ionization
1364 subroutine add_hypertc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1367 integer,
intent(in) :: ixi^
l,ixo^
l
1368 double precision,
intent(in) :: qdt
1369 double precision,
dimension(ixI^S,1:ndim),
intent(in) :: x
1370 double precision,
dimension(ixI^S,1:nw),
intent(in) :: wct,wctprim
1371 double precision,
dimension(ixI^S,1:nw),
intent(inout) :: w
1373 double precision,
dimension(ixI^S) :: te
1374 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1377 te(ixi^s)=wctprim(ixi^s,
p_)/wct(ixi^s,
rho_)
1381 do ix2=ixomin2,ixomax2
1382 do ix1=ixomin1,ixomax1
1389 sigma_t7=sigma_t5*te(ix^
d)
1393 sigma_t7=sigma_t5*te(ix^
d)
1395 sigmat5_bgradt=sigma_t5*(&
1396 block%B0(ix^
d,1,0)*((8.d0*(te(ix1+1,ix2)-te(ix1-1,ix2))-te(ix1+2,ix2)+te(ix1-2,ix2))/12.d0)/
block%ds(ix^
d,1)&
1397 +
block%B0(ix^
d,2,0)*((8.d0*(te(ix1,ix2+1)-te(ix1,ix2-1))-te(ix1,ix2+2)+te(ix1,ix2-2))/12.d0)/
block%ds(ix^
d,2))
1400 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^
d,
rho_)*(wctprim(ix^
d,
p_)/wct(ix^
d,
rho_))**1.5d0))
1402 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^
d,
q_))/tau
1404 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(sigmat5_bgradt+wct(ix^
d,
q_))/&
1411 do ix3=ixomin3,ixomax3
1412 do ix2=ixomin2,ixomax2
1413 do ix1=ixomin1,ixomax1
1420 sigma_t7=sigma_t5*te(ix^
d)
1424 sigma_t7=sigma_t5*te(ix^
d)
1426 sigmat5_bgradt=sigma_t5*(&
1427 block%B0(ix^
d,1,0)*((8.d0*(te(ix1+1,ix2,ix3)-te(ix1-1,ix2,ix3))-te(ix1+2,ix2,ix3)+te(ix1-2,ix2,ix3))/12.d0)/
block%ds(ix^
d,1)&
1428 +
block%B0(ix^
d,2,0)*((8.d0*(te(ix1,ix2+1,ix3)-te(ix1,ix2-1,ix3))-te(ix1,ix2+2,ix3)+te(ix1,ix2-2,ix3))/12.d0)/
block%ds(ix^
d,2)&
1429 +
block%B0(ix^
d,3,0)*((8.d0*(te(ix1,ix2,ix3+1)-te(ix1,ix2,ix3-1))-te(ix1,ix2,ix3+2)+te(ix1,ix2,ix3-2))/12.d0)/
block%ds(ix^
d,3))
1432 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^
d,
rho_)*(wctprim(ix^
d,
p_)/wct(ix^
d,
rho_))**1.5d0))
1434 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^
d,
q_))/tau
1436 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(sigmat5_bgradt+wct(ix^
d,
q_))/&
1443 end subroutine add_hypertc_source
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.
subroutine add_convert_method(phys_convert_vars, nwc, dataset_names, file_suffix)
Frozen-field hydrodynamics module.
integer, public, protected te_
Indices of temperature.
integer, public, protected ffhd_trac_type
Which TRAC method is used.
double precision, public hypertc_kappa
The thermal conductivity kappa in hyperbolic thermal conduction.
integer, public, protected e_
Index of the energy density (-1 if not present)
double precision, public ffhd_gamma
The adiabatic index.
logical, public, protected eq_state_units
double precision, public ffhd_adiab
The adiabatic constant.
procedure(sub_get_pthermal), pointer, public ffhd_get_temperature
logical, public, protected ffhd_hyperbolic_thermal_conduction
Whether hyperbolic type thermal conduction is used.
type(rc_fluid), allocatable, public rc_fl
type of fluid for radiative cooling
integer, public, protected rho_
Index of the density (in the w array)
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
procedure(sub_get_pthermal), pointer, public ffhd_get_pthermal
logical, public, protected ffhd_htc_sat
Wheterh saturation is considered for hyperbolic TC.
subroutine, public ffhd_get_rho(w, x, ixil, ixol, rho)
logical, public, protected ffhd_partial_ionization
Whether plasma is partially ionized.
procedure(sub_convert), pointer, public ffhd_to_conserved
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)
procedure(fun_kin_en), pointer, public ffhd_kin_en
subroutine, public ffhd_phys_init()
procedure(sub_get_v), pointer, public ffhd_get_v
subroutine, public ffhd_get_csound2(w, x, ixil, ixol, csound2)
logical, public, protected ffhd_energy
Whether an energy equation is used.
double precision, public, protected rr
type(tc_fluid), allocatable, public tc_fl
type of fluid for thermal conduction
type(te_fluid), allocatable, public te_fl_ffhd
type of fluid for thermal emission synthesis
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
logical, public, protected ffhd_radiative_cooling
Whether radiative cooling is added.
subroutine, public ffhd_get_v_idim(w, x, ixil, ixol, idim, v)
integer, public, protected q_
procedure(sub_convert), pointer, public ffhd_to_primitive
integer, public, protected tweight_
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
logical, public, protected ffhd_trac
Whether TRAC method is used.
logical, public, protected ffhd_thermal_conduction
Whether thermal conduction is used.
subroutine, public ffhd_ei_to_e(ixil, ixol, w, x)
logical, public, protected ffhd_gravity
Whether gravity is added.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
double precision, public, protected ffhd_trac_mask
Height of the mask used in the TRAC method.
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
integer, public, protected ffhd_trac_finegrid
Distance between two adjacent traced magnetic field lines (in finest cell size)
subroutine, public ffhd_e_to_ei(ixil, ixol, w, x)
Module for flux conservation near refinement boundaries.
Module with geometry-related routines (e.g., divergence, curl)
subroutine divvector(qvec, ixil, ixol, divq, nth_in)
subroutine gradient(q, ixil, ixol, idir, gradq, nth_in)
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.
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 unit_mass
Physical scaling factor for mass.
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision phys_trac_mask
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 cmax_global
global fastest wave speed needed in fd scheme and glm method
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
double precision dt
global time step
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision courantpar
The Courant (CFL) number used for the simulation.
double precision unit_velocity
Physical scaling factor for velocity.
logical b0field
split magnetic field as background B0 field
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 need_global_cmax
need global maximal wave speed
logical fix_small_values
fix small values with average or replace methods
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
double precision small_density
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer phys_trac_finegrid
integer, parameter unitconvert
integer number_equi_vars
number of equilibrium set variables, besides the mag field
Module for including gravity in (magneto)hydrodynamics simulations.
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
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.
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 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_mhd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (mhd implementation) Note: for multi-D MHD (1D MHD will use HD f...
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_mhd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
anisotropic thermal conduction with slope limited symmetric scheme Sharma 2007 Journal of Computation...
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_equi_vars), pointer usr_set_equi_vars