56 integer,
public,
protected ::
rho_
59 integer,
allocatable,
public,
protected ::
mom(:)
62 integer,
public,
protected ::
e_
65 integer,
public,
protected ::
p_
68 integer,
public,
protected ::
te_
73 integer,
public,
protected ::
q_
88 double precision,
public,
protected ::
h_ion_fr=1d0
91 double precision,
public,
protected ::
he_ion_fr=1d0
98 double precision,
public,
protected ::
rr=1d0
105 double precision :: gamma_1, inv_gamma_1
110 function fun_kin_en(w, ixI^L, ixO^L, inv_rho)
result(ke)
112 integer,
intent(in) :: ixi^
l, ixo^
l
113 double precision,
intent(in) :: w(ixi^s, nw)
114 double precision :: ke(ixo^s)
115 double precision,
intent(in),
optional :: inv_rho(ixo^s)
116 end function fun_kin_en
144 subroutine ffhd_read_params(files)
146 character(len=*),
intent(in) :: files(:)
154 do n = 1,
size(files)
155 open(
unitpar, file=trim(files(n)), status=
"old")
156 read(
unitpar, ffhd_list,
end=111)
159 end subroutine ffhd_read_params
162 subroutine ffhd_write_info(fh)
164 integer,
intent(in) :: fh
165 integer,
parameter :: n_par = 1
166 double precision :: values(n_par)
167 character(len=name_len) :: names(n_par)
168 integer,
dimension(MPI_STATUS_SIZE) :: st
171 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
175 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
176 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
177 end subroutine ffhd_write_info
195 if(
mype==0)
write(*,*)
'WARNING: set ffhd_thermal_conduction=F when ffhd_energy=F'
199 if(
mype==0)
write(*,*)
'WARNING: set ffhd_hyperbolic_tc=F when ffhd_energy=F'
203 if(
mype==0)
write(*,*)
'WARNING: set ffhd_radiative_cooling=F when ffhd_energy=F'
207 if(
mype==0)
write(*,*)
'WARNING: set ffhd_trac=F when ffhd_energy=F'
211 if(
mype==0)
write(*,*)
'WARNING: set ffhd_partial_ionization=F when ffhd_energy=F'
217 if(
mype==0)
write(*,*)
'WARNING: set ffhd_partial_ionization=F when eq_state_units=F'
223 if(
mype==0)
write(*,*)
'WARNING: turn off parabolic TC when using hyperbolic TC'
226 physics_type =
"ffhd"
228 phys_internal_e=.false.
240 if(
mype==0)
write(*,*)
'WARNING: reset ffhd_trac_type=1 for 1D simulation'
245 if(
mype==0)
write(*,*)
'WARNING: set ffhd_trac_mask==bigdouble for global TRAC method'
249 allocate(start_indices(number_species),stop_indices(number_species))
255 mom(:) = var_set_momentum(1)
259 e_ = var_set_energy()
267 q_ = var_set_fluxvar(
'q',
'q', need_bc=.false.)
277 stop_indices(1)=nwflux
281 te_ = var_set_auxvar(
'Te',
'Te')
302 if(.not.
allocated(flux_type))
then
303 allocate(flux_type(
ndir, nwflux))
304 flux_type = flux_default
305 else if(any(shape(flux_type) /= [
ndir, nwflux]))
then
306 call mpistop(
"phys_check error: flux_type has wrong shape")
309 phys_get_dt => ffhd_get_dt
310 phys_get_cmax => ffhd_get_cmax_origin
311 phys_get_tcutoff => ffhd_get_tcutoff
312 phys_get_cbounds => ffhd_get_cbounds
313 phys_to_primitive => ffhd_to_primitive_origin
315 phys_to_conserved => ffhd_to_conserved_origin
317 phys_get_flux => ffhd_get_flux
318 phys_get_v => ffhd_get_v_origin
322 phys_add_source_geom => ffhd_add_source_geom
323 phys_add_source => ffhd_add_source
324 phys_check_params => ffhd_check_params
325 phys_write_info => ffhd_write_info
326 phys_handle_small_values => ffhd_handle_small_values_origin
327 ffhd_handle_small_values => ffhd_handle_small_values_origin
328 phys_check_w => ffhd_check_w_origin
331 phys_get_pthermal => ffhd_get_pthermal_iso
334 phys_get_pthermal => ffhd_get_pthermal_origin
340 ffhd_get_rfactor=>rfactor_from_temperature_ionization
341 phys_update_temperature => ffhd_update_temperature
345 ffhd_get_rfactor=>rfactor_from_constant_ionization
355 call ffhd_physical_units()
367 call mpistop(
"thermal conduction needs ffhd_energy=T")
370 call mpistop(
"hyperbolic thermal conduction needs ffhd_energy=T")
373 call mpistop(
"radiative cooling needs ffhd_energy=T")
383 call add_sts_method(ffhd_get_tc_dt_ffhd,ffhd_sts_set_source_tc_ffhd,
e_,1,
e_,1,.false.)
384 tc_fl%get_temperature_from_conserved => ffhd_get_temperature_from_etot
385 tc_fl%get_temperature_from_eint => ffhd_get_temperature_from_eint
400 rc_fl%get_var_Rfactor => ffhd_get_rfactor
403 rc_fl%has_equi = .false.
404 rc_fl%subtract_equi = .false.
411 te_fl_ffhd%get_var_Rfactor => ffhd_get_rfactor
412 phys_te_images => ffhd_te_images
425 subroutine ffhd_te_images
430 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
432 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
434 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
436 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
439 call mpistop(
"Error in synthesize emission: Unknown convert_type")
441 end subroutine ffhd_te_images
444 subroutine ffhd_sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
448 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
449 double precision,
intent(in) :: x(ixi^s,1:
ndim)
450 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
451 double precision,
intent(in) :: my_dt
452 logical,
intent(in) :: fix_conserve_at_step
454 end subroutine ffhd_sts_set_source_tc_ffhd
456 function ffhd_get_tc_dt_ffhd(w,ixI^L,ixO^L,dx^D,x)
result(dtnew)
463 integer,
intent(in) :: ixi^
l, ixo^
l
464 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
465 double precision,
intent(in) :: w(ixi^s,1:nw)
466 double precision :: dtnew
469 end function ffhd_get_tc_dt_ffhd
471 subroutine ffhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
474 integer,
intent(in) :: ixi^
l,ixo^
l
475 double precision,
intent(inout) :: w(ixi^s,1:nw)
476 double precision,
intent(in) :: x(ixi^s,1:
ndim)
477 integer,
intent(in) :: step
478 character(len=140) :: error_msg
480 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
481 call ffhd_handle_small_ei(w,x,ixi^
l,ixo^
l,
e_,error_msg)
482 end subroutine ffhd_tc_handle_small_e
484 subroutine tc_params_read_ffhd(fl)
486 type(tc_fluid),
intent(inout) :: fl
489 logical :: tc_saturate=.false.
490 double precision :: tc_k_para=0d0
491 character(len=std_len) :: tc_slope_limiter=
"MC"
493 namelist /tc_list/ tc_saturate, tc_slope_limiter, tc_k_para
497 read(
unitpar, tc_list,
end=111)
501 fl%tc_saturate = tc_saturate
502 fl%tc_k_para = tc_k_para
503 select case(tc_slope_limiter)
505 fl%tc_slope_limiter = 0
508 fl%tc_slope_limiter = 1
511 fl%tc_slope_limiter = 2
514 fl%tc_slope_limiter = 3
517 fl%tc_slope_limiter = 4
519 call mpistop(
"Unknown tc_slope_limiter, choose MC, minmod")
521 end subroutine tc_params_read_ffhd
523 subroutine rc_params_read(fl)
526 type(rc_fluid),
intent(inout) :: fl
528 integer :: ncool = 4000
531 character(len=std_len) :: coolcurve=
'JCcorona'
534 logical :: tfix=.false.
540 logical :: rc_split=.false.
541 logical :: rad_cut=.false.
542 double precision :: rad_cut_hgt=0.5d0
543 double precision :: rad_cut_dey=0.15d0
545 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split, rad_cut, rad_cut_hgt, rad_cut_dey
549 read(
unitpar, rc_list,
end=111)
554 fl%coolcurve=coolcurve
559 fl%rad_cut_hgt=rad_cut_hgt
560 fl%rad_cut_dey=rad_cut_dey
561 end subroutine rc_params_read
563 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")
588 write(*,*)
'====FFHD run with settings===================='
589 write(*,*)
'Using mod_ffhd_phys with settings:'
591 write(*,*)
'Dimensionality :',
ndim
592 write(*,*)
'vector components:',
ndir
594 write(*,*)
'number of variables nw=',nw
595 write(*,*)
' start index iwstart=',iwstart
596 write(*,*)
'number of vector variables=',nvector
597 write(*,*)
'number of stagger variables nws=',nws
598 write(*,*)
'number of variables with BCs=',nwgc
599 write(*,*)
'number of vars with fluxes=',nwflux
600 write(*,*)
'number of vars with flux + BC=',nwfluxbc
601 write(*,*)
'number of auxiliary variables=',nwaux
602 write(*,*)
'number of extra vars without flux=',nwextra
603 write(*,*)
'number of extra vars for wextra=',nw_extra
604 write(*,*)
'number of auxiliary I/O variables=',
nwauxio
611 write(*,*)
'number due to phys_wider_stencil=',phys_wider_stencil
612 write(*,*)
'==========================================='
614 end subroutine ffhd_check_params
616 subroutine ffhd_physical_units()
618 double precision :: mp,kb
619 double precision :: a,b
706 end subroutine ffhd_physical_units
708 subroutine ffhd_check_w_origin(primitive,ixI^L,ixO^L,w,flag)
710 logical,
intent(in) :: primitive
711 integer,
intent(in) :: ixi^
l, ixo^
l
712 double precision,
intent(in) :: w(ixi^s,nw)
713 double precision :: tmp(ixi^s)
714 logical,
intent(inout) :: flag(ixi^s,1:nw)
724 where(tmp(ixo^s) <
small_e) flag(ixo^s,
e_) = .true.
727 end subroutine ffhd_check_w_origin
729 subroutine ffhd_to_conserved_origin(ixI^L,ixO^L,w,x)
731 integer,
intent(in) :: ixi^
l, ixo^
l
732 double precision,
intent(inout) :: w(ixi^s, nw)
733 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
736 w(ixo^s,
e_)=w(ixo^s,
p_)*inv_gamma_1+half*w(ixo^s,
mom(1))**2*w(ixo^s,
rho_)
739 end subroutine ffhd_to_conserved_origin
741 subroutine ffhd_to_primitive_origin(ixI^L,ixO^L,w,x)
743 integer,
intent(in) :: ixi^
l, ixo^
l
744 double precision,
intent(inout) :: w(ixi^s, nw)
745 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
748 call ffhd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'ffhd_to_primitive_origin')
751 w(ixo^s,
mom(1)) = w(ixo^s,
mom(1))/w(ixo^s,
rho_)
753 w(ixo^s,
p_)=gamma_1*(w(ixo^s,
e_)-half*w(ixo^s,
rho_)*w(ixo^s,
mom(1))**2)
755 end subroutine ffhd_to_primitive_origin
759 integer,
intent(in) :: ixi^
l, ixo^
l
760 double precision,
intent(inout) :: w(ixi^s, nw)
761 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
768 integer,
intent(in) :: ixi^
l, ixo^
l
769 double precision,
intent(inout) :: w(ixi^s, nw)
770 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
774 call ffhd_handle_small_ei(w,x,ixi^
l,ixi^
l,
e_,
'ffhd_e_to_ei')
778 subroutine ffhd_handle_small_values_origin(primitive, w, x, ixI^L, ixO^L, subname)
781 logical,
intent(in) :: primitive
782 integer,
intent(in) :: ixi^
l,ixo^
l
783 double precision,
intent(inout) :: w(ixi^s,1:nw)
784 double precision,
intent(in) :: x(ixi^s,1:
ndim)
785 character(len=*),
intent(in) :: subname
787 logical :: flag(ixi^s,1:nw)
788 double precision :: tmp2(ixi^s)
790 call phys_check_w(primitive, ixi^
l, ixi^
l, w, flag)
797 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(1)) = 0.0d0
803 where(flag(ixo^s,
e_))
820 if(.not.primitive)
then
829 end subroutine ffhd_handle_small_values_origin
831 subroutine ffhd_get_v_origin(w,x,ixI^L,ixO^L,v)
833 integer,
intent(in) :: ixi^
l, ixo^
l
834 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
835 double precision,
intent(out) :: v(ixi^s,
ndir)
836 double precision :: rho(ixi^s)
840 rho(ixo^s)=1.d0/rho(ixo^s)
842 v(ixo^s,
ndir) = w(ixo^s,
mom(1))*
block%B0(ixo^s,idir,0)*rho(ixo^s)
844 end subroutine ffhd_get_v_origin
848 integer,
intent(in) :: ixi^
l, ixo^
l, idim
849 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
850 double precision,
intent(out) :: v(ixi^s)
851 double precision :: rho(ixi^s)
854 v(ixo^s) = (w(ixo^s,
mom(1))*
block%B0(ixo^s,idim,0)) / rho(ixo^s)
857 subroutine ffhd_get_cmax_origin(wprim,x,ixI^L,ixO^L,idim,cmax)
859 integer,
intent(in) :: ixi^
l, ixo^
l, idim
861 double precision,
intent(in) :: wprim(ixi^s, nw), x(ixi^s,1:
ndim)
862 double precision,
intent(inout) :: cmax(ixi^s)
869 cmax(ixo^s)=dabs(wprim(ixo^s,
mom(1))*
block%B0(ixo^s,idim,0))+cmax(ixo^s)
871 end subroutine ffhd_get_cmax_origin
873 subroutine ffhd_get_tcutoff(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
876 integer,
intent(in) :: ixi^
l,ixo^
l
877 double precision,
intent(in) :: x(ixi^s,1:
ndim)
878 double precision,
intent(inout) :: w(ixi^s,1:nw)
879 double precision,
intent(out) :: tco_local,tmax_local
880 double precision,
parameter :: trac_delta=0.25d0
881 double precision :: te(ixi^s),lts(ixi^s)
882 double precision,
dimension(ixI^S,1:ndim) :: gradt
883 double precision :: bdir(
ndim)
884 double precision :: ltrc,ltrp,altr
885 integer :: idims,ix^
d,jxo^
l,hxo^
l,ixa^
d,ixb^
d
886 integer :: jxp^
l,hxp^
l,ixp^
l,ixq^
l
887 logical :: lrlt(ixi^s)
891 tmax_local=maxval(te(ixo^s))
901 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
903 where(lts(ixo^s) > trac_delta)
906 if(any(lrlt(ixo^s)))
then
907 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
918 lts(ixp^s)=0.5d0*dabs(te(jxp^s)-te(hxp^s))/te(ixp^s)
919 lts(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
920 lts(ixo^s)=0.25d0*(lts(jxo^s)+two*lts(ixo^s)+lts(hxo^s))
921 block%wextra(ixo^s,
tcoff_)=te(ixo^s)*lts(ixo^s)**0.4d0
923 call mpistop(
"ffhd_trac_type not allowed for 1D simulation")
938 call gradient(te,ixi^
l,ixo^
l,idims,gradt(ixi^s,idims))
944 ixb^
d=(ixomin^
d+ixomax^
d-1)/2+ixa^
d;
947 if(sum(bdir(:)**2) .gt. zero)
then
948 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
950 block%special_values(3:ndim+2)=bdir(1:ndim)
952 {
do ix^db=ixomin^db,ixomax^db\}
955 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2))*&
956 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
960 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2),block%ds(ix^d,3))*&
961 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
963 if(lts(ix^d)>trac_delta)
then
964 block%special_values(1)=max(block%special_values(1),te(ix^d))
967 block%special_values(2)=tmax_local
985 call gradient(te,ixi^l,ixq^l,idims,gradt(ixi^s,idims))
986 call gradientf(te,x,ixi^l,hxp^l,idims,gradt(ixi^s,idims),nghostcells,.true.)
987 call gradientf(te,x,ixi^l,jxp^l,idims,gradt(ixi^s,idims),nghostcells,.false.)
989 {
do ix^db=ixpmin^db,ixpmax^db\}
991 lts(ix^d)=abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
993 lts(ix^d)=min(^d&block%ds({ix^d},^d))*lts(ix^d)
994 lts(ix^d)=max(one,(exp(lts(ix^d))/ltrc)**ltrp)
999 {
do ix^db=ixpmin^db,ixpmax^db\}
1001 altr=0.25d0*((lts(ix1-1,ix2)+two*lts(ix^d)+lts(ix1+1,ix2))*block%B0(ix^d,1,0)**2+&
1002 (lts(ix1,ix2-1)+two*lts(ix^d)+lts(ix1,ix2+1))*block%B0(ix^d,2,0)**2)
1003 block%wextra(ix^d,
tcoff_)=te(ix^d)*altr**0.4d0
1006 altr=0.25d0*((lts(ix1-1,ix2,ix3)+two*lts(ix^d)+lts(ix1+1,ix2,ix3))*block%B0(ix^d,1,0)**2+&
1007 (lts(ix1,ix2-1,ix3)+two*lts(ix^d)+lts(ix1,ix2+1,ix3))*block%B0(ix^d,2,0)**2+&
1008 (lts(ix1,ix2,ix3-1)+two*lts(ix^d)+lts(ix1,ix2,ix3+1))*block%B0(ix^d,3,0)**2)
1009 block%wextra(ix^d,
tcoff_)=te(ix^d)*altr**0.4d0
1015 call mpistop(
"unknown ffhd_trac_type")
1018 end subroutine ffhd_get_tcutoff
1020 subroutine ffhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
1022 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1023 double precision,
intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1024 double precision,
intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1025 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1026 double precision,
intent(inout) :: cmax(ixi^s,1:number_species)
1027 double precision,
intent(inout),
optional :: cmin(ixi^s,1:number_species)
1028 double precision,
intent(in) :: hspeed(ixi^s,1:number_species)
1029 double precision :: wmean(ixi^s,nw)
1030 double precision,
dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1036 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1037 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1038 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1039 umean(ixo^s)=(wlp(ixo^s,
mom(1))*tmp1(ixo^s)&
1040 +wrp(ixo^s,
mom(1))*tmp2(ixo^s))*tmp3(ixo^s)
1041 umean(ixo^s)=umean(ixo^s)*
block%B0(ixo^s,idim,idim)
1044 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1045 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1046 ((wrp(ixo^s,
mom(1))-wlp(ixo^s,
mom(1)))*
block%B0(ixo^s,idim,idim))**2
1047 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1048 if(
present(cmin))
then
1049 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1050 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1052 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1055 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1056 tmp1(ixo^s)=wmean(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim)/wmean(ixo^s,
rho_)
1058 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1059 if(
present(cmin))
then
1060 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1061 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1063 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1069 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1070 if(
present(cmin))
then
1071 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1072 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))-csoundl(ixo^s)
1073 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1074 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1076 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1077 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1080 end subroutine ffhd_get_cbounds
1082 subroutine ffhd_get_pthermal_iso(w,x,ixI^L,ixO^L,pth)
1085 integer,
intent(in) :: ixi^
l, ixo^
l
1086 double precision,
intent(in) :: w(ixi^s,nw)
1087 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1088 double precision,
intent(out):: pth(ixi^s)
1092 end subroutine ffhd_get_pthermal_iso
1094 subroutine ffhd_get_pthermal_origin(w,x,ixI^L,ixO^L,pth)
1097 integer,
intent(in) :: ixi^
l, ixo^
l
1098 double precision,
intent(in) :: w(ixi^s,nw)
1099 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1100 double precision,
intent(out):: pth(ixi^s)
1105 {
do ix^db= ixo^lim^db\}
1110 elseif(check_small_values)
then
1111 {
do ix^db= ixo^lim^db\}
1112 if(pth(ix^d)<small_pressure)
then
1113 write(*,*)
"Error: small value of gas pressure",pth(ix^d),&
1114 " encountered when call ffhd_get_pthermal"
1115 write(*,*)
"Iteration: ", it,
" Time: ", global_time
1116 write(*,*)
"Location: ", x(ix^d,:)
1117 write(*,*)
"Cell number: ", ix^d
1119 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^d,iw)
1122 if(trace_small_values)
write(*,*) dsqrt(pth(ix^d)-bigdouble)
1123 write(*,*)
"Saving status at the previous time step"
1128 end subroutine ffhd_get_pthermal_origin
1130 subroutine ffhd_get_temperature_from_te(w, x, ixI^L, ixO^L, res)
1132 integer,
intent(in) :: ixi^
l, ixo^
l
1133 double precision,
intent(in) :: w(ixi^s, 1:nw)
1134 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1135 double precision,
intent(out):: res(ixi^s)
1137 res(ixo^s) = w(ixo^s,
te_)
1138 end subroutine ffhd_get_temperature_from_te
1140 subroutine ffhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1142 integer,
intent(in) :: ixi^
l, ixo^
l
1143 double precision,
intent(in) :: w(ixi^s, 1:nw)
1144 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1145 double precision,
intent(out):: res(ixi^s)
1146 double precision :: r(ixi^s)
1148 call ffhd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1149 res(ixo^s) = gamma_1 * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1150 end subroutine ffhd_get_temperature_from_eint
1152 subroutine ffhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1154 integer,
intent(in) :: ixi^
l, ixo^
l
1155 double precision,
intent(in) :: w(ixi^s, 1:nw)
1156 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1157 double precision,
intent(out):: res(ixi^s)
1159 double precision :: r(ixi^s)
1161 call ffhd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1163 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1164 end subroutine ffhd_get_temperature_from_etot
1168 integer,
intent(in) :: ixi^
l, ixo^
l
1169 double precision,
intent(in) :: w(ixi^s,nw)
1170 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1171 double precision,
intent(out) :: csound2(ixi^s)
1172 double precision :: rho(ixi^s)
1177 csound2(ixo^s)=
ffhd_gamma*csound2(ixo^s)/rho(ixo^s)
1183 subroutine ffhd_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
1186 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1188 double precision,
intent(in) :: wc(ixi^s,nw)
1190 double precision,
intent(in) :: w(ixi^s,nw)
1191 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1192 double precision,
intent(out) :: f(ixi^s,nwflux)
1193 double precision :: ptotal(ixo^s)
1198 ptotal(ixo^s)=w(ixo^s,
p_)
1204 f(ixo^s,
mom(1))=(wc(ixo^s,
mom(1))*w(ixo^s,
mom(1))+ptotal(ixo^s))*
block%B0(ixo^s,idim,idim)
1208 f(ixo^s,
e_)=w(ixo^s,
mom(1))*(wc(ixo^s,
e_)+ptotal(ixo^s))*
block%B0(ixo^s,idim,idim)
1210 f(ixo^s,
e_)=f(ixo^s,
e_)+w(ixo^s,
q_)*
block%B0(ixo^s,idim,idim)
1214 end subroutine ffhd_get_flux
1216 subroutine ffhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1220 integer,
intent(in) :: ixi^
l, ixo^
l
1221 double precision,
intent(in) :: qdt,dtfactor
1222 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:
ndim)
1223 double precision,
intent(inout) :: w(ixi^s,1:nw)
1224 logical,
intent(in) :: qsourcesplit
1225 logical,
intent(inout) :: active
1227 if (.not. qsourcesplit)
then
1229 call add_punitb(qdt,ixi^
l,ixo^
l,wct,w,x,wctprim)
1231 call add_hyperbolic_tc_source(qdt,ixi^
l,ixo^
l,wct,w,x,wctprim)
1237 w,x,qsourcesplit,active,
rc_fl)
1247 if(.not.qsourcesplit)
then
1249 call ffhd_update_temperature(ixi^
l,ixo^
l,wct,w,x)
1252 end subroutine ffhd_add_source
1254 subroutine add_punitb(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1257 integer,
intent(in) :: ixi^
l,ixo^
l
1258 double precision,
intent(in) :: qdt
1259 double precision,
intent(in) :: wct(ixi^s,1:nw),x(ixi^s,1:
ndim)
1260 double precision,
intent(in) :: wctprim(ixi^s,1:nw)
1261 double precision,
intent(inout) :: w(ixi^s,1:nw)
1263 integer :: idims,hxo^
l
1264 double precision :: divb(ixi^s)
1265 double precision :: rhovpar(ixi^s),gradrhov(ixi^s)
1270 hxo^
l=ixo^
l-
kr(idims,^
d);
1271 divb(ixo^s)=divb(ixo^s)+(
block%B0(ixo^s,idims,idims)-
block%B0(hxo^s,idims,idims))/
dxlevel(idims)
1276 w(ixo^s,
mom(1))=w(ixo^s,
mom(1))+qdt*wctprim(ixo^s,
p_)*divb(ixo^s)
1277 end subroutine add_punitb
1281 integer,
intent(in) :: ixi^
l, ixo^
l
1282 double precision,
intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:
ndim)
1283 double precision,
intent(out) :: rho(ixi^s)
1285 rho(ixo^s) = w(ixo^s,
rho_)
1288 subroutine ffhd_handle_small_ei(w, x, ixI^L, ixO^L, ie, subname)
1291 integer,
intent(in) :: ixi^
l,ixo^
l, ie
1292 double precision,
intent(inout) :: w(ixi^s,1:nw)
1293 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1294 character(len=*),
intent(in) :: subname
1296 logical :: flag(ixi^s,1:nw)
1297 double precision :: rho(ixi^s)
1300 where(w(ixo^s,ie)<
small_e) flag(ixo^s,ie)=.true.
1301 if(any(flag(ixo^s,ie)))
then
1304 where(flag(ixo^s,ie)) w(ixo^s,ie)=
small_e
1308 w(ixo^s,
e_)=w(ixo^s,
e_)*gamma_1
1310 w(ixo^s,
mom(1)) = w(ixo^s,
mom(1))/rho(ixo^s)
1314 end subroutine ffhd_handle_small_ei
1316 subroutine ffhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1319 integer,
intent(in) :: ixi^
l, ixo^
l
1320 double precision,
intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:
ndim)
1321 double precision,
intent(inout) :: w(ixi^s,1:nw)
1322 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1328 end subroutine ffhd_update_temperature
1330 subroutine ffhd_get_dt(wprim,ixI^L,ixO^L,dtnew,dx^D,x)
1334 integer,
intent(in) :: ixi^
l, ixo^
l
1335 double precision,
intent(inout) :: dtnew
1336 double precision,
intent(in) ::
dx^
d
1337 double precision,
intent(in) :: wprim(ixi^s,1:nw)
1338 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1345 end subroutine ffhd_get_dt
1347 subroutine ffhd_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
1349 integer,
intent(in) :: ixi^
l, ixo^
l
1350 double precision,
intent(in) :: qdt, dtfactor,x(ixi^s,1:
ndim)
1351 double precision,
intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw),w(ixi^s,1:nw)
1354 end subroutine ffhd_add_source_geom
1356 function ffhd_kin_en_origin(w, ixI^L, ixO^L, inv_rho)
result(ke)
1358 integer,
intent(in) :: ixi^
l, ixo^
l
1359 double precision,
intent(in) :: w(ixi^s, nw)
1360 double precision :: ke(ixo^s)
1361 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1363 if(
present(inv_rho))
then
1364 ke(ixo^s)=0.5d0*w(ixo^s,
mom(1))**2*inv_rho(ixo^s)
1366 ke(ixo^s)=0.5d0*w(ixo^s,
mom(1))**2/w(ixo^s,
rho_)
1368 end function ffhd_kin_en_origin
1370 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1373 integer,
intent(in) :: ixi^
l, ixo^
l
1374 double precision,
intent(in) :: w(ixi^s,1:nw)
1375 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1376 double precision,
intent(out):: rfactor(ixi^s)
1377 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1380 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)
1381 end subroutine rfactor_from_temperature_ionization
1383 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1385 integer,
intent(in) :: ixi^
l, ixo^
l
1386 double precision,
intent(in) :: w(ixi^s,1:nw)
1387 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1388 double precision,
intent(out):: rfactor(ixi^s)
1391 end subroutine rfactor_from_constant_ionization
1393 subroutine add_hyperbolic_tc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1396 integer,
intent(in) :: ixi^
l,ixo^
l
1397 double precision,
intent(in) :: qdt
1398 double precision,
dimension(ixI^S,1:ndim),
intent(in) :: x
1399 double precision,
dimension(ixI^S,1:nw),
intent(in) :: wct,wctprim
1400 double precision,
dimension(ixI^S,1:nw),
intent(inout) :: w
1402 double precision,
dimension(ixI^S) :: te
1403 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1406 te(ixi^s)=wctprim(ixi^s,
p_)/wct(ixi^s,
rho_)
1410 do ix2=ixomin2,ixomax2
1411 do ix1=ixomin1,ixomax1
1418 sigma_t7=sigma_t5*te(ix^
d)
1422 sigma_t7=sigma_t5*te(ix^
d)
1424 sigmat5_bgradt=sigma_t5*(&
1425 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)&
1426 +
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))
1429 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^
d,
rho_)*(wctprim(ix^
d,
p_)/wct(ix^
d,
rho_))**1.5d0))
1431 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^
d,
q_))/tau
1433 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(sigmat5_bgradt+wct(ix^
d,
q_))/&
1440 do ix3=ixomin3,ixomax3
1441 do ix2=ixomin2,ixomax2
1442 do ix1=ixomin1,ixomax1
1449 sigma_t7=sigma_t5*te(ix^
d)
1453 sigma_t7=sigma_t5*te(ix^
d)
1455 sigmat5_bgradt=sigma_t5*(&
1456 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)&
1457 +
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)&
1458 +
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))
1461 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^
d,
rho_)*(wctprim(ix^
d,
p_)/wct(ix^
d,
rho_))**1.5d0))
1463 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^
d,
q_))/tau
1465 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(sigmat5_bgradt+wct(ix^
d,
q_))/&
1472 end subroutine add_hyperbolic_tc_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.
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 hyperbolic_tc_kappa
The thermal conductivity kappa in hyperbolic thermal conduction.
double precision, public ffhd_adiab
The adiabatic constant.
procedure(sub_get_pthermal), pointer, public ffhd_get_temperature
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
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.
logical, public, protected ffhd_hyperbolic_tc_use_perp
Whether the perpendicular hyperbolic-TC channel is enabled.
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_hyperbolic_tc
Whether hyperbolic type thermal conduction is used.
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
logical, public, protected ffhd_hyperbolic_tc_sat
Whether saturation is considered for hyperbolic TC.
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.
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.
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
integer nghostcells
Number of ghost cells surrounding a grid.
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(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
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()
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