57 integer,
public,
protected ::
rho_
60 integer,
allocatable,
public,
protected ::
mom(:)
63 integer,
public,
protected ::
e_
66 integer,
public,
protected ::
p_
69 integer,
public,
protected ::
te_
74 integer,
public,
protected ::
q_
83 double precision,
protected :: small_e
92 double precision,
public,
protected ::
h_ion_fr=1d0
95 double precision,
public,
protected ::
he_ion_fr=1d0
102 double precision,
public,
protected ::
rr=1d0
109 double precision :: gamma_1, inv_gamma_1
114 function fun_kin_en(w, ixI^L, ixO^L, inv_rho)
result(ke)
116 integer,
intent(in) :: ixi^
l, ixo^
l
117 double precision,
intent(in) :: w(ixi^s, nw)
118 double precision :: ke(ixo^s)
119 double precision,
intent(in),
optional :: inv_rho(ixo^s)
120 end function fun_kin_en
148 subroutine ffhd_read_params(files)
151 character(len=*),
intent(in) :: files(:)
159 do n = 1,
size(files)
160 open(
unitpar, file=trim(files(n)), status=
"old")
161 read(
unitpar, ffhd_list,
end=111)
164 end subroutine ffhd_read_params
167 subroutine ffhd_write_info(fh)
169 integer,
intent(in) :: fh
170 integer,
parameter :: n_par = 1
171 double precision :: values(n_par)
172 character(len=name_len) :: names(n_par)
173 integer,
dimension(MPI_STATUS_SIZE) :: st
176 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
180 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
181 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
182 end subroutine ffhd_write_info
201 if(
mype==0)
write(*,*)
'WARNING: set ffhd_thermal_conduction=F when ffhd_energy=F'
205 if(
mype==0)
write(*,*)
'WARNING: set ffhd_hyperbolic_thermal_conduction=F when ffhd_energy=F'
209 if(
mype==0)
write(*,*)
'WARNING: set ffhd_radiative_cooling=F when ffhd_energy=F'
213 if(
mype==0)
write(*,*)
'WARNING: set ffhd_trac=F when ffhd_energy=F'
217 if(
mype==0)
write(*,*)
'WARNING: set ffhd_partial_ionization=F when ffhd_energy=F'
223 if(
mype==0)
write(*,*)
'WARNING: set ffhd_partial_ionization=F when eq_state_units=F'
229 if(
mype==0)
write(*,*)
'WARNING: turn off parabolic TC when using hyperbolic TC'
232 physics_type =
"ffhd"
234 phys_internal_e=.false.
246 if(
mype==0)
write(*,*)
'WARNING: reset ffhd_trac_type=1 for 1D simulation'
251 if(
mype==0)
write(*,*)
'WARNING: set ffhd_trac_mask==bigdouble for global TRAC method'
255 allocate(start_indices(number_species),stop_indices(number_species))
261 mom(:) = var_set_momentum(1)
265 e_ = var_set_energy()
283 stop_indices(1)=nwflux
287 te_ = var_set_auxvar(
'Te',
'Te')
308 if(.not.
allocated(flux_type))
then
309 allocate(flux_type(
ndir, nwflux))
310 flux_type = flux_default
311 else if(any(shape(flux_type) /= [
ndir, nwflux]))
then
312 call mpistop(
"phys_check error: flux_type has wrong shape")
315 phys_get_dt => ffhd_get_dt
316 phys_get_cmax => ffhd_get_cmax_origin
317 phys_get_a2max => ffhd_get_a2max
318 phys_get_tcutoff => ffhd_get_tcutoff
319 phys_get_cbounds => ffhd_get_cbounds
320 phys_to_primitive => ffhd_to_primitive_origin
322 phys_to_conserved => ffhd_to_conserved_origin
324 phys_get_flux => ffhd_get_flux
325 phys_get_v => ffhd_get_v_origin
329 phys_add_source_geom => ffhd_add_source_geom
330 phys_add_source => ffhd_add_source
331 phys_check_params => ffhd_check_params
332 phys_write_info => ffhd_write_info
333 phys_handle_small_values => ffhd_handle_small_values_origin
334 ffhd_handle_small_values => ffhd_handle_small_values_origin
335 phys_check_w => ffhd_check_w_origin
338 phys_get_pthermal => ffhd_get_pthermal_iso
341 phys_get_pthermal => ffhd_get_pthermal_origin
347 ffhd_get_rfactor=>rfactor_from_temperature_ionization
348 phys_update_temperature => ffhd_update_temperature
352 ffhd_get_rfactor=>rfactor_from_constant_ionization
362 call ffhd_physical_units()
368 call mpistop(
"thermal conduction needs ffhd_energy=T")
371 call mpistop(
"hyperbolic thermal conduction needs ffhd_energy=T")
374 call mpistop(
"radiative cooling needs ffhd_energy=T")
384 call add_sts_method(ffhd_get_tc_dt_ffhd,ffhd_sts_set_source_tc_ffhd,
e_,1,
e_,1,.false.)
385 tc_fl%get_temperature_from_conserved => ffhd_get_temperature_from_etot
386 tc_fl%get_temperature_from_eint => ffhd_get_temperature_from_eint
401 rc_fl%get_var_Rfactor => ffhd_get_rfactor
404 rc_fl%has_equi = .false.
409 te_fl_ffhd%get_var_Rfactor => ffhd_get_rfactor
411 phys_te_images => ffhd_te_images
426 subroutine ffhd_te_images
431 case(
'EIvtiCCmpi',
'EIvtuCCmpi')
433 case(
'ESvtiCCmpi',
'ESvtuCCmpi')
435 case(
'SIvtiCCmpi',
'SIvtuCCmpi')
437 case(
'WIvtiCCmpi',
'WIvtuCCmpi')
440 call mpistop(
"Error in synthesize emission: Unknown convert_type")
442 end subroutine ffhd_te_images
445 subroutine ffhd_sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
449 integer,
intent(in) :: ixi^
l, ixo^
l, igrid, nflux
450 double precision,
intent(in) :: x(ixi^s,1:
ndim)
451 double precision,
intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
452 double precision,
intent(in) :: my_dt
453 logical,
intent(in) :: fix_conserve_at_step
455 end subroutine ffhd_sts_set_source_tc_ffhd
457 function ffhd_get_tc_dt_ffhd(w,ixI^L,ixO^L,dx^D,x)
result(dtnew)
464 integer,
intent(in) :: ixi^
l, ixo^
l
465 double precision,
intent(in) ::
dx^
d, x(ixi^s,1:
ndim)
466 double precision,
intent(in) :: w(ixi^s,1:nw)
467 double precision :: dtnew
470 end function ffhd_get_tc_dt_ffhd
472 subroutine ffhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
475 integer,
intent(in) :: ixi^
l,ixo^
l
476 double precision,
intent(inout) :: w(ixi^s,1:nw)
477 double precision,
intent(in) :: x(ixi^s,1:
ndim)
478 integer,
intent(in) :: step
479 character(len=140) :: error_msg
481 write(error_msg,
"(a,i3)")
"Thermal conduction step ", step
482 call ffhd_handle_small_ei(w,x,ixi^
l,ixo^
l,
e_,error_msg)
483 end subroutine ffhd_tc_handle_small_e
485 subroutine tc_params_read_ffhd(fl)
487 type(tc_fluid),
intent(inout) :: fl
490 logical :: tc_saturate=.false.
491 double precision :: tc_k_para=0d0
492 character(len=std_len) :: tc_slope_limiter=
"MC"
494 namelist /tc_list/ tc_saturate, tc_slope_limiter, tc_k_para
498 read(
unitpar, tc_list,
end=111)
502 fl%tc_saturate = tc_saturate
503 fl%tc_k_para = tc_k_para
504 select case(tc_slope_limiter)
506 fl%tc_slope_limiter = 0
509 fl%tc_slope_limiter = 1
512 fl%tc_slope_limiter = 2
515 fl%tc_slope_limiter = 3
518 fl%tc_slope_limiter = 4
520 call mpistop(
"Unknown tc_slope_limiter, choose MC, minmod")
522 end subroutine tc_params_read_ffhd
524 subroutine rc_params_read(fl)
527 type(rc_fluid),
intent(inout) :: fl
529 integer :: ncool = 4000
530 double precision :: cfrac=0.1d0
533 character(len=std_len) :: coolcurve=
'JCcorona'
536 character(len=std_len) :: coolmethod=
'exact'
539 logical :: tfix=.false.
545 logical :: rc_split=.false.
546 logical :: rad_cut=.false.
547 double precision :: rad_cut_hgt=0.5d0
548 double precision :: rad_cut_dey=0.15d0
550 namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split, rad_cut, rad_cut_hgt, rad_cut_dey
554 read(
unitpar, rc_list,
end=111)
559 fl%coolcurve=coolcurve
560 fl%coolmethod=coolmethod
566 fl%rad_cut_hgt=rad_cut_hgt
567 fl%rad_cut_dey=rad_cut_dey
568 end subroutine rc_params_read
570 subroutine ffhd_check_params
577 if (
ffhd_gamma <= 0.0d0)
call mpistop (
"Error: ffhd_gamma <= 0")
578 if (
ffhd_adiab < 0.0d0)
call mpistop (
"Error: ffhd_adiab < 0")
582 call mpistop (
"Error: ffhd_gamma <= 0 or ffhd_gamma == 1")
583 inv_gamma_1=1.d0/gamma_1
588 call mpistop(
"usr_set_equi_vars has to be implemented in the user file")
590 end subroutine ffhd_check_params
592 subroutine ffhd_physical_units()
594 double precision :: mp,kb
595 double precision :: a,b
682 end subroutine ffhd_physical_units
684 subroutine ffhd_check_w_origin(primitive,ixI^L,ixO^L,w,flag)
686 logical,
intent(in) :: primitive
687 integer,
intent(in) :: ixi^
l, ixo^
l
688 double precision,
intent(in) :: w(ixi^s,nw)
689 double precision :: tmp(ixi^s)
690 logical,
intent(inout) :: flag(ixi^s,1:nw)
700 where(tmp(ixo^s) < small_e) flag(ixo^s,
e_) = .true.
703 end subroutine ffhd_check_w_origin
705 subroutine ffhd_to_conserved_origin(ixI^L,ixO^L,w,x)
707 integer,
intent(in) :: ixi^
l, ixo^
l
708 double precision,
intent(inout) :: w(ixi^s, nw)
709 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
712 w(ixo^s,
e_)=w(ixo^s,
p_)*inv_gamma_1+half*w(ixo^s,
mom(1))**2*w(ixo^s,
rho_)
715 end subroutine ffhd_to_conserved_origin
717 subroutine ffhd_to_primitive_origin(ixI^L,ixO^L,w,x)
719 integer,
intent(in) :: ixi^
l, ixo^
l
720 double precision,
intent(inout) :: w(ixi^s, nw)
721 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
724 call ffhd_handle_small_values(.false., w, x, ixi^
l, ixo^
l,
'ffhd_to_primitive_origin')
727 w(ixo^s,
mom(1)) = w(ixo^s,
mom(1))/w(ixo^s,
rho_)
729 w(ixo^s,
p_)=gamma_1*(w(ixo^s,
e_)-half*w(ixo^s,
rho_)*w(ixo^s,
mom(1))**2)
731 end subroutine ffhd_to_primitive_origin
735 integer,
intent(in) :: ixi^
l, ixo^
l
736 double precision,
intent(inout) :: w(ixi^s, nw)
737 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
744 integer,
intent(in) :: ixi^
l, ixo^
l
745 double precision,
intent(inout) :: w(ixi^s, nw)
746 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
750 call ffhd_handle_small_ei(w,x,ixi^
l,ixi^
l,
e_,
'ffhd_e_to_ei')
754 subroutine ffhd_handle_small_values_origin(primitive, w, x, ixI^L, ixO^L, subname)
757 logical,
intent(in) :: primitive
758 integer,
intent(in) :: ixi^
l,ixo^
l
759 double precision,
intent(inout) :: w(ixi^s,1:nw)
760 double precision,
intent(in) :: x(ixi^s,1:
ndim)
761 character(len=*),
intent(in) :: subname
763 logical :: flag(ixi^s,1:nw)
764 double precision :: tmp2(ixi^s)
766 call phys_check_w(primitive, ixi^
l, ixi^
l, w, flag)
773 where(flag(ixo^s,
rho_)) w(ixo^s,
mom(1)) = 0.0d0
779 where(flag(ixo^s,
e_))
796 if(.not.primitive)
then
805 end subroutine ffhd_handle_small_values_origin
807 subroutine ffhd_get_v_origin(w,x,ixI^L,ixO^L,v)
809 integer,
intent(in) :: ixi^
l, ixo^
l
810 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
811 double precision,
intent(out) :: v(ixi^s,
ndir)
812 double precision :: rho(ixi^s)
816 rho(ixo^s)=1.d0/rho(ixo^s)
818 v(ixo^s,
ndir) = w(ixo^s,
mom(1))*
block%B0(ixo^s,idir,0)*rho(ixo^s)
820 end subroutine ffhd_get_v_origin
824 integer,
intent(in) :: ixi^
l, ixo^
l, idim
825 double precision,
intent(in) :: w(ixi^s,nw), x(ixi^s,1:
ndim)
826 double precision,
intent(out) :: v(ixi^s)
827 double precision :: rho(ixi^s)
830 v(ixo^s) = (w(ixo^s,
mom(1))*
block%B0(ixo^s,idim,0)) / rho(ixo^s)
833 subroutine ffhd_get_cmax_origin(wprim,x,ixI^L,ixO^L,idim,cmax)
835 integer,
intent(in) :: ixi^
l, ixo^
l, idim
837 double precision,
intent(in) :: wprim(ixi^s, nw), x(ixi^s,1:
ndim)
838 double precision,
intent(inout) :: cmax(ixi^s)
845 cmax(ixo^s)=dabs(wprim(ixo^s,
mom(1))*
block%B0(ixo^s,idim,0))+cmax(ixo^s)
847 end subroutine ffhd_get_cmax_origin
849 subroutine ffhd_get_a2max(w,x,ixI^L,ixO^L,a2max)
852 integer,
intent(in) :: ixi^
l, ixo^
l
853 double precision,
intent(in) :: w(ixi^s, nw), x(ixi^s,1:
ndim)
854 double precision,
intent(inout) :: a2max(
ndim)
855 double precision :: a2(ixi^s,
ndim,nw)
856 integer :: gxo^
l,hxo^
l,jxo^
l,kxo^
l,i,j
859 call mpistop(
"subroutine get_a2max in mod_ffhd_phys adopts cartesian setting")
864 hxo^
l=ixo^
l-
kr(i,^
d);
865 gxo^
l=hxo^
l-
kr(i,^
d);
866 jxo^
l=ixo^
l+
kr(i,^
d);
867 kxo^
l=jxo^
l+
kr(i,^
d);
868 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
869 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
870 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/
dxlevel(i)**2
872 end subroutine ffhd_get_a2max
874 subroutine ffhd_get_tcutoff(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
877 integer,
intent(in) :: ixi^
l,ixo^
l
878 double precision,
intent(in) :: x(ixi^s,1:
ndim)
879 double precision,
intent(inout) :: w(ixi^s,1:nw)
880 double precision,
intent(out) :: tco_local,tmax_local
881 double precision,
parameter :: trac_delta=0.25d0
882 double precision :: te(ixi^s),lts(ixi^s)
883 double precision,
dimension(ixI^S,1:ndim) :: gradt
884 double precision :: bdir(
ndim)
885 double precision :: ltrc,ltrp,altr
886 integer :: idims,ix^
d,jxo^
l,hxo^
l,ixa^
d,ixb^
d
887 integer :: jxp^
l,hxp^
l,ixp^
l,ixq^
l
888 logical :: lrlt(ixi^s)
892 tmax_local=maxval(te(ixo^s))
902 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
904 where(lts(ixo^s) > trac_delta)
907 if(any(lrlt(ixo^s)))
then
908 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
919 lts(ixp^s)=0.5d0*dabs(te(jxp^s)-te(hxp^s))/te(ixp^s)
920 lts(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
921 lts(ixo^s)=0.25d0*(lts(jxo^s)+two*lts(ixo^s)+lts(hxo^s))
922 block%wextra(ixo^s,
tcoff_)=te(ixo^s)*lts(ixo^s)**0.4d0
924 call mpistop(
"ffhd_trac_type not allowed for 1D simulation")
939 call gradient(te,ixi^
l,ixo^
l,idims,gradt(ixi^s,idims))
945 ixb^
d=(ixomin^
d+ixomax^
d-1)/2+ixa^
d;
948 if(sum(bdir(:)**2) .gt. zero)
then
949 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
951 block%special_values(3:ndim+2)=bdir(1:ndim)
953 {
do ix^db=ixomin^db,ixomax^db\}
956 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2))*&
957 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
961 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2),block%ds(ix^d,3))*&
962 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
964 if(lts(ix^d)>trac_delta)
then
965 block%special_values(1)=max(block%special_values(1),te(ix^d))
968 block%special_values(2)=tmax_local
986 call gradient(te,ixi^l,ixq^l,idims,gradt(ixi^s,idims))
987 call gradientf(te,x,ixi^l,hxp^l,idims,gradt(ixi^s,idims),nghostcells,.true.)
988 call gradientf(te,x,ixi^l,jxp^l,idims,gradt(ixi^s,idims),nghostcells,.false.)
990 {
do ix^db=ixpmin^db,ixpmax^db\}
992 lts(ix^d)=abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
994 lts(ix^d)=min(^d&block%ds({ix^d},^d))*lts(ix^d)
995 lts(ix^d)=max(one,(exp(lts(ix^d))/ltrc)**ltrp)
1000 {
do ix^db=ixpmin^db,ixpmax^db\}
1002 altr=0.25d0*((lts(ix1-1,ix2)+two*lts(ix^d)+lts(ix1+1,ix2))*block%B0(ix^d,1,0)**2+&
1003 (lts(ix1,ix2-1)+two*lts(ix^d)+lts(ix1,ix2+1))*block%B0(ix^d,2,0)**2)
1004 block%wextra(ix^d,
tcoff_)=te(ix^d)*altr**0.4d0
1007 altr=0.25d0*((lts(ix1-1,ix2,ix3)+two*lts(ix^d)+lts(ix1+1,ix2,ix3))*block%B0(ix^d,1,0)**2+&
1008 (lts(ix1,ix2-1,ix3)+two*lts(ix^d)+lts(ix1,ix2+1,ix3))*block%B0(ix^d,2,0)**2+&
1009 (lts(ix1,ix2,ix3-1)+two*lts(ix^d)+lts(ix1,ix2,ix3+1))*block%B0(ix^d,3,0)**2)
1010 block%wextra(ix^d,
tcoff_)=te(ix^d)*altr**0.4d0
1016 call mpistop(
"unknown ffhd_trac_type")
1019 end subroutine ffhd_get_tcutoff
1021 subroutine ffhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
1023 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1024 double precision,
intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1025 double precision,
intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1026 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1027 double precision,
intent(inout) :: cmax(ixi^s,1:number_species)
1028 double precision,
intent(inout),
optional :: cmin(ixi^s,1:number_species)
1029 double precision,
intent(in) :: hspeed(ixi^s,1:number_species)
1030 double precision :: wmean(ixi^s,nw)
1031 double precision,
dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1037 tmp1(ixo^s)=dsqrt(wlp(ixo^s,
rho_))
1038 tmp2(ixo^s)=dsqrt(wrp(ixo^s,
rho_))
1039 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1040 umean(ixo^s)=(wlp(ixo^s,
mom(1))*tmp1(ixo^s)&
1041 +wrp(ixo^s,
mom(1))*tmp2(ixo^s))*tmp3(ixo^s)
1042 umean(ixo^s)=umean(ixo^s)*
block%B0(ixo^s,idim,idim)
1045 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1046 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1047 ((wrp(ixo^s,
mom(1))-wlp(ixo^s,
mom(1)))*
block%B0(ixo^s,idim,idim))**2
1048 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1049 if(
present(cmin))
then
1050 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1051 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1053 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1056 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1057 tmp1(ixo^s)=wmean(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim)/wmean(ixo^s,
rho_)
1059 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1060 if(
present(cmin))
then
1061 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1062 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1064 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1070 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1071 if(
present(cmin))
then
1072 cmin(ixo^s,1)=min(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1073 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))-csoundl(ixo^s)
1074 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1075 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1077 cmax(ixo^s,1)=max(wlp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim),&
1078 wrp(ixo^s,
mom(1))*
block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1081 end subroutine ffhd_get_cbounds
1083 subroutine ffhd_get_pthermal_iso(w,x,ixI^L,ixO^L,pth)
1086 integer,
intent(in) :: ixi^
l, ixo^
l
1087 double precision,
intent(in) :: w(ixi^s,nw)
1088 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1089 double precision,
intent(out):: pth(ixi^s)
1093 end subroutine ffhd_get_pthermal_iso
1095 subroutine ffhd_get_pthermal_origin(w,x,ixI^L,ixO^L,pth)
1098 integer,
intent(in) :: ixi^
l, ixo^
l
1099 double precision,
intent(in) :: w(ixi^s,nw)
1100 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1101 double precision,
intent(out):: pth(ixi^s)
1106 {
do ix^db= ixo^lim^db\}
1111 elseif(check_small_values)
then
1112 {
do ix^db= ixo^lim^db\}
1113 if(pth(ix^d)<small_pressure)
then
1114 write(*,*)
"Error: small value of gas pressure",pth(ix^d),&
1115 " encountered when call ffhd_get_pthermal"
1116 write(*,*)
"Iteration: ", it,
" Time: ", global_time
1117 write(*,*)
"Location: ", x(ix^d,:)
1118 write(*,*)
"Cell number: ", ix^d
1120 write(*,*) trim(cons_wnames(iw)),
": ",w(ix^d,iw)
1123 if(trace_small_values)
write(*,*) dsqrt(pth(ix^d)-bigdouble)
1124 write(*,*)
"Saving status at the previous time step"
1129 end subroutine ffhd_get_pthermal_origin
1131 subroutine ffhd_get_temperature_from_te(w, x, ixI^L, ixO^L, res)
1133 integer,
intent(in) :: ixi^
l, ixo^
l
1134 double precision,
intent(in) :: w(ixi^s, 1:nw)
1135 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1136 double precision,
intent(out):: res(ixi^s)
1138 res(ixo^s) = w(ixo^s,
te_)
1139 end subroutine ffhd_get_temperature_from_te
1141 subroutine ffhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1143 integer,
intent(in) :: ixi^
l, ixo^
l
1144 double precision,
intent(in) :: w(ixi^s, 1:nw)
1145 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1146 double precision,
intent(out):: res(ixi^s)
1147 double precision :: r(ixi^s)
1149 call ffhd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1150 res(ixo^s) = gamma_1 * w(ixo^s,
e_)/(w(ixo^s,
rho_)*r(ixo^s))
1151 end subroutine ffhd_get_temperature_from_eint
1153 subroutine ffhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1155 integer,
intent(in) :: ixi^
l, ixo^
l
1156 double precision,
intent(in) :: w(ixi^s, 1:nw)
1157 double precision,
intent(in) :: x(ixi^s, 1:
ndim)
1158 double precision,
intent(out):: res(ixi^s)
1160 double precision :: r(ixi^s)
1162 call ffhd_get_rfactor(w,x,ixi^
l,ixo^
l,r)
1164 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,
rho_))
1165 end subroutine ffhd_get_temperature_from_etot
1169 integer,
intent(in) :: ixi^
l, ixo^
l
1170 double precision,
intent(in) :: w(ixi^s,nw)
1171 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1172 double precision,
intent(out) :: csound2(ixi^s)
1173 double precision :: rho(ixi^s)
1178 csound2(ixo^s)=
ffhd_gamma*csound2(ixo^s)/rho(ixo^s)
1184 subroutine ffhd_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
1187 integer,
intent(in) :: ixi^
l, ixo^
l, idim
1189 double precision,
intent(in) :: wc(ixi^s,nw)
1191 double precision,
intent(in) :: w(ixi^s,nw)
1192 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1193 double precision,
intent(out) :: f(ixi^s,nwflux)
1194 double precision :: ptotal(ixo^s)
1199 ptotal(ixo^s)=w(ixo^s,
p_)
1205 f(ixo^s,
mom(1))=(wc(ixo^s,
mom(1))*w(ixo^s,
mom(1))+ptotal(ixo^s))*
block%B0(ixo^s,idim,idim)
1209 f(ixo^s,
e_)=w(ixo^s,
mom(1))*(wc(ixo^s,
e_)+ptotal(ixo^s))*
block%B0(ixo^s,idim,idim)
1211 f(ixo^s,
e_)=f(ixo^s,
e_)+w(ixo^s,
q_)*
block%B0(ixo^s,idim,idim)
1215 end subroutine ffhd_get_flux
1217 subroutine ffhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1222 integer,
intent(in) :: ixi^
l, ixo^
l
1223 double precision,
intent(in) :: qdt,dtfactor
1224 double precision,
intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:
ndim)
1225 double precision,
intent(inout) :: w(ixi^s,1:nw)
1226 logical,
intent(in) :: qsourcesplit
1227 logical,
intent(inout) :: active
1229 if (.not. qsourcesplit)
then
1231 call add_punitb(qdt,ixi^
l,ixo^
l,wct,w,x,wctprim)
1233 call add_hypertc_source(qdt,ixi^
l,ixo^
l,wct,w,x,wctprim)
1239 w,x,qsourcesplit,active,
rc_fl)
1254 if(.not.qsourcesplit)
then
1256 call ffhd_update_temperature(ixi^
l,ixo^
l,wct,w,x)
1259 end subroutine ffhd_add_source
1261 subroutine add_punitb(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1264 integer,
intent(in) :: ixi^
l,ixo^
l
1265 double precision,
intent(in) :: qdt
1266 double precision,
intent(in) :: wct(ixi^s,1:nw),x(ixi^s,1:
ndim)
1267 double precision,
intent(in) :: wctprim(ixi^s,1:nw)
1268 double precision,
intent(inout) :: w(ixi^s,1:nw)
1270 integer :: idims,hxo^
l
1271 double precision :: divb(ixi^s)
1272 double precision :: rhovpar(ixi^s),gradrhov(ixi^s)
1277 hxo^
l=ixo^
l-
kr(idims,^
d);
1278 divb(ixo^s)=divb(ixo^s)+(
block%B0(ixo^s,idims,idims)-
block%B0(hxo^s,idims,idims))/
dxlevel(idims)
1283 w(ixo^s,
mom(1))=w(ixo^s,
mom(1))+qdt*wctprim(ixo^s,
p_)*divb(ixo^s)
1284 end subroutine add_punitb
1288 integer,
intent(in) :: ixi^
l, ixo^
l
1289 double precision,
intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:
ndim)
1290 double precision,
intent(out) :: rho(ixi^s)
1292 rho(ixo^s) = w(ixo^s,
rho_)
1295 subroutine ffhd_handle_small_ei(w, x, ixI^L, ixO^L, ie, subname)
1298 integer,
intent(in) :: ixi^
l,ixo^
l, ie
1299 double precision,
intent(inout) :: w(ixi^s,1:nw)
1300 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1301 character(len=*),
intent(in) :: subname
1303 logical :: flag(ixi^s,1:nw)
1304 double precision :: rho(ixi^s)
1307 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
1308 if(any(flag(ixo^s,ie)))
then
1311 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
1315 w(ixo^s,
e_)=w(ixo^s,
e_)*gamma_1
1317 w(ixo^s,
mom(1)) = w(ixo^s,
mom(1))/rho(ixo^s)
1321 end subroutine ffhd_handle_small_ei
1323 subroutine ffhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1326 integer,
intent(in) :: ixi^
l, ixo^
l
1327 double precision,
intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:
ndim)
1328 double precision,
intent(inout) :: w(ixi^s,1:nw)
1329 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1335 end subroutine ffhd_update_temperature
1337 subroutine ffhd_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
1343 integer,
intent(in) :: ixi^
l, ixo^
l
1344 double precision,
intent(inout) :: dtnew
1345 double precision,
intent(in) ::
dx^
d
1346 double precision,
intent(in) :: w(ixi^s,1:nw)
1347 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1362 end subroutine ffhd_get_dt
1364 subroutine ffhd_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
1366 integer,
intent(in) :: ixi^
l, ixo^
l
1367 double precision,
intent(in) :: qdt, dtfactor,x(ixi^s,1:
ndim)
1368 double precision,
intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw),w(ixi^s,1:nw)
1371 end subroutine ffhd_add_source_geom
1373 function ffhd_kin_en_origin(w, ixI^L, ixO^L, inv_rho)
result(ke)
1375 integer,
intent(in) :: ixi^
l, ixo^
l
1376 double precision,
intent(in) :: w(ixi^s, nw)
1377 double precision :: ke(ixo^s)
1378 double precision,
intent(in),
optional :: inv_rho(ixo^s)
1380 if(
present(inv_rho))
then
1381 ke(ixo^s)=0.5d0*w(ixo^s,
mom(1))**2*inv_rho(ixo^s)
1383 ke(ixo^s)=0.5d0*w(ixo^s,
mom(1))**2/w(ixo^s,
rho_)
1385 end function ffhd_kin_en_origin
1387 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1390 integer,
intent(in) :: ixi^
l, ixo^
l
1391 double precision,
intent(in) :: w(ixi^s,1:nw)
1392 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1393 double precision,
intent(out):: rfactor(ixi^s)
1394 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1397 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)
1398 end subroutine rfactor_from_temperature_ionization
1400 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1402 integer,
intent(in) :: ixi^
l, ixo^
l
1403 double precision,
intent(in) :: w(ixi^s,1:nw)
1404 double precision,
intent(in) :: x(ixi^s,1:
ndim)
1405 double precision,
intent(out):: rfactor(ixi^s)
1408 end subroutine rfactor_from_constant_ionization
1410 subroutine add_hypertc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1412 integer,
intent(in) :: ixi^
l,ixo^
l
1413 double precision,
intent(in) :: qdt
1414 double precision,
dimension(ixI^S,1:ndim),
intent(in) :: x
1415 double precision,
dimension(ixI^S,1:nw),
intent(in) :: wct,wctprim
1416 double precision,
dimension(ixI^S,1:nw),
intent(inout) :: w
1418 double precision,
dimension(ixI^S) :: te
1419 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1422 te(ixi^s)=wctprim(ixi^s,
p_)/wct(ixi^s,
rho_)
1426 do ix2=ixomin2,ixomax2
1427 do ix1=ixomin1,ixomax1
1434 sigma_t7=sigma_t5*te(ix^
d)
1438 sigma_t7=sigma_t5*te(ix^
d)
1440 sigmat5_bgradt=sigma_t5*(&
1441 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)&
1442 +
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))
1444 f_sat=one/(one+abs(sigmat5_bgradt))/(1.5d0*wct(ix^
d,
rho_)*(
ffhd_gamma*wctprim(ix^
d,
p_)/wct(ix^
d,
rho_))**1.5d0)
1446 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^
d,
q_))/tau
1448 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(sigmat5_bgradt+wct(ix^
d,
q_))/&
1455 do ix3=ixomin3,ixomax3
1456 do ix2=ixomin2,ixomax2
1457 do ix1=ixomin1,ixomax1
1464 sigma_t7=sigma_t5*te(ix^
d)
1468 sigma_t7=sigma_t5*te(ix^
d)
1470 sigmat5_bgradt=sigma_t5*(&
1471 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)&
1472 +
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)&
1473 +
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))
1475 f_sat=one/(one+abs(sigmat5_bgradt))/(1.5d0*wct(ix^
d,
rho_)*(
ffhd_gamma*wctprim(ix^
d,
p_)/wct(ix^
d,
rho_))**1.5d0)
1477 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^
d,
q_))/tau
1479 w(ix^
d,
q_)=w(ix^
d,
q_)-qdt*(sigmat5_bgradt+wct(ix^
d,
q_))/&
1486 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_viscosity
Whether viscosity is added.
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, rhov, 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 for HD and MHD
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine cooling_get_dt(w, ixil, ixol, dtnew, dxd, x, fl)
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 1) in amrvac.par in sts_list set the following parameters which have...
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 explicut timestep for the TC (mhd implementation)
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
The module add viscous source terms and check time step.
subroutine viscosity_add_source(qdt, ixil, ixol, wct, w, x, energy, qsourcesplit, active)
subroutine viscosity_init(phys_wider_stencil)
Initialize the module.
subroutine viscosity_get_dt(w, ixil, ixol, dtnew, dxd, x)