MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
mod_thermal_conduction.t
Go to the documentation of this file.
1 !> Thermal conduction for HD and MHD or RHD and RMHD or twofl (plasma-neutral) module
2 !> Adaptation of mod_thermal_conduction for the mod_supertimestepping
3 !>
4 !> The TC is set by calling
5 !> tc_init_params()
6 !>
7 !> Organized such that it can call either isotropic (HD) or anisotropic (MHD) variants
8 !> it adds a heat conduction source to each energy equation
9 !> and can be recycled within a multi-fluid context (such as plasma-neutral twofl module)
10 !>
11 !>
12 !> 10.07.2011 developed by Chun Xia and Rony Keppens
13 !> 01.09.2012 moved to modules folder by Oliver Porth
14 !> 13.10.2013 optimized further by Chun Xia
15 !> 12.03.2014 implemented RKL2 super timestepping scheme to reduce iterations
16 !> and improve stability and accuracy up to second order in time by Chun Xia.
17 !> 23.08.2014 implemented saturation and perpendicular TC by Chun Xia
18 !> 12.01.2017 modulized by Chun Xia
19 !> adapted by Beatrice Popescu to twofluid settings
20 !> 06.09.2024 cleaned up for use in rhd and rmhd modules (Nishant Narechania and Rony Keppens)
21 !>
22 !> PURPOSE:
23 !> IN MHD ADD THE HEAT CONDUCTION SOURCE TO THE ENERGY EQUATION
24 !> S=DIV(KAPPA_i,j . GRAD_j T)
25 !> where KAPPA_i,j = tc_k_para b_i b_j + tc_k_perp (I - b_i b_j)
26 !> b_i b_j = B_i B_j / B**2, I is the unit matrix, and i, j= 1, 2, 3 for 3D
27 !> IN HD ADD THE HEAT CONDUCTION SOURCE TO THE ENERGY EQUATION
28 !> S=DIV(tc_k_para . GRAD T)
29 !> USAGE:
30 !> 1. in mod_usr.t -> subroutine usr_init(), add
31 !> unit_length=your length unit
32 !> unit_numberdensity=your number density unit
33 !> unit_velocity=your velocity unit
34 !> unit_temperature=your temperature unit
35 !> before call (m)hd_activate()
36 !> 2. to switch on thermal conduction in the (r)(m)hd_list of amrvac.par add:
37 !> (r)(m)hd_thermal_conduction=.true.
38 !> 3. in the tc_list of amrvac.par :
39 !> tc_perpendicular=.true. ! (default .false.) turn on thermal conduction perpendicular to magnetic field
40 !> tc_saturate=.true. ! (default .false. ) turn on thermal conduction saturate effect
41 !> tc_slope_limiter='MC' ! choose limiter for slope-limited anisotropic thermal conduction in MHD
42 !> note: twofl_list incorporates instances for charges and neutrals
43 
45  use mod_global_parameters, only: std_len
46  use mod_geometry
47  use mod_comm_lib, only: mpistop
48  implicit none
49 
50  !> The adiabatic index
51  double precision :: tc_gamma_1
52 
53  abstract interface
54  subroutine get_var_subr(w,x,ixI^L,ixO^L,res)
56  integer, intent(in) :: ixI^L, ixO^L
57  double precision, intent(in) :: w(ixI^S,nw)
58  double precision, intent(in) :: x(ixI^S,1:ndim)
59  double precision, intent(out):: res(ixI^S)
60  end subroutine get_var_subr
61  end interface
62 
63  type tc_fluid
64 
65  procedure(get_var_subr), pointer, nopass :: get_rho => null()
66  procedure(get_var_subr), pointer, nopass :: get_rho_equi => null()
67  procedure(get_var_subr), pointer,nopass :: get_temperature_from_eint => null()
68  procedure(get_var_subr), pointer,nopass :: get_temperature_from_conserved => null()
69  procedure(get_var_subr), pointer,nopass :: get_temperature_equi => null()
70  !> Indices of the variables
71  integer :: e_=-1
72  !> Index of cut off temperature for TRAC
73  integer :: tcoff_
74  ! if has_equi = .true. get_temperature_equi and get_rho_equi have to be set
75  logical :: has_equi=.false.
76 
77  ! BEGIN the following are read from param file or set in tc_read_hd_params or tc_read_mhd_params
78  !> Coefficient of thermal conductivity (parallel to magnetic field)
79  double precision :: tc_k_para
80 
81  !> Coefficient of thermal conductivity perpendicular to magnetic field
82  double precision :: tc_k_perp
83 
84  !> Name of slope limiter for transverse component of thermal flux
85  integer :: tc_slope_limiter
86 
87  !> Logical switch for test constant conductivity
88  logical :: tc_constant=.false.
89 
90  !> Calculate thermal conduction perpendicular to magnetic field (.true.) or not (.false.)
91  logical :: tc_perpendicular=.false.
92 
93  !> Consider thermal conduction saturation effect (.true.) or not (.false.)
94  logical :: tc_saturate=.false.
95  ! END the following are read from param file or set in tc_read_hd_params or tc_read_mhd_params
96  end type tc_fluid
97 
98 
99  public :: tc_get_mhd_params
100  public :: tc_get_hd_params
101  public :: tc_get_ffhd_params
102  public :: get_tc_dt_mhd
103  public :: get_tc_dt_hd
104  public :: get_tc_dt_ffhd
105  public :: sts_set_source_tc_mhd
106  public :: sts_set_source_tc_hd
107  public :: sts_set_source_tc_ffhd
108 
109 contains
110 
111  subroutine tc_init_params(phys_gamma)
113  double precision, intent(in) :: phys_gamma
114 
115  tc_gamma_1=phys_gamma-1d0
116  end subroutine tc_init_params
117 
118  !> Init TC coefficients: MHD case
119  subroutine tc_get_mhd_params(fl,read_mhd_params)
121 
122  interface
123  subroutine read_mhd_params(fl)
125  import tc_fluid
126  type(tc_fluid), intent(inout) :: fl
127 
128  end subroutine read_mhd_params
129  end interface
130  type(tc_fluid), intent(inout) :: fl
131 
132  fl%tc_slope_limiter=1
133  fl%tc_k_para=0.d0
134  fl%tc_k_perp=0.d0
135 
136  !> Read tc module parameters from par file: MHD case
137  call read_mhd_params(fl)
138 
139  if(fl%tc_k_para==0.d0 .and. fl%tc_k_perp==0.d0) then
140  if(si_unit) then
141  ! Spitzer thermal conductivity with SI units
142  fl%tc_k_para=8.d-12*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
143  ! thermal conductivity perpendicular to magnetic field
144  fl%tc_k_perp=4.d-30*unit_numberdensity**2/unit_magneticfield**2/unit_temperature**3*fl%tc_k_para
145  else
146  ! Spitzer thermal conductivity with cgs units
147  fl%tc_k_para=8.d-7*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
148  ! thermal conductivity perpendicular to magnetic field
149  fl%tc_k_perp=4.d-10*unit_numberdensity**2/unit_magneticfield**2/unit_temperature**3*fl%tc_k_para
150  end if
151  if(mype .eq. 0) print*, "Spitzer MHD: par: ",fl%tc_k_para, &
152  " ,perp: ",fl%tc_k_perp
153  else
154  fl%tc_constant=.true.
155  end if
156 
157  end subroutine tc_get_mhd_params
158 
159  !> Init TC coefficients: HD case
160  subroutine tc_get_hd_params(fl,read_hd_params)
162 
163  interface
164  subroutine read_hd_params(fl)
166  import tc_fluid
167  type(tc_fluid), intent(inout) :: fl
168 
169  end subroutine read_hd_params
170  end interface
171  type(tc_fluid), intent(inout) :: fl
172 
173  fl%tc_k_para=0.d0
174 
175  !> Read tc parameters from par file: HD case
176  call read_hd_params(fl)
177 
178  if(fl%tc_k_para==0.d0 ) then
179  if(si_unit) then
180  ! Spitzer thermal conductivity with SI units
181  fl%tc_k_para=8.d-12*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
182  else
183  ! Spitzer thermal conductivity with cgs units
184  fl%tc_k_para=8.d-7*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
185  end if
186  if(mype .eq. 0) print*, "Spitzer HD par: ",fl%tc_k_para
187  end if
188 
189  end subroutine tc_get_hd_params
190 
191  subroutine tc_get_ffhd_params(fl,read_ffhd_params)
193 
194  interface
195  subroutine read_ffhd_params(fl)
197  import tc_fluid
198  type(tc_fluid), intent(inout) :: fl
199  end subroutine read_ffhd_params
200  end interface
201 
202  type(tc_fluid), intent(inout) :: fl
203 
204  fl%tc_k_para=0.d0
205  !> Read tc parameters from par file: ffHD case
206  call read_ffhd_params(fl)
207  if(fl%tc_k_para==0.d0 ) then
208  if(si_unit) then
209  ! Spitzer thermal conductivity with SI units
210  fl%tc_k_para=8.d-12*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
211  else
212  ! Spitzer thermal conductivity with cgs units
213  fl%tc_k_para=8.d-7*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
214  endif
215  if(mype .eq. 0) print*, "Spitzer ffHD par: ",fl%tc_k_para
216  else
217  fl%tc_constant=.true.
218  endif
219  end subroutine tc_get_ffhd_params
220 
221  !> Get the explicut timestep for the TC (mhd implementation)
222  function get_tc_dt_mhd(w,ixI^L,ixO^L,dx^D,x,fl) result(dtnew)
223  !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
224  !where tc_k_para_i=tc_k_para*B_i**2/B**2
225  !and T=p/rho
227 
228  type(tc_fluid), intent(in) :: fl
229  integer, intent(in) :: ixi^l, ixo^l
230  double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
231  double precision, intent(in) :: w(ixi^s,1:nw)
232  double precision :: dtnew
233 
234  double precision :: mf(ixo^s,1:ndim),te(ixi^s),b2(ixi^s),gradt(ixi^s)
235  double precision :: tmp2(ixo^s),tmp(ixo^s),hfs(ixo^s)
236  double precision :: dtdiff_tcond,maxtmp2
237  integer :: idims
238 
239  !temperature
240  call fl%get_temperature_from_conserved(w,x,ixi^l,ixi^l,te)
241 
242  ! B
243  if(b0field) then
244  mf(ixo^s,1:ndim)=w(ixo^s,iw_mag(1:ndim))+block%B0(ixo^s,1:ndim,0)
245  else
246  mf(ixo^s,1:ndim)=w(ixo^s,iw_mag(1:ndim))
247  end if
248  ! Bsquared
249  b2(ixo^s)=sum(mf(ixo^s,1:ndim)**2,dim=ndim+1)
250  ! B_i**2/B**2
251  where(b2(ixo^s)/=0.d0)
252  ^d&mf(ixo^s,^d)=mf(ixo^s,^d)**2/b2(ixo^s);
253  elsewhere
254  ^d&mf(ixo^s,^d)=1.d0;
255  end where
256 
257  dtnew=bigdouble
258  ! B2 is now density
259  call fl%get_rho(w,x,ixi^l,ixo^l,b2)
260 
261  !tc_k_para_i
262  if(fl%tc_constant) then
263  tmp(ixo^s)=fl%tc_k_para
264  else
265  if(fl%tc_saturate) then
266  ! Kannan 2016 MN 458, 410
267  ! 3^1.5*kB^2/(4*sqrt(pi)*e^4)
268  ! l_mfpe=3.d0**1.5d0*kB_cgs**2/(4.d0*sqrt(dpi)*e_cgs**4*37.d0)=7093.9239487765044d0
269  tmp2(ixo^s)=te(ixo^s)**2/b2(ixo^s)*7093.9239487765044d0*unit_temperature**2/(unit_numberdensity*unit_length)
270  hfs=0.d0
271  do idims=1,ndim
272  call gradient(te,ixi^l,ixo^l,idims,gradt)
273  hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)*sqrt(mf(ixo^s,idims))
274  end do
275  ! kappa=kappa_Spizer/(1+4.2*l_mfpe/(T/|gradT.b|))
276  tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)/(1.d0+4.2d0*tmp2(ixo^s)*dabs(hfs(ixo^s))/te(ixo^s))
277  else
278  ! kappa=kappa_Spizer
279  tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)
280  end if
281  end if
282 
283  if(slab_uniform) then
284  do idims=1,ndim
285  ! approximate thermal conduction flux: tc_k_para_i/rho/dx*B_i**2/B**2
286  tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*dxlevel(idims))
287  maxtmp2=maxval(tmp2(ixo^s))
288  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
289  dtdiff_tcond=dxlevel(idims)/(tc_gamma_1*maxtmp2+smalldouble)
290  ! limit the time step
291  dtnew=min(dtnew,dtdiff_tcond)
292  end do
293  else
294  do idims=1,ndim
295  ! approximate thermal conduction flux: tc_k_para_i/rho/dx*B_i**2/B**2
296  tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*block%ds(ixo^s,idims))
297  maxtmp2=maxval(tmp2(ixo^s)/block%ds(ixo^s,idims))
298  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
299  dtdiff_tcond=1.d0/(tc_gamma_1*maxtmp2+smalldouble)
300  ! limit the time step
301  dtnew=min(dtnew,dtdiff_tcond)
302  end do
303  end if
304  dtnew=dtnew/dble(ndim)
305  end function get_tc_dt_mhd
306 
307  !> anisotropic thermal conduction with slope limited symmetric scheme
308  !> Sharma 2007 Journal of Computational Physics 227, 123
309  subroutine sts_set_source_tc_mhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,fl)
311  use mod_fix_conserve
312  integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
313  double precision, intent(in) :: x(ixi^s,1:ndim)
314  double precision, intent(in) :: w(ixi^s,1:nw)
315  double precision, intent(inout) :: wres(ixi^s,1:nw)
316  double precision, intent(in) :: my_dt
317  logical, intent(in) :: fix_conserve_at_step
318  type(tc_fluid), intent(in) :: fl
319 
320  !! qd store the heat conduction energy changing rate
321  double precision :: qd(ixi^s)
322  double precision :: rho(ixi^s),te(ixi^s)
323  double precision :: qvec(ixi^s,1:ndim)
324  double precision, allocatable, dimension(:^D&,:) :: qvec_equi
325 
326  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
327  double precision :: alpha,dxinv(ndim)
328  integer :: idims,idir,ix^d,ix^l,ixc^l,ixa^l,ixb^l
329 
330  ! coefficient of limiting on normal component
331  if(ndim<3) then
332  alpha=0.75d0
333  else
334  alpha=0.85d0
335  end if
336  ix^l=ixo^l^ladd1;
337 
338  dxinv=1.d0/dxlevel
339 
340  call fl%get_temperature_from_eint(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
341  call fl%get_rho(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
342  call set_source_tc_mhd(ixi^l,ixo^l,w,x,fl,qvec,rho,te,alpha)
343  if(fl%has_equi) then
344  allocate(qvec_equi(ixi^s,1:ndim))
345  call fl%get_temperature_equi(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
346  call fl%get_rho_equi(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
347  call set_source_tc_mhd(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te,alpha)
348  do idims=1,ndim
349  qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
350  end do
351  deallocate(qvec_equi)
352  endif
353 
354  qd=0.d0
355  if(slab_uniform) then
356  do idims=1,ndim
357  qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
358  ixb^l=ixo^l-kr(idims,^d);
359  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
360  end do
361  else
362  do idims=1,ndim
363  qvec(ix^s,idims)=qvec(ix^s,idims)*block%surfaceC(ix^s,idims)
364  ixb^l=ixo^l-kr(idims,^d);
365  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
366  end do
367  qd(ixo^s)=qd(ixo^s)/block%dvolume(ixo^s)
368  end if
369 
370  if(fix_conserve_at_step) then
371  allocate(fluxall(ixi^s,1,1:ndim))
372  fluxall(ixi^s,1,1:ndim)=my_dt*qvec(ixi^s,1:ndim)
373  call store_flux(igrid,fluxall,1,ndim,nflux)
374  deallocate(fluxall)
375  end if
376 
377  wres(ixo^s,fl%e_)=qd(ixo^s)
378  end subroutine sts_set_source_tc_mhd
379 
380  subroutine set_source_tc_mhd(ixI^L,ixO^L,w,x,fl,qvec,rho,Te,alpha)
382  integer, intent(in) :: ixI^L, ixO^L
383  double precision, intent(in) :: x(ixI^S,1:ndim)
384  double precision, intent(in) :: w(ixI^S,1:nw)
385  type(tc_fluid), intent(in) :: fl
386  double precision, intent(in) :: rho(ixI^S),Te(ixI^S)
387  double precision, intent(in) :: alpha
388  double precision, intent(out) :: qvec(ixI^S,1:ndim)
389  !! qd store the heat conduction energy changing rate
390  double precision :: qd(ixI^S)
391  double precision, dimension(ixI^S,1:ndim) :: mf,Bc,Bcf
392  double precision, dimension(ixI^S,1:ndim) :: gradT
393  double precision, dimension(ixI^S) :: ka,kaf,ke,kef,qdd,qe,Binv,minq,maxq,Bnorm
394  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
395  integer, dimension(ndim) :: lowindex
396  integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L
397 
398  ix^l=ixo^l^ladd1;
399 
400  ! T gradient at cell faces
401  ! B vector
402  if(b0field) then
403  mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))+block%B0(ixi^s,1:ndim,0)
404  else
405  mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))
406  end if
407  ! |B|
408  binv(ix^s)=dsqrt(sum(mf(ix^s,1:ndim)**2,dim=ndim+1))
409  where(binv(ix^s)/=0.d0)
410  binv(ix^s)=1.d0/binv(ix^s)
411  elsewhere
412  binv(ix^s)=bigdouble
413  end where
414  ! b unit vector: magnetic field direction vector
415  do idims=1,ndim
416  mf(ix^s,idims)=mf(ix^s,idims)*binv(ix^s)
417  end do
418  ! ixC is cell-corner index
419  ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
420  ! b unit vector at cell corner
421  bc=0.d0
422  {do ix^db=0,1\}
423  ixamin^d=ixcmin^d+ix^d;
424  ixamax^d=ixcmax^d+ix^d;
425  bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)+mf(ixa^s,1:ndim)
426  {end do\}
427  bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)*0.5d0**ndim
428  ! T gradient at cell faces
429  do idims=1,ndim
430  ixbmin^d=ixmin^d;
431  ixbmax^d=ixmax^d-kr(idims,^d);
432  call gradientc(te,ixi^l,ixb^l,idims,minq)
433  gradt(ixb^s,idims)=minq(ixb^s)
434  end do
435  if(fl%tc_constant) then
436  if(fl%tc_perpendicular) then
437  ka(ixc^s)=fl%tc_k_para-fl%tc_k_perp
438  ke(ixc^s)=fl%tc_k_perp
439  else
440  ka(ixc^s)=fl%tc_k_para
441  end if
442  else
443  ! conductivity at cell center
444  if(phys_trac) then
445  minq(ix^s)=te(ix^s)
446  where(minq(ix^s) < block%wextra(ix^s,fl%Tcoff_))
447  minq(ix^s)=block%wextra(ix^s,fl%Tcoff_)
448  end where
449  minq(ix^s)=fl%tc_k_para*sqrt(minq(ix^s)**5)
450  else
451  minq(ix^s)=fl%tc_k_para*sqrt(te(ix^s)**5)
452  end if
453  ka=0.d0
454  {do ix^db=0,1\}
455  ixbmin^d=ixcmin^d+ix^d;
456  ixbmax^d=ixcmax^d+ix^d;
457  ka(ixc^s)=ka(ixc^s)+minq(ixb^s)
458  {end do\}
459  ! cell corner conductivity
460  ka(ixc^s)=0.5d0**ndim*ka(ixc^s)
461  ! compensate with perpendicular conductivity
462  if(fl%tc_perpendicular) then
463  minq(ix^s)=fl%tc_k_perp*rho(ix^s)**2*binv(ix^s)**2/dsqrt(te(ix^s))
464  ke=0.d0
465  {do ix^db=0,1\}
466  ixbmin^d=ixcmin^d+ix^d;
467  ixbmax^d=ixcmax^d+ix^d;
468  ke(ixc^s)=ke(ixc^s)+minq(ixb^s)
469  {end do\}
470  ! cell corner conductivity: k_parallel-k_perpendicular
471  ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
472  where(ke(ixc^s)<ka(ixc^s))
473  ka(ixc^s)=ka(ixc^s)-ke(ixc^s)
474  elsewhere
475  ke(ixc^s)=ka(ixc^s)
476  ka(ixc^s)=0.d0
477  end where
478  end if
479  end if
480  if(fl%tc_slope_limiter==0) then
481  ! calculate thermal conduction flux with symmetric scheme
482  do idims=1,ndim
483  !qd corner values
484  qd=0.d0
485  {do ix^db=0,1 \}
486  if({ ix^d==0 .and. ^d==idims | .or.}) then
487  ixbmin^d=ixcmin^d+ix^d;
488  ixbmax^d=ixcmax^d+ix^d;
489  qd(ixc^s)=qd(ixc^s)+gradt(ixb^s,idims)
490  end if
491  {end do\}
492  ! temperature gradient at cell corner
493  qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
494  end do
495  ! b grad T at cell corner
496  qd(ixc^s)=sum(qvec(ixc^s,1:ndim)*bc(ixc^s,1:ndim),dim=ndim+1)
497  do idims=1,ndim
498  ! TC flux at cell corner
499  gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qd(ixc^s)
500  if(fl%tc_perpendicular) gradt(ixc^s,idims)=gradt(ixc^s,idims)+ke(ixc^s)*qvec(ixc^s,idims)
501  end do
502  ! TC flux at cell face
503  qvec=0.d0
504  do idims=1,ndim
505  ixb^l=ixo^l-kr(idims,^d);
506  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
507  {do ix^db=0,1 \}
508  if({ ix^d==0 .and. ^d==idims | .or.}) then
509  ixbmin^d=ixamin^d-ix^d;
510  ixbmax^d=ixamax^d-ix^d;
511  qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
512  end if
513  {end do\}
514  qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
515  if(fl%tc_saturate) then
516  ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
517  ! unsigned saturated TC flux = 5 phi rho c**3, c is isothermal sound speed
518  bcf=0.d0
519  {do ix^db=0,1 \}
520  if({ ix^d==0 .and. ^d==idims | .or.}) then
521  ixbmin^d=ixamin^d-ix^d;
522  ixbmax^d=ixamax^d-ix^d;
523  bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
524  end if
525  {end do\}
526  ! averaged b at face centers
527  bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
528  ixb^l=ixa^l+kr(idims,^d);
529  qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bcf(ixa^s,idims))
530  {do ix^db=ixamin^db,ixamax^db\}
531  if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
532  qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
533  end if
534  {end do\}
535  end if
536  end do
537  else
538  ! calculate thermal conduction flux with slope-limited symmetric scheme
539  qvec=0.d0
540  do idims=1,ndim
541  ixb^l=ixo^l-kr(idims,^d);
542  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
543  ! calculate normal of magnetic field
544  ixb^l=ixa^l+kr(idims,^d);
545  bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
546  bcf=0.d0
547  kaf=0.d0
548  kef=0.d0
549  {do ix^db=0,1 \}
550  if({ ix^d==0 .and. ^d==idims | .or.}) then
551  ixbmin^d=ixamin^d-ix^d;
552  ixbmax^d=ixamax^d-ix^d;
553  bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)+bc(ixb^s,1:ndim)
554  kaf(ixa^s)=kaf(ixa^s)+ka(ixb^s)
555  if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)+ke(ixb^s)
556  end if
557  {end do\}
558  ! averaged b at face centers
559  bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)*0.5d0**(ndim-1)
560  ! averaged thermal conductivity at face centers
561  kaf(ixa^s)=kaf(ixa^s)*0.5d0**(ndim-1)
562  if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)*0.5d0**(ndim-1)
563  ! limited normal component
564  minq(ixa^s)=min(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
565  maxq(ixa^s)=max(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
566  ! eq (19)
567  !corner values of gradT
568  qdd=0.d0
569  {do ix^db=0,1 \}
570  if({ ix^d==0 .and. ^d==idims | .or.}) then
571  ixbmin^d=ixcmin^d+ix^d;
572  ixbmax^d=ixcmax^d+ix^d;
573  qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
574  end if
575  {end do\}
576  ! temperature gradient at cell corner
577  qdd(ixc^s)=qdd(ixc^s)*0.5d0**(ndim-1)
578  ! eq (21)
579  qe=0.d0
580  {do ix^db=0,1 \}
581  qd(ixc^s)=qdd(ixc^s)
582  if({ ix^d==0 .and. ^d==idims | .or.}) then
583  ixbmin^d=ixamin^d-ix^d;
584  ixbmax^d=ixamax^d-ix^d;
585  where(qd(ixb^s)<=minq(ixa^s))
586  qd(ixb^s)=minq(ixa^s)
587  elsewhere(qd(ixb^s)>=maxq(ixa^s))
588  qd(ixb^s)=maxq(ixa^s)
589  end where
590  qvec(ixa^s,idims)=qvec(ixa^s,idims)+bc(ixb^s,idims)**2*qd(ixb^s)
591  if(fl%tc_perpendicular) qe(ixa^s)=qe(ixa^s)+qd(ixb^s)
592  end if
593  {end do\}
594  qvec(ixa^s,idims)=kaf(ixa^s)*qvec(ixa^s,idims)*0.5d0**(ndim-1)
595  ! add normal flux from perpendicular conduction
596  if(fl%tc_perpendicular) qvec(ixa^s,idims)=qvec(ixa^s,idims)+kef(ixa^s)*qe(ixa^s)*0.5d0**(ndim-1)
597  ! limited transverse component, eq (17)
598  ixbmin^d=ixamin^d;
599  ixbmax^d=ixamax^d+kr(idims,^d);
600  do idir=1,ndim
601  if(idir==idims) cycle
602  qd(ixi^s)=slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
603  qd(ixi^s)=slope_limiter(qd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
604  qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qd(ixa^s)
605  end do
606 
607  if(fl%tc_saturate) then
608  ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
609  ! unsigned saturated TC flux = 5 phi rho c**3, c=sqrt(p/rho) is isothermal sound speed, phi=1.1
610  ixb^l=ixa^l+kr(idims,^d);
611  qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bnorm(ixa^s))
612  {do ix^db=ixamin^db,ixamax^db\}
613  if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
614  !write(*,*) 'it',it,qvec(ix^D,idims),qd(ix^D),' TC saturated at ',&
615  !x(ix^D,:),' rho',rho(ix^D),' Te',Te(ix^D)
616  qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
617  end if
618  {end do\}
619  end if
620  end do
621  end if
622  end subroutine set_source_tc_mhd
623 
624  function slope_limiter(f,ixI^L,ixO^L,idims,pm,tc_slope_limiter) result(lf)
626  integer, intent(in) :: ixi^l, ixo^l, idims, pm
627  double precision, intent(in) :: f(ixi^s)
628  double precision :: lf(ixi^s)
629  integer, intent(in) :: tc_slope_limiter
630 
631  double precision :: signf(ixi^s)
632  integer :: ixb^l
633 
634  ixb^l=ixo^l+pm*kr(idims,^d);
635  signf(ixo^s)=sign(1.d0,f(ixo^s))
636  select case(tc_slope_limiter)
637  case(1)
638  ! 'MC' montonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
639  lf(ixo^s)=two*signf(ixo^s)* &
640  max(zero,min(dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s),&
641  signf(ixo^s)*quarter*(f(ixb^s)+f(ixo^s))))
642  case(2)
643  ! 'minmod' limiter
644  lf(ixo^s)=signf(ixo^s)*max(0.d0,min(abs(f(ixo^s)),signf(ixo^s)*f(ixb^s)))
645  case(3)
646  ! 'superbee' Roes superbee limiter (eq.3.51i)
647  lf(ixo^s)=signf(ixo^s)* &
648  max(zero,min(two*dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s)),&
649  min(dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s)))
650  case(4)
651  ! 'koren' Barry Koren Right variant
652  lf(ixo^s)=signf(ixo^s)* &
653  max(zero,min(two*dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s),&
654  (two*f(ixb^s)*signf(ixo^s)+dabs(f(ixo^s)))*third))
655  case default
656  call mpistop("Unknown slope limiter for thermal conduction")
657  end select
658  end function slope_limiter
659 
660  !> Calculate gradient of a scalar q at cell interfaces in direction idir
661  subroutine gradientc(q,ixI^L,ixO^L,idir,gradq)
663  integer, intent(in) :: ixI^L, ixO^L, idir
664  double precision, intent(in) :: q(ixI^S)
665  double precision, intent(inout) :: gradq(ixI^S)
666  integer :: jxO^L
667 
668  associate(x=>block%x)
669 
670  jxo^l=ixo^l+kr(idir,^d);
671  select case(coordinate)
672  case(cartesian)
673  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/dxlevel(idir)
674  case(cartesian_stretched)
675  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
676  case(spherical)
677  select case(idir)
678  case(1)
679  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,1)-x(ixo^s,1))
680  {^nooned
681  case(2)
682  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/( (x(jxo^s,2)-x(ixo^s,2))*x(ixo^s,1) )
683  }
684  {^ifthreed
685  case(3)
686  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/( (x(jxo^s,3)-x(ixo^s,3))*x(ixo^s,1)*dsin(x(ixo^s,2)) )
687  }
688  end select
689  case(cylindrical)
690  if(idir==phi_) then
691  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/((x(jxo^s,phi_)-x(ixo^s,phi_))*x(ixo^s,r_))
692  else
693  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
694  end if
695  case default
696  call mpistop('Unknown geometry')
697  end select
698 
699  end associate
700  end subroutine gradientc
701 
702  !> Get the explicit timestep for the TC (hd implementation)
703  function get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x,fl) result(dtnew)
704  ! Check diffusion time limit dt < dx_i**2 / ((gamma-1)*tc_k_para_i/rho)
706 
707  integer, intent(in) :: ixi^l, ixo^l
708  double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
709  double precision, intent(in) :: w(ixi^s,1:nw)
710  type(tc_fluid), intent(in) :: fl
711  double precision :: dtnew
712 
713  double precision :: tmp(ixo^s),tmp2(ixo^s),te(ixi^s),rho(ixi^s),hfs(ixo^s),gradt(ixi^s)
714  double precision :: dtdiff_tcond,maxtmp2
715  integer :: idim
716 
717  call fl%get_temperature_from_conserved(w,x,ixi^l,ixi^l,te)
718  call fl%get_rho(w,x,ixi^l,ixo^l,rho)
719 
720  if(fl%tc_saturate) then
721  ! Kannan 2016 MN 458, 410
722  ! 3^1.5*kB^2/(4*sqrt(pi)*e^4)
723  ! l_mfpe=3.d0**1.5d0*kB_cgs**2/(4.d0*sqrt(dpi)*e_cgs**4*37.d0)=7093.9239487765044d0
724  tmp2(ixo^s)=te(ixo^s)**2/rho(ixo^s)*7093.9239487765044d0*unit_temperature**2/(unit_numberdensity*unit_length)
725  hfs=0.d0
726  do idim=1,ndim
727  call gradient(te,ixi^l,ixo^l,idim,gradt)
728  hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)**2
729  end do
730  ! kappa=kappa_Spizer/(1+4.2*l_mfpe/(T/|gradT|))
731  tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/(rho(ixo^s)*(1.d0+4.2d0*tmp2(ixo^s)*sqrt(hfs(ixo^s))/te(ixo^s)))
732  else
733  tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/rho(ixo^s)
734  end if
735  dtnew = bigdouble
736 
737  if(slab_uniform) then
738  do idim=1,ndim
739  ! approximate thermal conduction flux: tc_k_para_i/rho/dx
740  tmp2(ixo^s)=tmp(ixo^s)/dxlevel(idim)
741  maxtmp2=maxval(tmp2(ixo^s))
742  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho)
743  dtdiff_tcond=dxlevel(idim)/(tc_gamma_1*maxtmp2+smalldouble)
744  ! limit the time step
745  dtnew=min(dtnew,dtdiff_tcond)
746  end do
747  else
748  do idim=1,ndim
749  ! approximate thermal conduction flux: tc_k_para_i/rho/dx
750  tmp2(ixo^s)=tmp(ixo^s)/block%ds(ixo^s,idim)
751  maxtmp2=maxval(tmp2(ixo^s)/block%ds(ixo^s,idim))
752  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
753  dtdiff_tcond=1.d0/(tc_gamma_1*maxtmp2+smalldouble)
754  ! limit the time step
755  dtnew=min(dtnew,dtdiff_tcond)
756  end do
757  end if
758  dtnew=dtnew/dble(ndim)
759  end function get_tc_dt_hd
760 
761  subroutine sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,fl)
763  use mod_fix_conserve
764 
765  integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
766  double precision, intent(in) :: x(ixi^s,1:ndim)
767  double precision, intent(in) :: w(ixi^s,1:nw)
768  double precision, intent(inout) :: wres(ixi^s,1:nw)
769  double precision, intent(in) :: my_dt
770  logical, intent(in) :: fix_conserve_at_step
771  type(tc_fluid), intent(in) :: fl
772 
773  double precision :: te(ixi^s),rho(ixi^s)
774  double precision :: qvec(ixi^s,1:ndim),qd(ixi^s)
775  double precision, allocatable, dimension(:^D&,:) :: qvec_equi
776  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
777 
778  double precision :: dxinv(ndim)
779  integer :: idims,ix^l,ixb^l
780 
781  ix^l=ixo^l^ladd1;
782 
783  dxinv=1.d0/dxlevel
784 
785  !calculate Te in whole domain (+ghosts)
786  call fl%get_temperature_from_eint(w, x, ixi^l, ixi^l, te)
787  call fl%get_rho(w, x, ixi^l, ixi^l, rho)
788  call set_source_tc_hd(ixi^l,ixo^l,w,x,fl,qvec,rho,te)
789  if(fl%has_equi) then
790  allocate(qvec_equi(ixi^s,1:ndim))
791  call fl%get_temperature_equi(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
792  call fl%get_rho_equi(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
793  call set_source_tc_hd(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te)
794  do idims=1,ndim
795  qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
796  end do
797  deallocate(qvec_equi)
798  endif
799 
800  qd=0.d0
801  if(slab_uniform) then
802  do idims=1,ndim
803  qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
804  ixb^l=ixo^l-kr(idims,^d);
805  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
806  end do
807  else
808  do idims=1,ndim
809  qvec(ix^s,idims)=qvec(ix^s,idims)*block%surfaceC(ix^s,idims)
810  ixb^l=ixo^l-kr(idims,^d);
811  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
812  end do
813  qd(ixo^s)=qd(ixo^s)/block%dvolume(ixo^s)
814  end if
815 
816  if(fix_conserve_at_step) then
817  allocate(fluxall(ixi^s,1,1:ndim))
818  fluxall(ixi^s,1,1:ndim)=my_dt*qvec(ixi^s,1:ndim)
819  call store_flux(igrid,fluxall,1,ndim,nflux)
820  deallocate(fluxall)
821  end if
822  wres(ixo^s,fl%e_)=qd(ixo^s)
823  end subroutine sts_set_source_tc_hd
824 
825  subroutine set_source_tc_hd(ixI^L,ixO^L,w,x,fl,qvec,rho,Te)
827  integer, intent(in) :: ixI^L, ixO^L
828  double precision, intent(in) :: x(ixI^S,1:ndim)
829  double precision, intent(in) :: w(ixI^S,1:nw)
830  type(tc_fluid), intent(in) :: fl
831  double precision, intent(in) :: Te(ixI^S),rho(ixI^S)
832  double precision, intent(out) :: qvec(ixI^S,1:ndim)
833  double precision :: gradT(ixI^S,1:ndim),ke(ixI^S),qd(ixI^S)
834  integer :: idims,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixD^L
835 
836  ix^l=ixo^l^ladd1;
837  ! ixC is cell-corner index
838  ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
839 
840  ! calculate thermal conduction flux with symmetric scheme
841  ! T gradient (central difference) at cell corners
842  do idims=1,ndim
843  ixbmin^d=ixmin^d;
844  ixbmax^d=ixmax^d-kr(idims,^d);
845  call gradientc(te,ixi^l,ixb^l,idims,ke)
846  qd=0.d0
847  {do ix^db=0,1 \}
848  if({ix^d==0 .and. ^d==idims |.or. }) then
849  ixbmin^d=ixcmin^d+ix^d;
850  ixbmax^d=ixcmax^d+ix^d;
851  qd(ixc^s)=qd(ixc^s)+ke(ixb^s)
852  end if
853  {end do\}
854  ! temperature gradient at cell corner
855  qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
856  end do
857  ! conductivity at cell center
858  if(phys_trac) then
859  ! transition region adaptive conduction
860  where(te(ix^s) < block%wextra(ix^s,fl%Tcoff_))
861  qd(ix^s)=fl%tc_k_para*dsqrt(block%wextra(ix^s,fl%Tcoff_))**5
862  else where
863  qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s))**5
864  end where
865  else
866  qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s))**5
867  end if
868  ke=0.d0
869  {do ix^db=0,1\}
870  ixbmin^d=ixcmin^d+ix^d;
871  ixbmax^d=ixcmax^d+ix^d;
872  ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
873  {end do\}
874  ! cell corner conductivity
875  ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
876  ! cell corner conduction flux
877  do idims=1,ndim
878  gradt(ixc^s,idims)=ke(ixc^s)*qvec(ixc^s,idims)
879  end do
880 
881  if(fl%tc_saturate) then
882  ! consider saturation with unsigned saturated TC flux = 5 phi rho c**3
883  ! saturation flux at cell center
884  qd(ix^s)=5.5d0*rho(ix^s)*dsqrt(te(ix^s)**3)
885  !cell corner values of qd in ke
886  ke=0.d0
887  {do ix^db=0,1\}
888  ixbmin^d=ixcmin^d+ix^d;
889  ixbmax^d=ixcmax^d+ix^d;
890  ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
891  {end do\}
892  ! cell corner saturation flux
893  ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
894  ! magnitude of cell corner conduction flux
895  qd(ixc^s)=norm2(gradt(ixc^s,:),dim=ndim+1)
896  {do ix^db=ixcmin^db,ixcmax^db\}
897  if(qd(ix^d)>ke(ix^d)) then
898  ke(ix^d)=ke(ix^d)/(qd(ix^d)+smalldouble)
899  do idims=1,ndim
900  gradt(ix^d,idims)=ke(ix^d)*gradt(ix^d,idims)
901  end do
902  end if
903  {end do\}
904  end if
905 
906  ! conduction flux at cell face
907  qvec=0.d0
908  do idims=1,ndim
909  ixb^l=ixo^l-kr(idims,^d);
910  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
911  {do ix^db=0,1 \}
912  if({ ix^d==0 .and. ^d==idims | .or.}) then
913  ixbmin^d=ixamin^d-ix^d;
914  ixbmax^d=ixamax^d-ix^d;
915  qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
916  end if
917  {end do\}
918  qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
919  end do
920  end subroutine set_source_tc_hd
921 
922  !> Get the explicut timestep for the TC (ffhd implementation)
923  function get_tc_dt_ffhd(w,ixI^L,ixO^L,dx^D,x,fl) result(dtnew)
924  !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
925  !where tc_k_para_i=tc_k_para*B_i**2/B**2
926  !and T=p/rho
928  type(tc_fluid), intent(in) :: fl
929  integer, intent(in) :: ixi^l, ixo^l
930  double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
931  double precision, intent(in) :: w(ixi^s,1:nw)
932  double precision :: dtnew
933  double precision :: mf(ixo^s,1:ndim),te(ixi^s),b2(ixi^s),gradt(ixi^s)
934  double precision :: tmp2(ixo^s),tmp(ixo^s),hfs(ixo^s)
935  double precision :: dtdiff_tcond,maxtmp2
936  integer :: idims
937 
938  !temperature
939  call fl%get_temperature_from_conserved(w,x,ixi^l,ixi^l,te)
940  mf(ixo^s,1:ndim)=(block%B0(ixo^s,1:ndim,0))**2
941 
942  dtnew=bigdouble
943  ! B2 is now density
944  call fl%get_rho(w,x,ixi^l,ixo^l,b2)
945 
946  !tc_k_para_i
947  if(fl%tc_constant) then
948  tmp(ixo^s)=fl%tc_k_para
949  else
950  if(fl%tc_saturate) then
951  ! Kannan 2016 MN 458, 410
952  ! 3^1.5*kB^2/(4*sqrt(pi)*e^4)
953  ! l_mfpe=3.d0**1.5d0*kB_cgs**2/(4.d0*sqrt(dpi)*e_cgs**4*37.d0)=7093.9239487765044d0
954  tmp2(ixo^s)=te(ixo^s)**2/b2(ixo^s)*7093.9239487765044d0*unit_temperature**2/(unit_numberdensity*unit_length)
955  hfs=0.d0
956  do idims=1,ndim
957  call gradient(te,ixi^l,ixo^l,idims,gradt)
958  hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)*sqrt(mf(ixo^s,idims))
959  end do
960  ! kappa=kappa_Spizer/(1+4.2*l_mfpe/(T/|gradT.b|))
961  tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)/(1.d0+4.2d0*tmp2(ixo^s)*dabs(hfs(ixo^s))/te(ixo^s))
962  else
963  ! kappa=kappa_Spizer
964  tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)
965  end if
966  end if
967 
968  if(slab_uniform) then
969  do idims=1,ndim
970  ! approximate thermal conduction flux: tc_k_para_i/rho/dx*B_i**2/B**2
971  tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*dxlevel(idims))
972  maxtmp2=maxval(tmp2(ixo^s))
973  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
974  dtdiff_tcond=dxlevel(idims)/(tc_gamma_1*maxtmp2+smalldouble)
975  ! limit the time step
976  dtnew=min(dtnew,dtdiff_tcond)
977  end do
978  else
979  do idims=1,ndim
980  ! approximate thermal conduction flux: tc_k_para_i/rho/dx*B_i**2/B**2
981  tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*block%ds(ixo^s,idims))
982  maxtmp2=maxval(tmp2(ixo^s)/block%ds(ixo^s,idims))
983  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
984  dtdiff_tcond=1.d0/(tc_gamma_1*maxtmp2+smalldouble)
985  ! limit the time step
986  dtnew=min(dtnew,dtdiff_tcond)
987  end do
988  end if
989  dtnew=dtnew/dble(ndim)
990  end function get_tc_dt_ffhd
991 
992  !> anisotropic thermal conduction with slope limited symmetric scheme
993  !> Sharma 2007 Journal of Computational Physics 227, 123
994  subroutine sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,fl)
996  use mod_fix_conserve
997  integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
998  double precision, intent(in) :: x(ixi^s,1:ndim)
999  double precision, intent(in) :: w(ixi^s,1:nw)
1000  double precision, intent(inout) :: wres(ixi^s,1:nw)
1001  double precision, intent(in) :: my_dt
1002  logical, intent(in) :: fix_conserve_at_step
1003  type(tc_fluid), intent(in) :: fl
1004  !! qd store the heat conduction energy changing rate
1005  double precision :: qd(ixi^s)
1006  double precision :: rho(ixi^s),te(ixi^s)
1007  double precision :: qvec(ixi^s,1:ndim)
1008  double precision, allocatable, dimension(:^D&,:) :: qvec_equi
1009  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
1010  double precision :: alpha,dxinv(ndim)
1011  integer :: idims,idir,ix^d,ix^l,ixc^l,ixa^l,ixb^l
1012 
1013  ! coefficient of limiting on normal component
1014  if(ndim<3) then
1015  alpha=0.75d0
1016  else
1017  alpha=0.85d0
1018  end if
1019  ix^l=ixo^l^ladd1;
1020  dxinv=1.d0/dxlevel
1021  call fl%get_temperature_from_eint(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
1022  call fl%get_rho(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
1023  call set_source_tc_ffhd(ixi^l,ixo^l,w,x,fl,qvec,rho,te,alpha)
1024 
1025  qd=0.d0
1026  if(slab_uniform) then
1027  do idims=1,ndim
1028  qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
1029  ixb^l=ixo^l-kr(idims,^d);
1030  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
1031  end do
1032  else
1033  do idims=1,ndim
1034  qvec(ix^s,idims)=qvec(ix^s,idims)*block%surfaceC(ix^s,idims)
1035  ixb^l=ixo^l-kr(idims,^d);
1036  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
1037  end do
1038  qd(ixo^s)=qd(ixo^s)/block%dvolume(ixo^s)
1039  end if
1040 
1041  if(fix_conserve_at_step) then
1042  allocate(fluxall(ixi^s,1,1:ndim))
1043  fluxall(ixi^s,1,1:ndim)=my_dt*qvec(ixi^s,1:ndim)
1044  call store_flux(igrid,fluxall,1,ndim,nflux)
1045  deallocate(fluxall)
1046  end if
1047  wres(ixo^s,fl%e_)=qd(ixo^s)
1048  end subroutine sts_set_source_tc_ffhd
1049 
1050  subroutine set_source_tc_ffhd(ixI^L,ixO^L,w,x,fl,qvec,rho,Te,alpha)
1052  integer, intent(in) :: ixI^L, ixO^L
1053  double precision, intent(in) :: x(ixI^S,1:ndim)
1054  double precision, intent(in) :: w(ixI^S,1:nw)
1055  type(tc_fluid), intent(in) :: fl
1056  double precision, intent(in) :: rho(ixI^S),Te(ixI^S)
1057  double precision, intent(in) :: alpha
1058  double precision, intent(out) :: qvec(ixI^S,1:ndim)
1059  !! qd store the heat conduction energy changing rate
1060  double precision :: qd(ixI^S)
1061  double precision, dimension(ixI^S,1:ndim) :: mf,Bc,Bcf
1062  double precision, dimension(ixI^S,1:ndim) :: gradT
1063  double precision, dimension(ixI^S) :: ka,kaf,ke,kef,qdd,qe,Binv,minq,maxq,Bnorm
1064  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
1065  integer, dimension(ndim) :: lowindex
1066  integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixA^D,ixB^D
1067 
1068  ix^l=ixo^l^ladd1;
1069 
1070  ! T gradient at cell faces
1071  ! B vector
1072  mf(ixi^s,1:ndim)=block%B0(ixi^s,1:ndim,0)
1073  ! ixC is cell-corner index
1074  ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
1075  ! b unit vector at cell corner
1076  bc=0.d0
1077  {do ix^db=0,1\}
1078  ixamin^d=ixcmin^d+ix^d;
1079  ixamax^d=ixcmax^d+ix^d;
1080  bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)+mf(ixa^s,1:ndim)
1081  {end do\}
1082  bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)*0.5d0**ndim
1083  ! T gradient at cell faces
1084  do idims=1,ndim
1085  ixbmin^d=ixmin^d;
1086  ixbmax^d=ixmax^d-kr(idims,^d);
1087  call gradientc(te,ixi^l,ixb^l,idims,minq)
1088  gradt(ixb^s,idims)=minq(ixb^s)
1089  end do
1090  if(fl%tc_constant) then
1091  ka(ixc^s)=fl%tc_k_para
1092  else
1093  ! conductivity at cell center
1094  if(phys_trac) then
1095  minq(ix^s)=te(ix^s)
1096  {do ix^db=ixmin^db,ixmax^db\}
1097  if(minq(ix^d) < block%wextra(ix^d,fl%Tcoff_)) then
1098  minq(ix^d)=block%wextra(ix^d,fl%Tcoff_)
1099  end if
1100  {end do\}
1101  minq(ix^s)=fl%tc_k_para*sqrt(minq(ix^s)**5)
1102  else
1103  minq(ix^s)=fl%tc_k_para*sqrt(te(ix^s)**5)
1104  end if
1105  ka=0.d0
1106  {do ix^db=0,1\}
1107  ixbmin^d=ixcmin^d+ix^d;
1108  ixbmax^d=ixcmax^d+ix^d;
1109  ka(ixc^s)=ka(ixc^s)+minq(ixb^s)
1110  {end do\}
1111  ! cell corner conductivity
1112  ka(ixc^s)=0.5d0**ndim*ka(ixc^s)
1113  end if
1114  if(fl%tc_slope_limiter==0) then
1115  ! calculate thermal conduction flux with symmetric scheme
1116  do idims=1,ndim
1117  !qd corner values
1118  qd=0.d0
1119  {do ix^db=0,1 \}
1120  if({ ix^d==0 .and. ^d==idims | .or.}) then
1121  ixbmin^d=ixcmin^d+ix^d;
1122  ixbmax^d=ixcmax^d+ix^d;
1123  qd(ixc^s)=qd(ixc^s)+gradt(ixb^s,idims)
1124  end if
1125  {end do\}
1126  ! temperature gradient at cell corner
1127  qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
1128  end do
1129  ! b grad T at cell corner
1130  qd(ixc^s)=0.d0
1131  do idims=1,ndim
1132  qd(ixc^s)=qvec(ixc^s,idims)*bc(ixc^s,idims)+qd(ixc^s)
1133  end do
1134  do idims=1,ndim
1135  ! TC flux at cell corner
1136  gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qd(ixc^s)
1137  end do
1138  ! TC flux at cell face
1139  qvec=0.d0
1140  do idims=1,ndim
1141  ixb^l=ixo^l-kr(idims,^d);
1142  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
1143  {do ix^db=0,1 \}
1144  if({ ix^d==0 .and. ^d==idims | .or.}) then
1145  ixbmin^d=ixamin^d-ix^d;
1146  ixbmax^d=ixamax^d-ix^d;
1147  qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
1148  end if
1149  {end do\}
1150  qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
1151  if(fl%tc_saturate) then
1152  ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
1153  ! unsigned saturated TC flux = 5 phi rho c**3, c is isothermal sound speed
1154  bcf=0.d0
1155  {do ix^db=0,1 \}
1156  if({ ix^d==0 .and. ^d==idims | .or.}) then
1157  ixbmin^d=ixamin^d-ix^d;
1158  ixbmax^d=ixamax^d-ix^d;
1159  bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
1160  end if
1161  {end do\}
1162  ! averaged b at face centers
1163  bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
1164  ixb^l=ixa^l+kr(idims,^d);
1165  qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bcf(ixa^s,idims))
1166  {do ix^db=ixamin^db,ixamax^db\}
1167  if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
1168  qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
1169  end if
1170  {end do\}
1171  end if
1172  end do
1173  else
1174  ! calculate thermal conduction flux with slope-limited symmetric scheme
1175  qvec=0.d0
1176  do idims=1,ndim
1177  ixb^l=ixo^l-kr(idims,^d);
1178  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
1179  ! calculate normal of magnetic field
1180  ixb^l=ixa^l+kr(idims,^d);
1181  bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
1182  bcf=0.d0
1183  kaf=0.d0
1184  kef=0.d0
1185  {do ix^db=0,1 \}
1186  if({ ix^d==0 .and. ^d==idims | .or.}) then
1187  ixbmin^d=ixamin^d-ix^d;
1188  ixbmax^d=ixamax^d-ix^d;
1189  bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)+bc(ixb^s,1:ndim)
1190  kaf(ixa^s)=kaf(ixa^s)+ka(ixb^s)
1191  end if
1192  {end do\}
1193  ! averaged b at face centers
1194  bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)*0.5d0**(ndim-1)
1195  ! averaged thermal conductivity at face centers
1196  kaf(ixa^s)=kaf(ixa^s)*0.5d0**(ndim-1)
1197  ! limited normal component
1198  minq(ixa^s)=min(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
1199  maxq(ixa^s)=max(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
1200  ! eq (19)
1201  !corner values of gradT
1202  qdd=0.d0
1203  {do ix^db=0,1 \}
1204  if({ ix^d==0 .and. ^d==idims | .or.}) then
1205  ixbmin^d=ixcmin^d+ix^d;
1206  ixbmax^d=ixcmax^d+ix^d;
1207  qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
1208  end if
1209  {end do\}
1210  ! temperature gradient at cell corner
1211  qdd(ixc^s)=qdd(ixc^s)*0.5d0**(ndim-1)
1212  ! eq (21)
1213  qe=0.d0
1214  {do ix^db=0,1 \}
1215  qd(ixc^s)=qdd(ixc^s)
1216  if({ ix^d==0 .and. ^d==idims | .or.}) then
1217  ixbmin^d=ixamin^d-ix^d;
1218  ixbmax^d=ixamax^d-ix^d;
1219  {do ixa^db=ixamin^db,ixamax^db
1220  ixb^db=ixa^db-ix^db\}
1221  if(qd(ixb^d)<=minq(ixa^d)) then
1222  qd(ixb^d)=minq(ixa^d)
1223  else if(qd(ixb^d)>=maxq(ixa^d)) then
1224  qd(ixb^d)=maxq(ixa^d)
1225  end if
1226  {end do\}
1227  qvec(ixa^s,idims)=qvec(ixa^s,idims)+bc(ixb^s,idims)**2*qd(ixb^s)
1228  end if
1229  {end do\}
1230  qvec(ixa^s,idims)=kaf(ixa^s)*qvec(ixa^s,idims)*0.5d0**(ndim-1)
1231  ! limited transverse component, eq (17)
1232  ixbmin^d=ixamin^d;
1233  ixbmax^d=ixamax^d+kr(idims,^d);
1234  do idir=1,ndim
1235  if(idir==idims) cycle
1236  qd(ixi^s)=slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
1237  qd(ixi^s)=slope_limiter(qd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
1238  qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qd(ixa^s)
1239  end do
1240  if(fl%tc_saturate) then
1241  ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
1242  ! unsigned saturated TC flux = 5 phi rho c**3, c=sqrt(p/rho) is isothermal sound speed, phi=1.1
1243  ixb^l=ixa^l+kr(idims,^d);
1244  qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bnorm(ixa^s))
1245  {do ix^db=ixamin^db,ixamax^db\}
1246  if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
1247  !write(*,*) 'it',it,qvec(ix^D,idims),qd(ix^D),' TC saturated at ',&
1248  !x(ix^D,:),' rho',rho(ix^D),' Te',Te(ix^D)
1249  qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
1250  end if
1251  {end do\}
1252  end if
1253  end do
1254  end if
1255  end subroutine set_source_tc_ffhd
1256 end module mod_thermal_conduction
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module for flux conservation near refinement boundaries.
subroutine, public store_flux(igrid, fC, idimLIM, nwfluxin)
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
subroutine gradient(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q within ixL in direction idir.
Definition: mod_geometry.t:321
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 unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision unit_numberdensity
Physical scaling factor for number density.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
integer, dimension(:), allocatable, parameter d
double precision unit_magneticfield
Physical scaling factor for magnetic field.
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
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
double precision, dimension(ndim) dxlevel
Thermal conduction for HD and MHD or RHD and RMHD or twofl (plasma-neutral) module Adaptation of mod_...
subroutine, public sts_set_source_tc_hd(ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
subroutine, public tc_get_hd_params(fl, read_hd_params)
Init TC coefficients: HD case.
subroutine, public sts_set_source_tc_ffhd(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 set_source_tc_mhd(ixIL, ixOL, w, x, fl, qvec, rho, Te, alpha)
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 gradientc(q, ixIL, ixOL, idir, gradq)
Calculate gradient of a scalar q at cell interfaces in direction idir.
double precision function, public get_tc_dt_hd(w, ixIL, ixOL, dxD, x, fl)
Get the explicit timestep for the TC (hd implementation)
subroutine set_source_tc_hd(ixIL, ixOL, w, x, fl, qvec, rho, Te)
double precision tc_gamma_1
The adiabatic index.
subroutine, public tc_get_mhd_params(fl, read_mhd_params)
Init TC coefficients: MHD case.
subroutine set_source_tc_ffhd(ixIL, ixOL, w, x, fl, qvec, rho, Te, alpha)
double precision function, public get_tc_dt_ffhd(w, ixIL, ixOL, dxD, x, fl)
Get the explicut timestep for the TC (ffhd implementation)
double precision function, dimension(ixi^s) slope_limiter(f, ixIL, ixOL, idims, pm, tc_slope_limiter)
subroutine, public tc_get_ffhd_params(fl, read_ffhd_params)
double precision function, public get_tc_dt_mhd(w, ixIL, ixOL, dxD, x, fl)
Get the explicut timestep for the TC (mhd implementation)