MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_thermal_conduction.t
Go to the documentation of this file.
1 !> Thermal conduction for HD and MHD
2 !> Adaptation of mod_thermal_conduction for the mod_supertimestepping
3 !> In order to use it set use_mhd_tc=1 (for the mhd impl) or 2 (for the hd impl) in mhd_list (for the mhd module both hd and mhd impl can be used)
4 !> or use_new_hd_tc in hd_list parameters to true
5 !> (for the hd module, hd implementation has to be used)
6 !> The TC is set by calling one
7 !> tc_init_hd_for_total_energy and tc_init_mhd_for_total_energy might
8 !> The second argument: ixArray has to be [rho_,e_,mag(1)] for mhd (Be aware that the other components of the mag field are assumed consecutive) and [rho_,e_] for hd
9 !> additionally when internal energy equation is solved, an additional element of this array is eaux_: the index of the internal energy variable.
10 !>
11 !> 10.07.2011 developed by Chun Xia and Rony Keppens
12 !> 01.09.2012 moved to modules folder by Oliver Porth
13 !> 13.10.2013 optimized further by Chun Xia
14 !> 12.03.2014 implemented RKL2 super timestepping scheme to reduce iterations
15 !> and improve stability and accuracy up to second order in time by Chun Xia.
16 !> 23.08.2014 implemented saturation and perpendicular TC by Chun Xia
17 !> 12.01.2017 modulized by Chun Xia
18 !>
19 !> PURPOSE:
20 !> IN MHD ADD THE HEAT CONDUCTION SOURCE TO THE ENERGY EQUATION
21 !> S=DIV(KAPPA_i,j . GRAD_j T)
22 !> where KAPPA_i,j = tc_k_para b_i b_j + tc_k_perp (I - b_i b_j)
23 !> b_i b_j = B_i B_j / B**2, I is the unit matrix, and i, j= 1, 2, 3 for 3D
24 !> IN HD ADD THE HEAT CONDUCTION SOURCE TO THE ENERGY EQUATION
25 !> S=DIV(tc_k_para . GRAD T)
26 !> USAGE:
27 !> 1. in mod_usr.t -> subroutine usr_init(), add
28 !> unit_length=your length unit
29 !> unit_numberdensity=your number density unit
30 !> unit_velocity=your velocity unit
31 !> unit_temperature=your temperature unit
32 !> before call (m)hd_activate()
33 !> 2. to switch on thermal conduction in the (m)hd_list of amrvac.par add:
34 !> (m)hd_thermal_conduction=.true.
35 !> 3. in the tc_list of amrvac.par :
36 !> tc_perpendicular=.true. ! (default .false.) turn on thermal conduction perpendicular to magnetic field
37 !> tc_saturate=.true. ! (default .false. ) turn on thermal conduction saturate effect
38 !> tc_slope_limiter='MC' ! choose limiter for slope-limited anisotropic thermal conduction in MHD
39 
41  use mod_global_parameters, only: std_len
42  use mod_geometry
43  use mod_comm_lib, only: mpistop
44  implicit none
45 
46  !> The adiabatic index
47  double precision :: tc_gamma_1
48 
49  abstract interface
50  subroutine get_var_subr(w,x,ixI^L,ixO^L,res)
52  integer, intent(in) :: ixI^L, ixO^L
53  double precision, intent(in) :: w(ixI^S,nw)
54  double precision, intent(in) :: x(ixI^S,1:ndim)
55  double precision, intent(out):: res(ixI^S)
56  end subroutine get_var_subr
57 
58 
59  end interface
60 
61  type tc_fluid
62 
63  procedure(get_var_subr), pointer, nopass :: get_rho => null()
64  procedure(get_var_subr), pointer, nopass :: get_rho_equi => null()
65  procedure(get_var_subr), pointer,nopass :: get_temperature_from_eint => null()
66  procedure(get_var_subr), pointer,nopass :: get_temperature_from_conserved => null()
67  procedure(get_var_subr), pointer,nopass :: get_temperature_equi => null()
68  !> Indices of the variables
69  integer :: e_=-1
70  !> Index of cut off temperature for TRAC
71  integer :: tcoff_
72  ! if has_equi = .true. get_temperature_equi and get_rho_equi have to be set
73  logical :: has_equi=.false.
74 
75  ! the following are read from param file or set in tc_read_hd_params or tc_read_mhd_params
76  !> Coefficient of thermal conductivity (parallel to magnetic field)
77  double precision :: tc_k_para
78 
79  !> Coefficient of thermal conductivity perpendicular to magnetic field
80  double precision :: tc_k_perp
81 
82  !> Name of slope limiter for transverse component of thermal flux
83  integer :: tc_slope_limiter
84 
85  !> Logical switch for test constant conductivity
86  logical :: tc_constant=.false.
87 
88  !> Calculate thermal conduction perpendicular to magnetic field (.true.) or not (.false.)
89  logical :: tc_perpendicular=.false.
90 
91  !> Consider thermal conduction saturation effect (.true.) or not (.false.)
92  logical :: tc_saturate=.false.
93  ! END the following are read from param file or set in tc_read_hd_params or tc_read_mhd_params
94  end type tc_fluid
95 
96 
97  public :: tc_get_mhd_params
98  public :: tc_get_hd_params
99  public :: get_tc_dt_mhd
100  public :: get_tc_dt_hd
101  public :: sts_set_source_tc_mhd
102  public :: sts_set_source_tc_hd
103 
104 contains
105 
106  subroutine tc_init_params(phys_gamma)
108  double precision, intent(in) :: phys_gamma
109 
110  tc_gamma_1=phys_gamma-1d0
111  end subroutine tc_init_params
112 
113  !> Init TC coeffiecients: MHD case
114  subroutine tc_get_mhd_params(fl,read_mhd_params)
116 
117  interface
118  subroutine read_mhd_params(fl)
120  import tc_fluid
121  type(tc_fluid), intent(inout) :: fl
122 
123  end subroutine read_mhd_params
124  end interface
125 
126  type(tc_fluid), intent(inout) :: fl
127 
128  fl%tc_slope_limiter=1
129 
130  fl%tc_k_para=0.d0
131 
132  fl%tc_k_perp=0.d0
133 
134  call read_mhd_params(fl)
135 
136  if(fl%tc_k_para==0.d0 .and. fl%tc_k_perp==0.d0) then
137  if(si_unit) then
138  ! Spitzer thermal conductivity with SI units
139  fl%tc_k_para=8.d-12*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
140  ! thermal conductivity perpendicular to magnetic field
141  fl%tc_k_perp=4.d-30*unit_numberdensity**2/unit_magneticfield**2/unit_temperature**3*fl%tc_k_para
142  else
143  ! Spitzer thermal conductivity with cgs units
144  fl%tc_k_para=8.d-7*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
145  ! thermal conductivity perpendicular to magnetic field
146  fl%tc_k_perp=4.d-10*unit_numberdensity**2/unit_magneticfield**2/unit_temperature**3*fl%tc_k_para
147  end if
148  if(mype .eq. 0) print*, "Spitzer MHD: par: ",fl%tc_k_para, &
149  " ,perp: ",fl%tc_k_perp
150  else
151  fl%tc_constant=.true.
152  end if
153 
154  contains
155 
156  !> Read tc module parameters from par file: MHD case
157 
158  end subroutine tc_get_mhd_params
159 
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  !> Read tc parameters from par file: HD case
175  call read_hd_params(fl)
176  if(fl%tc_k_para==0.d0 ) then
177  if(si_unit) then
178  ! Spitzer thermal conductivity with SI units
179  fl%tc_k_para=8.d-12*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
180  else
181  ! Spitzer thermal conductivity with cgs units
182  fl%tc_k_para=8.d-7*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
183  end if
184  if(mype .eq. 0) print*, "Spitzer HD par: ",fl%tc_k_para
185  end if
186  end subroutine tc_get_hd_params
187 
188  !> Get the explicut timestep for the TC (mhd implementation)
189  function get_tc_dt_mhd(w,ixI^L,ixO^L,dx^D,x,fl) result(dtnew)
190  !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
191  !where tc_k_para_i=tc_k_para*B_i**2/B**2
192  !and T=p/rho
194 
195  type(tc_fluid), intent(in) :: fl
196  integer, intent(in) :: ixi^l, ixo^l
197  double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
198  double precision, intent(in) :: w(ixi^s,1:nw)
199  double precision :: dtnew
200 
201  double precision :: mf(ixo^s,1:ndim),te(ixi^s),b2(ixi^s),gradt(ixi^s)
202  double precision :: tmp2(ixo^s),tmp(ixo^s),hfs(ixo^s)
203  double precision :: dtdiff_tcond,maxtmp2
204  integer :: idims
205 
206  !temperature
207  call fl%get_temperature_from_conserved(w,x,ixi^l,ixi^l,te)
208 
209  ! B
210  if(b0field) then
211  mf(ixo^s,1:ndim)=w(ixo^s,iw_mag(1:ndim))+block%B0(ixo^s,1:ndim,0)
212  else
213  mf(ixo^s,1:ndim)=w(ixo^s,iw_mag(1:ndim))
214  end if
215  ! B**2
216  tmp=sum(mf**2,dim=ndim+1)
217  ! B_i**2/B**2
218  where(tmp(ixo^s)/=0.d0)
219  ^d&mf(ixo^s,^d)=mf(ixo^s,^d)**2/tmp(ixo^s);
220  elsewhere
221  ^d&mf(ixo^s,^d)=1.d0;
222  end where
223 
224  dtnew=bigdouble
225  ! B2 is now density
226  call fl%get_rho(w,x,ixi^l,ixo^l,b2)
227 
228  !tc_k_para_i
229  if(fl%tc_constant) then
230  tmp(ixo^s)=fl%tc_k_para
231  else
232  if(fl%tc_saturate) then
233  ! Kannan 2016 MN 458, 410
234  ! 3^1.5*kB^2/(4*sqrt(pi)*e^4)
235  ! l_mfpe=3.d0**1.5d0*kB_cgs**2/(4.d0*sqrt(dpi)*e_cgs**4*37.d0)=7093.9239487765044d0
236  tmp2(ixo^s)=te(ixo^s)**2/b2(ixo^s)*7093.9239487765044d0*unit_temperature**2/(unit_numberdensity*unit_length)
237  hfs=0.d0
238  do idims=1,ndim
239  call gradient(te,ixi^l,ixo^l,idims,gradt)
240  hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)*sqrt(mf(ixo^s,idims))
241  end do
242  ! kappa=kappa_Spizer/(1+4.2*l_mfpe/(T/|gradT.b|))
243  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))
244  else
245  ! kappa=kappa_Spizer
246  tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)
247  end if
248  end if
249 
250  if(slab_uniform) then
251  do idims=1,ndim
252  ! approximate thermal conduction flux: tc_k_para_i/rho/dx*B_i**2/B**2
253  tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*dxlevel(idims))
254  maxtmp2=maxval(tmp2(ixo^s))
255  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
256  dtdiff_tcond=dxlevel(idims)/(tc_gamma_1*maxtmp2+smalldouble)
257  ! limit the time step
258  dtnew=min(dtnew,dtdiff_tcond)
259  end do
260  else
261  do idims=1,ndim
262  ! approximate thermal conduction flux: tc_k_para_i/rho/dx*B_i**2/B**2
263  tmp2(ixo^s)=tmp(ixo^s)*mf(ixo^s,idims)/(b2(ixo^s)*block%ds(ixo^s,idims))
264  maxtmp2=maxval(tmp2(ixo^s)/block%ds(ixo^s,idims))
265  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
266  dtdiff_tcond=1.d0/(tc_gamma_1*maxtmp2+smalldouble)
267  ! limit the time step
268  dtnew=min(dtnew,dtdiff_tcond)
269  end do
270  end if
271  dtnew=dtnew/dble(ndim)
272 
273  end function get_tc_dt_mhd
274 
275  !> anisotropic thermal conduction with slope limited symmetric scheme
276  !> Sharma 2007 Journal of Computational Physics 227, 123
277  subroutine sts_set_source_tc_mhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,fl)
279  use mod_fix_conserve
280  integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
281  double precision, intent(in) :: x(ixi^s,1:ndim)
282  double precision, intent(in) :: w(ixi^s,1:nw)
283  double precision, intent(inout) :: wres(ixi^s,1:nw)
284  double precision, intent(in) :: my_dt
285  logical, intent(in) :: fix_conserve_at_step
286  type(tc_fluid), intent(in) :: fl
287 
288  !! qd store the heat conduction energy changing rate
289  double precision :: qd(ixi^s)
290  double precision :: rho(ixi^s),te(ixi^s)
291  double precision :: qvec(ixi^s,1:ndim)
292  double precision, allocatable, dimension(:^D&,:) :: qvec_equi
293 
294  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
295  double precision :: alpha,dxinv(ndim)
296  integer :: idims,idir,ix^d,ix^l,ixc^l,ixa^l,ixb^l
297 
298  ! coefficient of limiting on normal component
299  if(ndim<3) then
300  alpha=0.75d0
301  else
302  alpha=0.85d0
303  end if
304  ix^l=ixo^l^ladd1;
305 
306  dxinv=1.d0/dxlevel
307 
308  call fl%get_temperature_from_eint(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
309  call fl%get_rho(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
310  call set_source_tc_mhd(ixi^l,ixo^l,w,x,fl,qvec,rho,te,alpha)
311  if(fl%has_equi) then
312  allocate(qvec_equi(ixi^s,1:ndim))
313  call fl%get_temperature_equi(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
314  call fl%get_rho_equi(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
315  call set_source_tc_mhd(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te,alpha)
316  do idims=1,ndim
317  qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
318  end do
319  deallocate(qvec_equi)
320  endif
321 
322  qd=0.d0
323  if(slab_uniform) then
324  do idims=1,ndim
325  qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
326  ixb^l=ixo^l-kr(idims,^d);
327  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
328  end do
329  else
330  do idims=1,ndim
331  qvec(ix^s,idims)=qvec(ix^s,idims)*block%surfaceC(ix^s,idims)
332  ixb^l=ixo^l-kr(idims,^d);
333  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
334  end do
335  qd(ixo^s)=qd(ixo^s)/block%dvolume(ixo^s)
336  end if
337 
338  if(fix_conserve_at_step) then
339  allocate(fluxall(ixi^s,1,1:ndim))
340  fluxall(ixi^s,1,1:ndim)=my_dt*qvec(ixi^s,1:ndim)
341  call store_flux(igrid,fluxall,1,ndim,nflux)
342  deallocate(fluxall)
343  end if
344 
345  wres(ixo^s,fl%e_)=qd(ixo^s)
346 
347  end subroutine sts_set_source_tc_mhd
348 
349  subroutine set_source_tc_mhd(ixI^L,ixO^L,w,x,fl,qvec,rho,Te,alpha)
351 
352  integer, intent(in) :: ixI^L, ixO^L
353  double precision, intent(in) :: x(ixI^S,1:ndim)
354  double precision, intent(in) :: w(ixI^S,1:nw)
355  type(tc_fluid), intent(in) :: fl
356  double precision, intent(in) :: rho(ixI^S),Te(ixI^S)
357  double precision, intent(in) :: alpha
358  double precision, intent(out) :: qvec(ixI^S,1:ndim)
359 
360  !! qd store the heat conduction energy changing rate
361  double precision :: qd(ixI^S)
362 
363  double precision, dimension(ixI^S,1:ndim) :: mf,Bc,Bcf
364  double precision, dimension(ixI^S,1:ndim) :: gradT
365  double precision, dimension(ixI^S) :: ka,kaf,ke,kef,qdd,qe,Binv,minq,maxq,Bnorm
366  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
367  integer, dimension(ndim) :: lowindex
368  integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixA^D,ixB^D
369 
370  ix^l=ixo^l^ladd1;
371 
372  ! T gradient at cell faces
373  ! B vector
374  if(b0field) then
375  mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))+block%B0(ixi^s,1:ndim,0)
376  else
377  mf(ixi^s,1:ndim)=w(ixi^s,iw_mag(1:ndim))
378  end if
379  ! |B|
380  binv=dsqrt(sum(mf**2,dim=ndim+1))
381  {do ix^db=ixmin^db,ixmax^db\}
382  if(binv(ix^d)/=0.d0) then
383  binv(ix^d)=1.d0/binv(ix^d)
384  else
385  binv(ix^d)=bigdouble
386  end if
387  {end do\}
388  ! b unit vector: magnetic field direction vector
389  do idims=1,ndim
390  mf(ix^s,idims)=mf(ix^s,idims)*binv(ix^s)
391  end do
392  ! ixC is cell-corner index
393  ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
394  ! b unit vector at cell corner
395  bc=0.d0
396  {do ix^db=0,1\}
397  ixamin^d=ixcmin^d+ix^d;
398  ixamax^d=ixcmax^d+ix^d;
399  bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)+mf(ixa^s,1:ndim)
400  {end do\}
401  bc(ixc^s,1:ndim)=bc(ixc^s,1:ndim)*0.5d0**ndim
402  ! T gradient at cell faces
403  do idims=1,ndim
404  ixbmin^d=ixmin^d;
405  ixbmax^d=ixmax^d-kr(idims,^d);
406  call gradientc(te,ixi^l,ixb^l,idims,minq)
407  gradt(ixb^s,idims)=minq(ixb^s)
408  end do
409  if(fl%tc_constant) then
410  if(fl%tc_perpendicular) then
411  ka(ixc^s)=fl%tc_k_para-fl%tc_k_perp
412  ke(ixc^s)=fl%tc_k_perp
413  else
414  ka(ixc^s)=fl%tc_k_para
415  end if
416  else
417  ! conductivity at cell center
418  if(phys_trac) then
419  minq(ix^s)=te(ix^s)
420  {do ix^db=ixmin^db,ixmax^db\}
421  if(minq(ix^d) < block%wextra(ix^d,fl%Tcoff_)) then
422  minq(ix^d)=block%wextra(ix^d,fl%Tcoff_)
423  end if
424  {end do\}
425  minq(ix^s)=fl%tc_k_para*sqrt(minq(ix^s)**5)
426  else
427  minq(ix^s)=fl%tc_k_para*sqrt(te(ix^s)**5)
428  end if
429  ka=0.d0
430  {do ix^db=0,1\}
431  ixbmin^d=ixcmin^d+ix^d;
432  ixbmax^d=ixcmax^d+ix^d;
433  ka(ixc^s)=ka(ixc^s)+minq(ixb^s)
434  {end do\}
435  ! cell corner conductivity
436  ka(ixc^s)=0.5d0**ndim*ka(ixc^s)
437  ! compensate with perpendicular conductivity
438  if(fl%tc_perpendicular) then
439  minq(ix^s)=fl%tc_k_perp*rho(ix^s)**2*binv(ix^s)**2/dsqrt(te(ix^s))
440  ke=0.d0
441  {do ix^db=0,1\}
442  ixbmin^d=ixcmin^d+ix^d;
443  ixbmax^d=ixcmax^d+ix^d;
444  ke(ixc^s)=ke(ixc^s)+minq(ixb^s)
445  {end do\}
446  ! cell corner conductivity: k_parallel-k_perpendicular
447  ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
448  {do ix^db=ixcmin^db,ixcmax^db\}
449  if(ke(ix^d)<ka(ix^d)) then
450  ka(ix^d)=ka(ix^d)-ke(ix^d)
451  else
452  ke(ix^d)=ka(ix^d)
453  ka(ix^d)=0.d0
454  end if
455  {end do\}
456  end if
457  end if
458  if(fl%tc_slope_limiter==0) then
459  ! calculate thermal conduction flux with symmetric scheme
460  do idims=1,ndim
461  !qd corner values
462  qd=0.d0
463  {do ix^db=0,1 \}
464  if({ ix^d==0 .and. ^d==idims | .or.}) then
465  ixbmin^d=ixcmin^d+ix^d;
466  ixbmax^d=ixcmax^d+ix^d;
467  qd(ixc^s)=qd(ixc^s)+gradt(ixb^s,idims)
468  end if
469  {end do\}
470  ! temperature gradient at cell corner
471  qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
472  end do
473  ! b grad T at cell corner
474  qd(ixc^s)=0.d0
475  do idims=1,ndim
476  qd(ixc^s)=qvec(ixc^s,idims)*bc(ixc^s,idims)+qd(ixc^s)
477  end do
478  do idims=1,ndim
479  ! TC flux at cell corner
480  gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qd(ixc^s)
481  if(fl%tc_perpendicular) gradt(ixc^s,idims)=gradt(ixc^s,idims)+ke(ixc^s)*qvec(ixc^s,idims)
482  end do
483  ! TC flux at cell face
484  qvec=0.d0
485  do idims=1,ndim
486  ixb^l=ixo^l-kr(idims,^d);
487  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
488  {do ix^db=0,1 \}
489  if({ ix^d==0 .and. ^d==idims | .or.}) then
490  ixbmin^d=ixamin^d-ix^d;
491  ixbmax^d=ixamax^d-ix^d;
492  qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
493  end if
494  {end do\}
495  qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
496  if(fl%tc_saturate) then
497  ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
498  ! unsigned saturated TC flux = 5 phi rho c**3, c is isothermal sound speed
499  bcf=0.d0
500  {do ix^db=0,1 \}
501  if({ ix^d==0 .and. ^d==idims | .or.}) then
502  ixbmin^d=ixamin^d-ix^d;
503  ixbmax^d=ixamax^d-ix^d;
504  bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
505  end if
506  {end do\}
507  ! averaged b at face centers
508  bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
509  ixb^l=ixa^l+kr(idims,^d);
510  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))
511  {do ix^db=ixamin^db,ixamax^db\}
512  if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
513  qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
514  end if
515  {end do\}
516  end if
517  end do
518  else
519  ! calculate thermal conduction flux with slope-limited symmetric scheme
520  qvec=0.d0
521  do idims=1,ndim
522  ixb^l=ixo^l-kr(idims,^d);
523  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
524  ! calculate normal of magnetic field
525  ixb^l=ixa^l+kr(idims,^d);
526  bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
527  bcf=0.d0
528  kaf=0.d0
529  kef=0.d0
530  {do ix^db=0,1 \}
531  if({ ix^d==0 .and. ^d==idims | .or.}) then
532  ixbmin^d=ixamin^d-ix^d;
533  ixbmax^d=ixamax^d-ix^d;
534  bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)+bc(ixb^s,1:ndim)
535  kaf(ixa^s)=kaf(ixa^s)+ka(ixb^s)
536  if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)+ke(ixb^s)
537  end if
538  {end do\}
539  ! averaged b at face centers
540  bcf(ixa^s,1:ndim)=bcf(ixa^s,1:ndim)*0.5d0**(ndim-1)
541  ! averaged thermal conductivity at face centers
542  kaf(ixa^s)=kaf(ixa^s)*0.5d0**(ndim-1)
543  if(fl%tc_perpendicular) kef(ixa^s)=kef(ixa^s)*0.5d0**(ndim-1)
544  ! limited normal component
545  minq(ixa^s)=min(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
546  maxq(ixa^s)=max(alpha*gradt(ixa^s,idims),gradt(ixa^s,idims)/alpha)
547  ! eq (19)
548  !corner values of gradT
549  qdd=0.d0
550  {do ix^db=0,1 \}
551  if({ ix^d==0 .and. ^d==idims | .or.}) then
552  ixbmin^d=ixcmin^d+ix^d;
553  ixbmax^d=ixcmax^d+ix^d;
554  qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
555  end if
556  {end do\}
557  ! temperature gradient at cell corner
558  qdd(ixc^s)=qdd(ixc^s)*0.5d0**(ndim-1)
559  ! eq (21)
560  qe=0.d0
561  {do ix^db=0,1 \}
562  qd(ixc^s)=qdd(ixc^s)
563  if({ ix^d==0 .and. ^d==idims | .or.}) then
564  ixbmin^d=ixamin^d-ix^d;
565  ixbmax^d=ixamax^d-ix^d;
566  {do ixa^db=ixamin^db,ixamax^db
567  ixb^db=ixa^db-ix^db\}
568  if(qd(ixb^d)<=minq(ixa^d)) then
569  qd(ixb^d)=minq(ixa^d)
570  else if(qd(ixb^d)>=maxq(ixa^d)) then
571  qd(ixb^d)=maxq(ixa^d)
572  end if
573  {end do\}
574  qvec(ixa^s,idims)=qvec(ixa^s,idims)+bc(ixb^s,idims)**2*qd(ixb^s)
575  if(fl%tc_perpendicular) qe(ixa^s)=qe(ixa^s)+qd(ixb^s)
576  end if
577  {end do\}
578  qvec(ixa^s,idims)=kaf(ixa^s)*qvec(ixa^s,idims)*0.5d0**(ndim-1)
579  ! add normal flux from perpendicular conduction
580  if(fl%tc_perpendicular) qvec(ixa^s,idims)=qvec(ixa^s,idims)+kef(ixa^s)*qe(ixa^s)*0.5d0**(ndim-1)
581  ! limited transverse component, eq (17)
582  ixbmin^d=ixamin^d;
583  ixbmax^d=ixamax^d+kr(idims,^d);
584  do idir=1,ndim
585  if(idir==idims) cycle
586  qd(ixi^s)=slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
587  qd(ixi^s)=slope_limiter(qd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
588  qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qd(ixa^s)
589  end do
590 
591  if(fl%tc_saturate) then
592  ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
593  ! unsigned saturated TC flux = 5 phi rho c**3, c=sqrt(p/rho) is isothermal sound speed, phi=1.1
594  ixb^l=ixa^l+kr(idims,^d);
595  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))
596  {do ix^db=ixamin^db,ixamax^db\}
597  if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
598  !write(*,*) 'it',it,qvec(ix^D,idims),qd(ix^D),' TC saturated at ',&
599  !x(ix^D,:),' rho',rho(ix^D),' Te',Te(ix^D)
600  qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
601  end if
602  {end do\}
603  end if
604  end do
605  end if
606 
607  end subroutine set_source_tc_mhd
608 
609  function slope_limiter(f,ixI^L,ixO^L,idims,pm,tc_slope_limiter) result(lf)
611  integer, intent(in) :: ixi^l, ixo^l, idims, pm
612  double precision, intent(in) :: f(ixi^s)
613  double precision :: lf(ixi^s)
614  integer, intent(in) :: tc_slope_limiter
615 
616  double precision :: signf(ixi^s)
617  integer :: ixb^l
618 
619  ixb^l=ixo^l+pm*kr(idims,^d);
620  signf(ixo^s)=sign(1.d0,f(ixo^s))
621  select case(tc_slope_limiter)
622  case(1)
623  ! 'MC' montonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
624  lf(ixo^s)=two*signf(ixo^s)* &
625  max(zero,min(dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s),&
626  signf(ixo^s)*quarter*(f(ixb^s)+f(ixo^s))))
627  case(2)
628  ! 'minmod' limiter
629  lf(ixo^s)=signf(ixo^s)*max(0.d0,min(abs(f(ixo^s)),signf(ixo^s)*f(ixb^s)))
630  case(3)
631  ! 'superbee' Roes superbee limiter (eq.3.51i)
632  lf(ixo^s)=signf(ixo^s)* &
633  max(zero,min(two*dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s)),&
634  min(dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s)))
635  case(4)
636  ! 'koren' Barry Koren Right variant
637  lf(ixo^s)=signf(ixo^s)* &
638  max(zero,min(two*dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s),&
639  (two*f(ixb^s)*signf(ixo^s)+dabs(f(ixo^s)))*third))
640  case default
641  call mpistop("Unknown slope limiter for thermal conduction")
642  end select
643 
644  end function slope_limiter
645 
646  !> Calculate gradient of a scalar q at cell interfaces in direction idir
647  subroutine gradientc(q,ixI^L,ixO^L,idir,gradq)
649 
650  integer, intent(in) :: ixI^L, ixO^L, idir
651  double precision, intent(in) :: q(ixI^S)
652  double precision, intent(inout) :: gradq(ixI^S)
653  integer :: jxO^L
654 
655  associate(x=>block%x)
656 
657  jxo^l=ixo^l+kr(idir,^d);
658  select case(coordinate)
659  case(cartesian)
660  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/dxlevel(idir)
661  case(cartesian_stretched)
662  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
663  case(spherical)
664  select case(idir)
665  case(1)
666  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,1)-x(ixo^s,1))
667  {^nooned
668  case(2)
669  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/( (x(jxo^s,2)-x(ixo^s,2))*x(ixo^s,1) )
670  }
671  {^ifthreed
672  case(3)
673  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)) )
674  }
675  end select
676  case(cylindrical)
677  if(idir==phi_) then
678  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/((x(jxo^s,phi_)-x(ixo^s,phi_))*x(ixo^s,r_))
679  else
680  gradq(ixo^s)=(q(jxo^s)-q(ixo^s))/(x(jxo^s,idir)-x(ixo^s,idir))
681  end if
682  case default
683  call mpistop('Unknown geometry')
684  end select
685 
686  end associate
687  end subroutine gradientc
688 
689  function get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x,fl) result(dtnew)
690  ! Check diffusion time limit dt < dx_i**2 / ((gamma-1)*tc_k_para_i/rho)
692 
693  integer, intent(in) :: ixi^l, ixo^l
694  double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
695  double precision, intent(in) :: w(ixi^s,1:nw)
696  type(tc_fluid), intent(in) :: fl
697  double precision :: dtnew
698 
699  double precision :: tmp(ixo^s),tmp2(ixo^s),te(ixi^s),rho(ixi^s),hfs(ixo^s),gradt(ixi^s)
700  double precision :: dtdiff_tcond,maxtmp2
701  integer :: idim
702 
703  call fl%get_temperature_from_conserved(w,x,ixi^l,ixi^l,te)
704  call fl%get_rho(w,x,ixi^l,ixo^l,rho)
705 
706  if(fl%tc_saturate) then
707  ! Kannan 2016 MN 458, 410
708  ! 3^1.5*kB^2/(4*sqrt(pi)*e^4)
709  ! l_mfpe=3.d0**1.5d0*kB_cgs**2/(4.d0*sqrt(dpi)*e_cgs**4*37.d0)=7093.9239487765044d0
710  tmp2(ixo^s)=te(ixo^s)**2/rho(ixo^s)*7093.9239487765044d0*unit_temperature**2/(unit_numberdensity*unit_length)
711  hfs=0.d0
712  do idim=1,ndim
713  call gradient(te,ixi^l,ixo^l,idim,gradt)
714  hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)**2
715  end do
716  ! kappa=kappa_Spizer/(1+4.2*l_mfpe/(T/|gradT|))
717  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)))
718  else
719  tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/rho(ixo^s)
720  end if
721  dtnew = bigdouble
722 
723  if(slab_uniform) then
724  do idim=1,ndim
725  ! approximate thermal conduction flux: tc_k_para_i/rho/dx
726  tmp2(ixo^s)=tmp(ixo^s)/dxlevel(idim)
727  maxtmp2=maxval(tmp2(ixo^s))
728  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho)
729  dtdiff_tcond=dxlevel(idim)/(tc_gamma_1*maxtmp2+smalldouble)
730  ! limit the time step
731  dtnew=min(dtnew,dtdiff_tcond)
732  end do
733  else
734  do idim=1,ndim
735  ! approximate thermal conduction flux: tc_k_para_i/rho/dx
736  tmp2(ixo^s)=tmp(ixo^s)/block%ds(ixo^s,idim)
737  maxtmp2=maxval(tmp2(ixo^s)/block%ds(ixo^s,idim))
738  ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
739  dtdiff_tcond=1.d0/(tc_gamma_1*maxtmp2+smalldouble)
740  ! limit the time step
741  dtnew=min(dtnew,dtdiff_tcond)
742  end do
743  end if
744  dtnew=dtnew/dble(ndim)
745 
746  end function get_tc_dt_hd
747 
748  subroutine sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,fl)
750  use mod_fix_conserve
751 
752  integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
753  double precision, intent(in) :: x(ixi^s,1:ndim)
754  double precision, intent(in) :: w(ixi^s,1:nw)
755  double precision, intent(inout) :: wres(ixi^s,1:nw)
756  double precision, intent(in) :: my_dt
757  logical, intent(in) :: fix_conserve_at_step
758  type(tc_fluid), intent(in) :: fl
759 
760  double precision :: te(ixi^s),rho(ixi^s)
761  double precision :: qvec(ixi^s,1:ndim),qd(ixi^s)
762  double precision, allocatable, dimension(:^D&,:) :: qvec_equi
763  double precision, allocatable, dimension(:^D&,:,:) :: fluxall
764 
765  double precision :: dxinv(ndim)
766  integer :: idims,ix^l,ixb^l
767 
768  ix^l=ixo^l^ladd1;
769 
770  dxinv=1.d0/dxlevel
771 
772  !calculate Te in whole domain (+ghosts)
773  call fl%get_temperature_from_eint(w, x, ixi^l, ixi^l, te)
774  call fl%get_rho(w, x, ixi^l, ixi^l, rho)
775  call set_source_tc_hd(ixi^l,ixo^l,w,x,fl,qvec,rho,te)
776  if(fl%has_equi) then
777  allocate(qvec_equi(ixi^s,1:ndim))
778  call fl%get_temperature_equi(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
779  call fl%get_rho_equi(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
780  call set_source_tc_hd(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te)
781  do idims=1,ndim
782  qvec(ix^s,idims)=qvec(ix^s,idims) - qvec_equi(ix^s,idims)
783  end do
784  deallocate(qvec_equi)
785  endif
786 
787 
788  qd=0.d0
789  if(slab_uniform) then
790  do idims=1,ndim
791  qvec(ix^s,idims)=dxinv(idims)*qvec(ix^s,idims)
792  ixb^l=ixo^l-kr(idims,^d);
793  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
794  end do
795  else
796  do idims=1,ndim
797  qvec(ix^s,idims)=qvec(ix^s,idims)*block%surfaceC(ix^s,idims)
798  ixb^l=ixo^l-kr(idims,^d);
799  qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixb^s,idims)
800  end do
801  qd(ixo^s)=qd(ixo^s)/block%dvolume(ixo^s)
802  end if
803 
804  if(fix_conserve_at_step) then
805  allocate(fluxall(ixi^s,1,1:ndim))
806  fluxall(ixi^s,1,1:ndim)=my_dt*qvec(ixi^s,1:ndim)
807  call store_flux(igrid,fluxall,1,ndim,nflux)
808  deallocate(fluxall)
809  end if
810  wres(ixo^s,fl%e_)=qd(ixo^s)
811 
812  end subroutine sts_set_source_tc_hd
813 
814  subroutine set_source_tc_hd(ixI^L,ixO^L,w,x,fl,qvec,rho,Te)
816 
817  integer, intent(in) :: ixI^L, ixO^L
818  double precision, intent(in) :: x(ixI^S,1:ndim)
819  double precision, intent(in) :: w(ixI^S,1:nw)
820  type(tc_fluid), intent(in) :: fl
821  double precision, intent(in) :: Te(ixI^S),rho(ixI^S)
822  double precision, intent(out) :: qvec(ixI^S,1:ndim)
823 
824 
825  double precision :: gradT(ixI^S,1:ndim),ke(ixI^S),qd(ixI^S)
826 
827  integer :: idims,ix^D,ix^L,ixC^L,ixA^L,ixB^L,ixD^L
828 
829  ix^l=ixo^l^ladd1;
830  ! ixC is cell-corner index
831  ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
832 
833  !{^IFONED
834  !! cell corner temperature in ke
835  !ke=0.d0
836  !ixAmax^D=ixmax^D; ixAmin^D=ixmin^D-1;
837  !{do ix^DB=0,1\}
838  ! ixBmin^D=ixAmin^D+ix^D;
839  ! ixBmax^D=ixAmax^D+ix^D;
840  ! ke(ixA^S)=ke(ixA^S)+Te(ixB^S)
841  !{end do\}
842  !ke(ixA^S)=0.5d0**ndim*ke(ixA^S)
843  !do idims=1,ndim
844  ! ixBmin^D=ixmin^D;
845  ! ixBmax^D=ixmax^D-kr(idims,^D);
846  ! call gradient(ke,ixI^L,ixB^L,idims,qd)
847  ! gradT(ixB^S,idims)=qd(ixB^S)
848  !end do
849  !! transition region adaptive conduction
850  !if(phys_trac) then
851  ! where(ke(ixI^S) < block%wextra(ixI^S,fl%Tcoff_))
852  ! ke(ixI^S)=block%wextra(ixI^S,fl%Tcoff_)
853  ! end where
854  !end if
855  !! cell corner conduction flux
856  !do idims=1,ndim
857  ! gradT(ixC^S,idims)=gradT(ixC^S,idims)*fl%tc_k_para*sqrt(ke(ixC^S)**5)
858  !end do
859  !}
860 
861  ! calculate thermal conduction flux with symmetric scheme
862  ! T gradient (central difference) at cell corners
863  do idims=1,ndim
864  ixbmin^d=ixmin^d;
865  ixbmax^d=ixmax^d-kr(idims,^d);
866  call gradientc(te,ixi^l,ixb^l,idims,ke)
867  qd=0.d0
868  {do ix^db=0,1 \}
869  if({ix^d==0 .and. ^d==idims |.or. }) then
870  ixbmin^d=ixcmin^d+ix^d;
871  ixbmax^d=ixcmax^d+ix^d;
872  qd(ixc^s)=qd(ixc^s)+ke(ixb^s)
873  end if
874  {end do\}
875  ! temperature gradient at cell corner
876  qvec(ixc^s,idims)=qd(ixc^s)*0.5d0**(ndim-1)
877  end do
878  ! conductivity at cell center
879  if(phys_trac) then
880  ! transition region adaptive conduction
881  {do ix^db=ixmin^db,ixmax^db\}
882  if(te(ix^d) < block%wextra(ix^d,fl%Tcoff_)) then
883  qd(ix^d)=fl%tc_k_para*dsqrt(block%wextra(ix^d,fl%Tcoff_))**5
884  else
885  qd(ix^d)=fl%tc_k_para*dsqrt(te(ix^d))**5
886  end if
887  {end do\}
888  else
889  qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s))**5
890  end if
891  ke=0.d0
892  {do ix^db=0,1\}
893  ixbmin^d=ixcmin^d+ix^d;
894  ixbmax^d=ixcmax^d+ix^d;
895  ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
896  {end do\}
897  ! cell corner conductivity
898  ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
899  ! cell corner conduction flux
900  do idims=1,ndim
901  gradt(ixc^s,idims)=ke(ixc^s)*qvec(ixc^s,idims)
902  end do
903 
904  if(fl%tc_saturate) then
905  ! consider saturation with unsigned saturated TC flux = 5 phi rho c**3
906  ! saturation flux at cell center
907  qd(ix^s)=5.5d0*rho(ix^s)*dsqrt(te(ix^s)**3)
908  !cell corner values of qd in ke
909  ke=0.d0
910  {do ix^db=0,1\}
911  ixbmin^d=ixcmin^d+ix^d;
912  ixbmax^d=ixcmax^d+ix^d;
913  ke(ixc^s)=ke(ixc^s)+qd(ixb^s)
914  {end do\}
915  ! cell corner saturation flux
916  ke(ixc^s)=0.5d0**ndim*ke(ixc^s)
917  ! magnitude of cell corner conduction flux
918  qd(ixc^s)=norm2(gradt(ixc^s,:),dim=ndim+1)
919  {do ix^db=ixcmin^db,ixcmax^db\}
920  if(qd(ix^d)>ke(ix^d)) then
921  ke(ix^d)=ke(ix^d)/(qd(ix^d)+smalldouble)
922  do idims=1,ndim
923  gradt(ix^d,idims)=ke(ix^d)*gradt(ix^d,idims)
924  end do
925  end if
926  {end do\}
927  end if
928 
929  ! conduction flux at cell face
930  qvec=0.d0
931  do idims=1,ndim
932  ixb^l=ixo^l-kr(idims,^d);
933  ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
934  {do ix^db=0,1 \}
935  if({ ix^d==0 .and. ^d==idims | .or.}) then
936  ixbmin^d=ixamin^d-ix^d;
937  ixbmax^d=ixamax^d-ix^d;
938  qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
939  end if
940  {end do\}
941  qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
942  end do
943 
944  end subroutine set_source_tc_hd
945 
946 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 Adaptation of mod_thermal_conduction for the mod_supertimestepping ...
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)
Read tc module parameters from par file: MHD case.
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)
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 coeffiecients: MHD case.
double precision function, dimension(ixi^s) slope_limiter(f, ixIL, ixOL, idims, pm, tc_slope_limiter)
double precision function, public get_tc_dt_mhd(w, ixIL, ixOL, dxD, x, fl)
Get the explicut timestep for the TC (mhd implementation)