MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
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!> 30.11.2025 Minor cleanup (for consistency between hd and mhd)
22!>
23!> PURPOSE:
24!> IN MHD ADD THE HEAT CONDUCTION SOURCE TO THE ENERGY EQUATION
25!> S=DIV(KAPPA_i,j . GRAD_j T)
26!> where KAPPA_i,j = tc_k_para b_i b_j + tc_k_perp (I - b_i b_j)
27!> b_i b_j = B_i B_j / B**2, I is the unit matrix, and i, j= 1, 2, 3 for 3D
28!> IN HD ADD THE HEAT CONDUCTION SOURCE TO THE ENERGY EQUATION
29!> S=DIV(tc_k_para . GRAD T)
30!> USAGE:
31!> 1. in mod_usr.t -> subroutine usr_init(), add
32!> unit_length=your length unit
33!> unit_numberdensity=your number density unit
34!> unit_velocity=your velocity unit
35!> unit_temperature=your temperature unit
36!> before call (m)hd_activate()
37!> 2. to switch on thermal conduction in the (r)(m)hd_list of amrvac.par add:
38!> (r)(m)hd_thermal_conduction=.true.
39!> 3. in the tc_list of amrvac.par :
40!> tc_perpendicular=.true. ! (default .false.) turn on thermal conduction perpendicular to magnetic field
41!> tc_saturate=.true. ! (default .false. ) turn on thermal conduction saturate effect
42!> tc_slope_limiter='MC' ! choose limiter for slope-limited anisotropic thermal conduction in MHD
43!> note: twofl_list incorporates instances for charges and neutrals
44
46 use mod_global_parameters, only: std_len
47 use mod_geometry
48 use mod_comm_lib, only: mpistop
49 implicit none
50
51 !> The adiabatic index
52 double precision :: tc_gamma_1
53
54 abstract interface
55 subroutine get_var_subr(w,x,ixI^L,ixO^L,res)
57 integer, intent(in) :: ixI^L, ixO^L
58 double precision, intent(in) :: w(ixI^S,nw)
59 double precision, intent(in) :: x(ixI^S,1:ndim)
60 double precision, intent(out):: res(ixI^S)
61 end subroutine get_var_subr
62 end interface
63
65
66 ! BEGIN the following are read from param file or set in tc_read_hd_params or tc_read_mhd_params
67 !> Coefficient of thermal conductivity (parallel to magnetic field)
68 double precision :: tc_k_para
69
70 !> Coefficient of thermal conductivity perpendicular to magnetic field
71 double precision :: tc_k_perp
72
73 !> Indices of the variables
74 integer :: e_=-1
75 !> Index of cut off temperature for TRAC
76 integer :: tcoff_
77 !> Name of slope limiter for transverse component of thermal flux
78 integer :: tc_slope_limiter
79
80 ! if has_equi = .true. get_temperature_equi and get_rho_equi have to be set
81 logical :: has_equi=.false.
82
83 !> Logical switch for test constant conductivity
84 logical :: tc_constant=.false.
85
86 !> Calculate thermal conduction perpendicular to magnetic field (.true.) or not (.false.)
87 logical :: tc_perpendicular=.false.
88
89 !> Consider thermal conduction saturation effect (.true.) or not (.false.)
90 logical :: tc_saturate=.false.
91 ! END the following are read from param file or set in tc_read_hd_params or tc_read_mhd_params
92 procedure(get_var_subr), pointer, nopass :: get_rho => null()
93 procedure(get_var_subr), pointer, nopass :: get_rho_equi => null()
94 procedure(get_var_subr), pointer,nopass :: get_temperature_from_eint => null()
95 procedure(get_var_subr), pointer,nopass :: get_temperature_from_conserved => null()
96 procedure(get_var_subr), pointer,nopass :: get_temperature_equi => null()
97 end type tc_fluid
98
99 public :: tc_get_mhd_params
100 public :: tc_get_hd_params
101 public :: get_tc_dt_mhd
102 public :: get_tc_dt_hd
103 public :: sts_set_source_tc_mhd
104 public :: sts_set_source_tc_hd
105
106contains
107
108 subroutine tc_init_params(phys_gamma)
110 double precision, intent(in) :: phys_gamma
111
112 tc_gamma_1=phys_gamma-1d0
113 end subroutine tc_init_params
114
115 !> Init TC coefficients: MHD case
116 subroutine tc_get_mhd_params(fl,read_mhd_params)
118
119 interface
120 subroutine read_mhd_params(fl)
122 import tc_fluid
123 type(tc_fluid), intent(inout) :: fl
124
125 end subroutine read_mhd_params
126 end interface
127 type(tc_fluid), intent(inout) :: fl
128
129 fl%tc_slope_limiter=1
130 fl%tc_k_para=0.d0
131 fl%tc_k_perp=0.d0
132
133 !> Read tc module parameters from par file: MHD case
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 if(mype .eq. 0) print*, "Constant thermal conduction coefficients with values: ",fl%tc_k_para,fl%tc_k_perp
153 end if
154 end subroutine tc_get_mhd_params
155
156 !> Init TC coefficients: HD case
157 subroutine tc_get_hd_params(fl,read_hd_params)
159
160 interface
161 subroutine read_hd_params(fl)
163 import tc_fluid
164 type(tc_fluid), intent(inout) :: fl
165
166 end subroutine read_hd_params
167 end interface
168 type(tc_fluid), intent(inout) :: fl
169
170 fl%tc_k_para=0.d0
171
172 !> Read tc parameters from par file: HD case
173 call read_hd_params(fl)
174
175 if(fl%tc_k_para==0.d0) then
176 if(si_unit) then
177 ! Spitzer thermal conductivity with SI units
178 fl%tc_k_para=8.d-12*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
179 else
180 ! Spitzer thermal conductivity with cgs units
181 fl%tc_k_para=8.d-7*unit_temperature**3.5d0/unit_length/unit_density/unit_velocity**3
182 end if
183 if(mype .eq. 0) print*, "Spitzer HD par: ",fl%tc_k_para
184 else
185 fl%tc_constant=.true.
186 if(mype .eq. 0) print*, "Constant thermal conduction coefficient at value: ",fl%tc_k_para
187 end if
188
189 end subroutine tc_get_hd_params
190
191 !> Get the explicit timestep for the TC (mhd implementation)
192 !> Note: for multi-D MHD (1D MHD will use HD fall-back)
193 function get_tc_dt_mhd(w,ixI^L,ixO^L,dx^D,x,fl) result(dtnew)
194 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
195 !where tc_k_para_i=tc_k_para*B_i**2/B**2
196 !and T=p/rho
198
199 type(tc_fluid), intent(in) :: fl
200 integer, intent(in) :: ixi^l, ixo^l
201 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
202 double precision, intent(in) :: w(ixi^s,1:nw)
203 double precision :: dtnew
204
205 double precision :: mf(ixo^s,1:ndim),te(ixi^s),rho(ixi^s),gradt(ixi^s)
206 double precision :: tmp(ixo^s),hfs(ixo^s),blocal(1:ndim)
207 double precision :: dtdiff_tcond,maxtmp2
208 integer :: idims,ix^d
209
210
211 ! B
212 if(allocated(iw_mag)) then
213 if(b0field) then
214 {do ix^db=ixomin^db,ixomax^db\}
215 ^d&blocal(^d)=w({ix^d},iw_mag(^d))+block%B0({ix^d},^d,0)\
216 {^iftwod
217 if(blocal(1)/=0.d0) then
218 mf(ix^d,1)=sign(1.d0,blocal(1))/dsqrt(1.d0+(blocal(2)/blocal(1))**2)
219 else
220 mf(ix^d,1)=0.d0
221 end if
222 if(blocal(2)/=0.d0) then
223 mf(ix^d,2)=sign(1.d0,blocal(2))/dsqrt(1.d0+(blocal(1)/blocal(2))**2)
224 else
225 mf(ix^d,2)=0.d0
226 end if
227 }
228 {^ifthreed
229 if(blocal(1)/=0.d0) then
230 mf(ix^d,1)=sign(1.d0,blocal(1))/dsqrt(1.d0+(blocal(2)/blocal(1))**2+(blocal(3)/blocal(1))**2)
231 else
232 mf(ix^d,1)=0.d0
233 end if
234 if(blocal(2)/=0.d0) then
235 mf(ix^d,2)=sign(1.d0,blocal(2))/dsqrt(1.d0+(blocal(1)/blocal(2))**2+(blocal(3)/blocal(2))**2)
236 else
237 mf(ix^d,2)=0.d0
238 end if
239 if(blocal(3)/=0.d0) then
240 mf(ix^d,3)=sign(1.d0,blocal(3))/dsqrt(1.d0+(blocal(1)/blocal(3))**2+(blocal(2)/blocal(3))**2)
241 else
242 mf(ix^d,3)=0.d0
243 end if
244 }
245 {end do\}
246 else
247 {do ix^db=ixomin^db,ixomax^db\}
248 {^iftwod
249 if(w(ix^d,iw_mag(1))/=0.d0) then
250 mf(ix^d,1)=sign(1.d0,w(ix^d,iw_mag(1)))/dsqrt(1.d0+(w(ix^d,iw_mag(2))/w(ix^d,iw_mag(1)))**2)
251 else
252 mf(ix^d,1)=0.d0
253 end if
254 if(w(ix^d,iw_mag(2))/=0.d0) then
255 mf(ix^d,2)=sign(1.d0,w(ix^d,iw_mag(2)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(2)))**2)
256 else
257 mf(ix^d,2)=0.d0
258 end if
259 }
260 {^ifthreed
261 if(w(ix^d,iw_mag(1))/=0.d0) then
262 mf(ix^d,1)=sign(1.d0,w(ix^d,iw_mag(1)))/dsqrt(1.d0+(w(ix^d,iw_mag(2))/w(ix^d,iw_mag(1)))**2+&
263 (w(ix^d,iw_mag(3))/w(ix^d,iw_mag(1)))**2)
264 else
265 mf(ix^d,1)=0.d0
266 end if
267 if(w(ix^d,iw_mag(2))/=0.d0) then
268 mf(ix^d,2)=sign(1.d0,w(ix^d,iw_mag(2)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(2)))**2+&
269 (w(ix^d,iw_mag(3))/w(ix^d,iw_mag(2)))**2)
270 else
271 mf(ix^d,2)=0.d0
272 end if
273 if(w(ix^d,iw_mag(3))/=0.d0) then
274 mf(ix^d,3)=sign(1.d0,w(ix^d,iw_mag(3)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(3)))**2+&
275 (w(ix^d,iw_mag(2))/w(ix^d,iw_mag(3)))**2)
276 else
277 mf(ix^d,3)=0.d0
278 end if
279 }
280 {end do\}
281 end if
282 else
283 ! this if for ffhd or other physics modules without B components
284 mf(ixo^s,1:ndim)=block%B0(ixo^s,1:ndim,0)
285 end if
286
287 !temperature
288 call fl%get_temperature_from_conserved(w,x,ixi^l,ixi^l,te)
289 call fl%get_rho(w,x,ixi^l,ixo^l,rho)
290
291 !tc_k_para_i
292 if(fl%tc_constant) then
293 tmp(ixo^s)=fl%tc_k_para
294 else
295 if(fl%tc_saturate) then
296 ! Kannan 2016 MN 458, 410
297 ! 3^1.5*kB^2/(4*sqrt(pi)*e^4)
298 ! l_mfpe=3.d0**1.5d0*kB_cgs**2/(4.d0*sqrt(dpi)*e_cgs**4*37.d0)=7093.9239487765044d0
299 tmp(ixo^s)=te(ixo^s)**2/rho(ixo^s)*7093.9239487765044d0*unit_temperature**2/(unit_numberdensity*unit_length)
300 do idims=1,ndim
301 call gradient(te,ixi^l,ixo^l,idims,gradt)
302 if(idims==1) then
303 hfs(ixo^s)=gradt(ixo^s)*mf(ixo^s,idims)
304 else
305 hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)*mf(ixo^s,idims)
306 end if
307 end do
308 ! kappa=kappa_Spitzer/(1+4.2*l_mfpe/(T/|gradT.b|))
309 tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)/(1.d0+4.2d0*tmp(ixo^s)*dabs(hfs(ixo^s))/te(ixo^s))
310 else
311 ! kappa=kappa_Spitzer
312 tmp(ixo^s)=fl%tc_k_para*dsqrt(te(ixo^s)**5)
313 end if
314 end if
315
316 dtnew=bigdouble
317 do idims=1,ndim
318 ! approximate thermal conduction flux: tc_k_para_i/rho/dx*B_i**2/B**2
319 maxtmp2=maxval(tmp(ixo^s)*mf(ixo^s,idims)**2/(rho(ixo^s)*block%ds(ixo^s,idims)**2))
320 ! dt< dx_idim**2/((gamma-1)*tc_k_para_i/rho*B_i**2/B**2)
321 dtdiff_tcond=1.d0/(tc_gamma_1*maxtmp2+smalldouble)
322 ! limit the time step
323 dtnew=min(dtnew,dtdiff_tcond)
324 end do
325 dtnew=dtnew/dble(ndim)
326 end function get_tc_dt_mhd
327
328 !> anisotropic thermal conduction with slope limited symmetric scheme
329 !> Sharma 2007 Journal of Computational Physics 227, 123
330 subroutine sts_set_source_tc_mhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,fl)
333 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
334 double precision, intent(in) :: x(ixi^s,1:ndim)
335 double precision, intent(in) :: w(ixi^s,1:nw)
336 double precision, intent(inout) :: wres(ixi^s,1:nw)
337 double precision, intent(in) :: my_dt
338 logical, intent(in) :: fix_conserve_at_step
339 type(tc_fluid), intent(in) :: fl
340
341 !! qd store the heat conduction energy changing rate
342 double precision :: qd(ixo^s)
343 double precision :: rho(ixi^s),te(ixi^s)
344 double precision :: qvec(ixi^s,1:ndim)
345 double precision :: fluxall(ixi^s,1,1:ndim)
346 double precision :: alpha,dxinv(ndim)
347 double precision, allocatable, dimension(:^D&,:) :: qvec_equi
348 integer :: idims,ixa^l
349
350 ! coefficient of limiting on normal component
351 if(ndim<3) then
352 alpha=0.75d0
353 else
354 alpha=0.85d0
355 end if
356
357 dxinv=1.d0/dxlevel
358
359 call fl%get_temperature_from_eint(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
360 call fl%get_rho(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
361 if(slab_uniform) then
362 call set_source_tc_mhd(ixi^l,ixo^l,w,x,fl,qvec,rho,te,alpha)
363 else
364 call set_source_tc_mhd_geo(ixi^l,ixo^l,w,x,fl,qvec,rho,te,alpha)
365 end if
366 if(fl%has_equi) then
367 allocate(qvec_equi(ixi^s,1:ndim))
368 call fl%get_temperature_equi(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
369 call fl%get_rho_equi(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
370 if(slab_uniform) then
371 call set_source_tc_mhd(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te,alpha)
372 else
373 call set_source_tc_mhd_geo(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te,alpha)
374 end if
375 do idims=1,ndim
376 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
377 qvec(ixa^s,idims)=qvec(ixa^s,idims)-qvec_equi(ixa^s,idims)
378 end do
379 deallocate(qvec_equi)
380 end if
381
382 if(slab_uniform) then
383 do idims=1,ndim
384 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
385 qvec(ixa^s,idims)=dxinv(idims)*qvec(ixa^s,idims)
386 ixa^l=ixo^l-kr(idims,^d);
387 if(idims==1) then
388 qd(ixo^s)=qvec(ixo^s,idims)-qvec(ixa^s,idims)
389 else
390 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixa^s,idims)
391 end if
392 end do
393 else
394 do idims=1,ndim
395 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
396 qvec(ixa^s,idims)=qvec(ixa^s,idims)*block%surfaceC(ixa^s,idims)
397 ixa^l=ixo^l-kr(idims,^d);
398 if(idims==1) then
399 qd(ixo^s)=qvec(ixo^s,idims)-qvec(ixa^s,idims)
400 else
401 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixa^s,idims)
402 end if
403 end do
404 qd(ixo^s)=qd(ixo^s)/block%dvolume(ixo^s)
405 end if
406
407 if(fix_conserve_at_step) then
408 do idims=1,ndim
409 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
410 fluxall(ixa^s,1,idims)=my_dt*qvec(ixa^s,idims)
411 end do
412 call store_flux(igrid,fluxall,1,ndim,nflux)
413 end if
414
415 wres(ixo^s,fl%e_)=qd(ixo^s)
416 end subroutine sts_set_source_tc_mhd
417
418 subroutine set_source_tc_mhd(ixI^L,ixO^L,w,x,fl,qvec,rho,Te,alpha)
420 integer, intent(in) :: ixI^L, ixO^L
421 double precision, intent(in) :: x(ixI^S,1:ndim)
422 double precision, intent(in) :: w(ixI^S,1:nw)
423 type(tc_fluid), intent(in) :: fl
424 double precision, intent(in) :: rho(ixI^S),Te(ixI^S)
425 double precision, intent(in) :: alpha
426 double precision, intent(out) :: qvec(ixI^S,1:ndim)
427
428 !! qdd store the heat conduction energy changing rate
429 double precision, dimension(ixI^S,1:ndim) :: mf,Bc,Bcf,gradT
430 double precision, dimension(ixI^S) :: ka,kaf,ke,kef,qdd,Bnorm
431 double precision :: minq,maxq,qd(ixI^S,2**(ndim-1)), blocal(ndim)
432 integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L
433
434 ix^l=ixo^l^ladd1;
435
436 ! T gradient at cell faces
437 ! b unit vector mf: magnetic field direction vector
438 if(allocated(iw_mag)) then
439 if(b0field) then
440 {do ix^db=ixmin^db,ixmax^db\}
441 ^d&blocal(^d)=w({ix^d},iw_mag(^d))+block%B0({ix^d},^d,0)\
442 {^iftwod
443 if(blocal(1)/=0.d0) then
444 mf(ix^d,1)=sign(1.d0,blocal(1))/dsqrt(1.d0+(blocal(2)/blocal(1))**2)
445 else
446 mf(ix^d,1)=0.d0
447 end if
448 if(blocal(2)/=0.d0) then
449 mf(ix^d,2)=sign(1.d0,blocal(2))/dsqrt(1.d0+(blocal(1)/blocal(2))**2)
450 else
451 mf(ix^d,2)=0.d0
452 end if
453 }
454 {^ifthreed
455 if(blocal(1)/=0.d0) then
456 mf(ix^d,1)=sign(1.d0,blocal(1))/dsqrt(1.d0+(blocal(2)/blocal(1))**2+(blocal(3)/blocal(1))**2)
457 else
458 mf(ix^d,1)=0.d0
459 end if
460 if(blocal(2)/=0.d0) then
461 mf(ix^d,2)=sign(1.d0,blocal(2))/dsqrt(1.d0+(blocal(1)/blocal(2))**2+(blocal(3)/blocal(2))**2)
462 else
463 mf(ix^d,2)=0.d0
464 end if
465 if(blocal(3)/=0.d0) then
466 mf(ix^d,3)=sign(1.d0,blocal(3))/dsqrt(1.d0+(blocal(1)/blocal(3))**2+(blocal(2)/blocal(3))**2)
467 else
468 mf(ix^d,3)=0.d0
469 end if
470 }
471 {end do\}
472 else
473 {do ix^db=ixmin^db,ixmax^db\}
474 {^iftwod
475 if(w(ix^d,iw_mag(1))/=0.d0) then
476 mf(ix^d,1)=sign(1.d0,w(ix^d,iw_mag(1)))/dsqrt(1.d0+(w(ix^d,iw_mag(2))/w(ix^d,iw_mag(1)))**2)
477 else
478 mf(ix^d,1)=0.d0
479 end if
480 if(w(ix^d,iw_mag(2))/=0.d0) then
481 mf(ix^d,2)=sign(1.d0,w(ix^d,iw_mag(2)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(2)))**2)
482 else
483 mf(ix^d,2)=0.d0
484 end if
485 }
486 {^ifthreed
487 if(w(ix^d,iw_mag(1))/=0.d0) then
488 mf(ix^d,1)=sign(1.d0,w(ix^d,iw_mag(1)))/dsqrt(1.d0+(w(ix^d,iw_mag(2))/w(ix^d,iw_mag(1)))**2+&
489 (w(ix^d,iw_mag(3))/w(ix^d,iw_mag(1)))**2)
490 else
491 mf(ix^d,1)=0.d0
492 end if
493 if(w(ix^d,iw_mag(2))/=0.d0) then
494 mf(ix^d,2)=sign(1.d0,w(ix^d,iw_mag(2)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(2)))**2+&
495 (w(ix^d,iw_mag(3))/w(ix^d,iw_mag(2)))**2)
496 else
497 mf(ix^d,2)=0.d0
498 end if
499 if(w(ix^d,iw_mag(3))/=0.d0) then
500 mf(ix^d,3)=sign(1.d0,w(ix^d,iw_mag(3)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(3)))**2+&
501 (w(ix^d,iw_mag(2))/w(ix^d,iw_mag(3)))**2)
502 else
503 mf(ix^d,3)=0.d0
504 end if
505 }
506 {end do\}
507 end if
508 else
509 mf(ix^s,1:ndim)=block%B0(ix^s,1:ndim,0)
510 endif
511 ! ixC is cell-corner index
512 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
513 ! b unit vector at cell corner
514 {^ifthreed
515 do idims=1,3
516 {do ix^db=ixcmin^db,ixcmax^db\}
517 bc(ix^d,idims)=0.125d0*(mf(ix1,ix2,ix3,idims)+mf(ix1+1,ix2,ix3,idims)&
518 +mf(ix1,ix2+1,ix3,idims)+mf(ix1+1,ix2+1,ix3,idims)&
519 +mf(ix1,ix2,ix3+1,idims)+mf(ix1+1,ix2,ix3+1,idims)&
520 +mf(ix1,ix2+1,ix3+1,idims)+mf(ix1+1,ix2+1,ix3+1,idims))
521 {end do\}
522 end do
523 }
524 {^iftwod
525 do idims=1,2
526 {do ix^db=ixcmin^db,ixcmax^db\}
527 bc(ix^d,idims)=0.25d0*(mf(ix1,ix2,idims)+mf(ix1+1,ix2,idims)&
528 +mf(ix1,ix2+1,idims)+mf(ix1+1,ix2+1,idims))
529 {end do\}
530 end do
531 }
532 ! T gradient at cell faces
533 do idims=1,ndim
534 ixbmin^d=ixmin^d;
535 ixbmax^d=ixmax^d-kr(idims,^d);
536 call gradientf(te,x,ixi^l,ixb^l,idims,gradt(ixi^s,idims))
537 end do
538 if(fl%tc_constant) then
539 if(fl%tc_perpendicular) then
540 ka(ixc^s)=fl%tc_k_para-fl%tc_k_perp
541 ke(ixc^s)=fl%tc_k_perp
542 else
543 ka(ixc^s)=fl%tc_k_para
544 end if
545 else
546 ! conductivity at cell center
547 if(phys_trac) then
548 {do ix^db=ixmin^db,ixmax^db\}
549 if(te(ix^d) < block%wextra(ix^d,fl%Tcoff_)) then
550 qdd(ix^d)=fl%tc_k_para*dsqrt(block%wextra(ix^d,fl%Tcoff_)**5)
551 else
552 qdd(ix^d)=fl%tc_k_para*dsqrt(te(ix^d)**5)
553 end if
554 {end do\}
555 else
556 qdd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s)**5)
557 end if
558 ! cell corner parallel conductivity in ka
559 {^ifthreed
560 {do ix^db=ixcmin^db,ixcmax^db\}
561 ka(ix^d)=0.125d0*(qdd(ix1,ix2,ix3)+qdd(ix1+1,ix2,ix3)&
562 +qdd(ix1,ix2+1,ix3)+qdd(ix1+1,ix2+1,ix3)&
563 +qdd(ix1,ix2,ix3+1)+qdd(ix1+1,ix2,ix3+1)&
564 +qdd(ix1,ix2+1,ix3+1)+qdd(ix1+1,ix2+1,ix3+1))
565 {end do\}
566 }
567 {^iftwod
568 {do ix^db=ixcmin^db,ixcmax^db\}
569 ka(ix^d)=0.25d0*(qdd(ix1,ix2)+qdd(ix1+1,ix2)&
570 +qdd(ix1,ix2+1)+qdd(ix1+1,ix2+1))
571 {end do\}
572 }
573 ! compensate with perpendicular conductivity
574 if(fl%tc_perpendicular) then
575 if(b0field) then
576 qdd(ix^s)=fl%tc_k_perp*rho(ix^s)**2/((^d&(w(ix^s,iw_mag(^d))+block%B0(ix^s,^d,0))**2+)*dsqrt(te(ix^s))+smalldouble)
577 else
578 qdd(ix^s)=fl%tc_k_perp*rho(ix^s)**2/((^d&w(ix^s,iw_mag(^d))**2+)*dsqrt(te(ix^s))+smalldouble)
579 end if
580 {^ifthreed
581 {do ix^db=ixcmin^db,ixcmax^db\}
582 ke(ix^d)=0.125d0*(qdd(ix1,ix2,ix3)+qdd(ix1+1,ix2,ix3)&
583 +qdd(ix1,ix2+1,ix3)+qdd(ix1+1,ix2+1,ix3)&
584 +qdd(ix1,ix2,ix3+1)+qdd(ix1+1,ix2,ix3+1)&
585 +qdd(ix1,ix2+1,ix3+1)+qdd(ix1+1,ix2+1,ix3+1))
586 if(ke(ix^d)<ka(ix^d)) then
587 ka(ix^d)=ka(ix^d)-ke(ix^d)
588 else
589 ke(ix^d)=ka(ix^d)
590 ka(ix^d)=0.d0
591 end if
592 {end do\}
593 }
594 {^iftwod
595 {do ix^db=ixcmin^db,ixcmax^db\}
596 ke(ix^d)=0.25d0*(qdd(ix1,ix2)+qdd(ix1+1,ix2)&
597 +qdd(ix1,ix2+1)+qdd(ix1+1,ix2+1))
598 if(ke(ix^d)<ka(ix^d)) then
599 ka(ix^d)=ka(ix^d)-ke(ix^d)
600 else
601 ke(ix^d)=ka(ix^d)
602 ka(ix^d)=0.d0
603 end if
604 {end do\}
605 }
606 end if
607 end if
608 if(fl%tc_slope_limiter==0) then
609 ! calculate thermal conduction flux with symmetric scheme
610 do idims=1,ndim
611 !qdd corner values
612 qdd=0.d0
613 {do ix^db=0,1 \}
614 if({ ix^d==0 .and. ^d==idims | .or.}) then
615 ixbmin^d=ixcmin^d+ix^d;
616 ixbmax^d=ixcmax^d+ix^d;
617 qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
618 end if
619 {end do\}
620 ! temperature gradient at cell corner
621 qvec(ixc^s,idims)=qdd(ixc^s)*0.5d0**(ndim-1)
622 end do
623 ! b grad T at cell corner
624 qdd(ixc^s)=sum(qvec(ixc^s,1:ndim)*bc(ixc^s,1:ndim),dim=ndim+1)
625 do idims=1,ndim
626 ! TC flux at cell corner
627 gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qdd(ixc^s)
628 if(fl%tc_perpendicular) gradt(ixc^s,idims)=gradt(ixc^s,idims)+ke(ixc^s)*qvec(ixc^s,idims)
629 end do
630 ! TC flux at cell face
631 qvec=0.d0
632 do idims=1,ndim
633 ixb^l=ixo^l-kr(idims,^d);
634 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
635 {do ix^db=0,1 \}
636 if({ ix^d==0 .and. ^d==idims | .or.}) then
637 ixbmin^d=ixamin^d-ix^d;
638 ixbmax^d=ixamax^d-ix^d;
639 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
640 end if
641 {end do\}
642 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
643 if(fl%tc_saturate) then
644 ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
645 ! unsigned saturated TC flux = 5 phi rho c**3, c is isothermal sound speed
646 bcf=0.d0
647 {do ix^db=0,1 \}
648 if({ ix^d==0 .and. ^d==idims | .or.}) then
649 ixbmin^d=ixamin^d-ix^d;
650 ixbmax^d=ixamax^d-ix^d;
651 bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
652 end if
653 {end do\}
654 ! averaged b at face centers
655 bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
656 ixb^l=ixa^l+kr(idims,^d);
657 qdd(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))
658 {do ix^db=ixamin^db,ixamax^db\}
659 if(dabs(qvec(ix^d,idims))>qdd(ix^d)) then
660 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qdd(ix^d)
661 end if
662 {end do\}
663 end if
664 end do
665 else
666 ! calculate thermal conduction flux with slope-limited symmetric scheme
667 do idims=1,ndim
668 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
669 {^ifthreed
670 if(idims==1) then
671 {do ix^db=ixamin^db,ixamax^db\}
672 ! averaged b at face centers
673 ^d&bcf({ix^d},^d)=0.25d0*(bc({ix^d},^d)+bc(ix1,ix2-1,ix3,^d)&
674 +bc(ix1,ix2,ix3-1,^d)+bc(ix1,ix2-1,ix3-1,^d))\
675 kaf(ix^d)=0.25d0*(ka(ix1,ix2,ix3)+ka(ix1,ix2-1,ix3)&
676 +ka(ix1,ix2,ix3-1)+ka(ix1,ix2-1,ix3-1))
677 ! averaged thermal conductivity at face centers
678 if(fl%tc_perpendicular) &
679 kef(ix^d)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1,ix2-1,ix3)&
680 +ke(ix1,ix2,ix3-1)+ke(ix1,ix2-1,ix3-1))
681 {end do\}
682 else if(idims==2) then
683 {do ix^db=ixamin^db,ixamax^db\}
684 ^d&bcf({ix^d},^d)=0.25d0*(bc({ix^d},^d)+bc(ix1-1,ix2,ix3,^d)&
685 +bc(ix1,ix2,ix3-1,^d)+bc(ix1-1,ix2,ix3-1,^d))\
686 kaf(ix^d)=0.25d0*(ka(ix1,ix2,ix3)+ka(ix1-1,ix2,ix3)&
687 +ka(ix1,ix2,ix3-1)+ka(ix1-1,ix2,ix3-1))
688 if(fl%tc_perpendicular) &
689 kef(ix^d)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1-1,ix2,ix3)&
690 +ke(ix1,ix2,ix3-1)+ke(ix1-1,ix2,ix3-1))
691 {end do\}
692 else
693 {do ix^db=ixamin^db,ixamax^db\}
694 ^d&bcf({ix^d},^d)=0.25d0*(bc({ix^d},^d)+bc(ix1,ix2-1,ix3,^d)&
695 +bc(ix1-1,ix2,ix3,^d)+bc(ix1-1,ix2-1,ix3,^d))\
696 kaf(ix^d)=0.25d0*(ka(ix1,ix2,ix3)+ka(ix1,ix2-1,ix3)&
697 +ka(ix1-1,ix2,ix3)+ka(ix1-1,ix2-1,ix3))
698 if(fl%tc_perpendicular) &
699 kef(ix^d)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1,ix2-1,ix3)&
700 +ke(ix1-1,ix2,ix3)+ke(ix1-1,ix2-1,ix3))
701 {end do\}
702 end if
703 }
704 {^iftwod
705 if(idims==1) then
706 {do ix^db=ixamin^db,ixamax^db\}
707 ^d&bcf({ix^d},^d)=0.5d0*(bc(ix1,ix2,^d)+bc(ix1,ix2-1,^d))\
708 kaf(ix^d)=0.5d0*(ka(ix1,ix2)+ka(ix1,ix2-1))
709 if(fl%tc_perpendicular) &
710 kef(ix^d)=0.5d0*(ke(ix1,ix2)+ke(ix1,ix2-1))
711 {end do\}
712 else
713 {do ix^db=ixamin^db,ixamax^db\}
714 ^d&bcf({ix^d},^d)=0.5d0*(bc(ix1,ix2,^d)+bc(ix1-1,ix2,^d))\
715 kaf(ix^d)=0.5d0*(ka(ix1,ix2)+ka(ix1-1,ix2))
716 if(fl%tc_perpendicular) &
717 kef(ix^d)=0.5d0*(ke(ix1,ix2)+ke(ix1-1,ix2))
718 {end do\}
719 end if
720 }
721 ! eq (19)
722 ! temperature gradient at cell corner
723 {^ifthreed
724 if(idims==1) then
725 {do ix^db=ixcmin^db,ixcmax^db\}
726 qdd(ix^d)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1,ix2+1,ix3,idims)&
727 +gradt(ix1,ix2,ix3+1,idims)+gradt(ix1,ix2+1,ix3+1,idims))
728 {end do\}
729 else if(idims==2) then
730 {do ix^db=ixcmin^db,ixcmax^db\}
731 qdd(ix^d)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1+1,ix2,ix3,idims)&
732 +gradt(ix1,ix2,ix3+1,idims)+gradt(ix1+1,ix2,ix3+1,idims))
733 {end do\}
734 else
735 {do ix^db=ixcmin^db,ixcmax^db\}
736 qdd(ix^d)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1+1,ix2,ix3,idims)&
737 +gradt(ix1,ix2+1,ix3,idims)+gradt(ix1+1,ix2+1,ix3,idims))
738 {end do\}
739 end if
740 }
741 {^iftwod
742 if(idims==1) then
743 {do ix^db=ixcmin^db,ixcmax^db\}
744 qdd(ix^d)=0.5d0*(gradt(ix1,ix2,idims)+gradt(ix1,ix2+1,idims))
745 {end do\}
746 else
747 {do ix^db=ixcmin^db,ixcmax^db\}
748 qdd(ix^d)=0.5d0*(gradt(ix1,ix2,idims)+gradt(ix1+1,ix2,idims))
749 {end do\}
750 end if
751 }
752 ! eq (21)
753 {^ifthreed
754 if(idims==1) then
755 {do ix^db=ixamin^db,ixamax^db\}
756 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
757 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
758 if(qdd(ix^d)<minq) then
759 qd(ix^d,1)=minq
760 else if(qdd(ix^d)>maxq) then
761 qd(ix^d,1)=maxq
762 else
763 qd(ix^d,1)=qdd(ix^d)
764 end if
765 if(qdd(ix1,ix2-1,ix3)<minq) then
766 qd(ix^d,2)=minq
767 else if(qdd(ix1,ix2-1,ix3)>maxq) then
768 qd(ix^d,2)=maxq
769 else
770 qd(ix^d,2)=qdd(ix1,ix2-1,ix3)
771 end if
772 if(qdd(ix1,ix2,ix3-1)<minq) then
773 qd(ix^d,3)=minq
774 else if(qdd(ix1,ix2,ix3-1)>maxq) then
775 qd(ix^d,3)=maxq
776 else
777 qd(ix^d,3)=qdd(ix1,ix2,ix3-1)
778 end if
779 if(qdd(ix1,ix2-1,ix3-1)<minq) then
780 qd(ix^d,4)=minq
781 else if(qdd(ix1,ix2-1,ix3-1)>maxq) then
782 qd(ix^d,4)=maxq
783 else
784 qd(ix^d,4)=qdd(ix1,ix2-1,ix3-1)
785 end if
786 qvec(ix^d,idims)=kaf(ix^d)*0.25d0*(bc(ix^d,idims)**2*qd(ix^d,1)+bc(ix1,ix2-1,ix3,idims)**2*qd(ix^d,2)&
787 +bc(ix1,ix2,ix3-1,idims)**2*qd(ix^d,3)+bc(ix1,ix2-1,ix3-1,idims)**2*qd(ix^d,4))
788 if(fl%tc_perpendicular) &
789 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.25d0*(qd(ix^d,1)+qd(ix^d,2)+qd(ix^d,3)+qd(ix^d,4))
790 {end do\}
791 else if(idims==2) then
792 {do ix^db=ixamin^db,ixamax^db\}
793 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
794 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
795 if(qdd(ix^d)<minq) then
796 qd(ix^d,1)=minq
797 else if(qdd(ix^d)>maxq) then
798 qd(ix^d,1)=maxq
799 else
800 qd(ix^d,1)=qdd(ix^d)
801 end if
802 if(qdd(ix1-1,ix2,ix3)<minq) then
803 qd(ix^d,2)=minq
804 else if(qdd(ix1-1,ix2,ix3)>maxq) then
805 qd(ix^d,2)=maxq
806 else
807 qd(ix^d,2)=qdd(ix1-1,ix2,ix3)
808 end if
809 if(qdd(ix1,ix2,ix3-1)<minq) then
810 qd(ix^d,3)=minq
811 else if(qdd(ix1,ix2,ix3-1)>maxq) then
812 qd(ix^d,3)=maxq
813 else
814 qd(ix^d,3)=qdd(ix1,ix2,ix3-1)
815 end if
816 if(qdd(ix1-1,ix2,ix3-1)<minq) then
817 qd(ix^d,4)=minq
818 else if(qdd(ix1-1,ix2,ix3-1)>maxq) then
819 qd(ix^d,4)=maxq
820 else
821 qd(ix^d,4)=qdd(ix1-1,ix2,ix3-1)
822 end if
823 qvec(ix^d,idims)=kaf(ix^d)*0.25d0*(bc(ix^d,idims)**2*qd(ix^d,1)+bc(ix1-1,ix2,ix3,idims)**2*qd(ix^d,2)&
824 +bc(ix1,ix2,ix3-1,idims)**2*qd(ix^d,3)+bc(ix1-1,ix2,ix3-1,idims)**2*qd(ix^d,4))
825 if(fl%tc_perpendicular) &
826 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.25d0*(qd(ix^d,1)+qd(ix^d,2)+qd(ix^d,3)+qd(ix^d,4))
827 {end do\}
828 else
829 {do ix^db=ixamin^db,ixamax^db\}
830 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
831 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
832 if(qdd(ix^d)<minq) then
833 qd(ix^d,1)=minq
834 else if(qdd(ix^d)>maxq) then
835 qd(ix^d,1)=maxq
836 else
837 qd(ix^d,1)=qdd(ix^d)
838 end if
839 if(qdd(ix1-1,ix2,ix3)<minq) then
840 qd(ix^d,2)=minq
841 else if(qdd(ix1-1,ix2,ix3)>maxq) then
842 qd(ix^d,2)=maxq
843 else
844 qd(ix^d,2)=qdd(ix1-1,ix2,ix3)
845 end if
846 if(qdd(ix1,ix2-1,ix3)<minq) then
847 qd(ix^d,3)=minq
848 else if(qdd(ix1,ix2-1,ix3)>maxq) then
849 qd(ix^d,3)=maxq
850 else
851 qd(ix^d,3)=qdd(ix1,ix2-1,ix3)
852 end if
853 if(qdd(ix1-1,ix2-1,ix3)<minq) then
854 qd(ix^d,4)=minq
855 else if(qdd(ix1-1,ix2-1,ix3)>maxq) then
856 qd(ix^d,4)=maxq
857 else
858 qd(ix^d,4)=qdd(ix1-1,ix2-1,ix3)
859 end if
860 qvec(ix^d,idims)=kaf(ix^d)*0.25d0*(bc(ix^d,idims)**2*qd(ix^d,1)+bc(ix1-1,ix2,ix3,idims)**2*qd(ix^d,2)&
861 +bc(ix1,ix2-1,ix3,idims)**2*qd(ix^d,3)+bc(ix1-1,ix2-1,ix3,idims)**2*qd(ix^d,4))
862 if(fl%tc_perpendicular) &
863 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.25d0*(qd(ix^d,1)+qd(ix^d,2)+qd(ix^d,3)+qd(ix^d,4))
864 {end do\}
865 end if
866 }
867 {^iftwod
868 if(idims==1) then
869 {do ix^db=ixamin^db,ixamax^db\}
870 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
871 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
872 if(qdd(ix^d)<minq) then
873 qd(ix^d,1)=minq
874 else if(qdd(ix^d)>maxq) then
875 qd(ix^d,1)=maxq
876 else
877 qd(ix^d,1)=qdd(ix^d)
878 end if
879 if(qdd(ix1,ix2-1)<minq) then
880 qd(ix^d,2)=minq
881 else if(qdd(ix1,ix2-1)>maxq) then
882 qd(ix^d,2)=maxq
883 else
884 qd(ix^d,2)=qdd(ix1,ix2-1)
885 end if
886 qvec(ix^d,idims)=kaf(ix^d)*0.5d0*(bc(ix1,ix2,idims)**2*qd(ix^d,1)+bc(ix1,ix2-1,idims)**2*qd(ix^d,2))
887 if(fl%tc_perpendicular) &
888 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.5d0*(qd(ix^d,1)+qd(ix^d,2))
889 {end do\}
890 else
891 {do ix^db=ixamin^db,ixamax^db\}
892 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
893 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
894 if(qdd(ix^d)<minq) then
895 qd(ix^d,1)=minq
896 else if(qdd(ix^d)>maxq) then
897 qd(ix^d,1)=maxq
898 else
899 qd(ix^d,1)=qdd(ix^d)
900 end if
901 if(qdd(ix1-1,ix2)<minq) then
902 qd(ix^d,2)=minq
903 else if(qdd(ix1-1,ix2)>maxq) then
904 qd(ix^d,2)=maxq
905 else
906 qd(ix^d,2)=qdd(ix1-1,ix2)
907 end if
908 qvec(ix^d,idims)=kaf(ix^d)*0.5d0*(bc(ix1,ix2,idims)**2*qd(ix^d,1)+bc(ix1-1,ix2,idims)**2*qd(ix^d,2))
909 if(fl%tc_perpendicular) &
910 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.5d0*(qd(ix^d,1)+qd(ix^d,2))
911 {end do\}
912 end if
913 }
914 ! calculate normal of magnetic field
915 ixb^l=ixa^l+kr(idims,^d);
916 bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
917 ! limited transverse component, eq (17)
918 ixbmin^d=ixamin^d;
919 ixbmax^d=ixamax^d+kr(idims,^d);
920 do idir=1,ndim
921 if(idir==idims) cycle
922 qdd(ixi^s)=slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
923 qdd(ixi^s)=slope_limiter(qdd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
924 qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qdd(ixa^s)
925 end do
926 if(fl%tc_saturate) then
927 ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
928 ! unsigned saturated TC flux = 5 phi rho c**3, c=sqrt(p/rho) is isothermal sound speed, phi=1.1
929 ixb^l=ixa^l+kr(idims,^d);
930 qdd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bnorm(ixa^s))
931 {do ix^db=ixamin^db,ixamax^db\}
932 if(dabs(qvec(ix^d,idims))>qdd(ix^d)) then
933 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qdd(ix^d)
934 end if
935 {end do\}
936 end if
937 end do
938 end if
939 end subroutine set_source_tc_mhd
940
941 subroutine set_source_tc_mhd_geo(ixI^L,ixO^L,w,x,fl,qvec,rho,Te,alpha)
943 integer, intent(in) :: ixI^L, ixO^L
944 double precision, intent(in) :: x(ixI^S,1:ndim)
945 double precision, intent(in) :: w(ixI^S,1:nw)
946 type(tc_fluid), intent(in) :: fl
947 double precision, intent(in) :: rho(ixI^S),Te(ixI^S)
948 double precision, intent(in) :: alpha
949 double precision, intent(out) :: qvec(ixI^S,1:ndim)
950
951 !! qdd store the heat conduction energy changing rate
952 double precision, dimension(ixI^S,1:ndim) :: mf,Bc,Bcf,gradT
953 double precision, dimension(ixI^S) :: ka,kaf,ke,kef,qdd,Bnorm
954 double precision :: minq,maxq,qd(ixI^S,2**(ndim-1)), blocal(ndim)
955 integer :: idims,idir,ix^D,ix^L,ixC^L,ixA^L,ixB^L
956
957 ix^l=ixo^l^ladd1;
958
959 ! T gradient at cell faces
960 ! b unit vector mf: magnetic field direction vector
961 if(allocated(iw_mag)) then
962 if(b0field) then
963 {do ix^db=ixmin^db,ixmax^db\}
964 ^d&blocal(^d)=w({ix^d},iw_mag(^d))+block%B0({ix^d},^d,0)\
965 {^iftwod
966 if(blocal(1)/=0.d0) then
967 mf(ix^d,1)=sign(1.d0,blocal(1))/dsqrt(1.d0+(blocal(2)/blocal(1))**2)
968 else
969 mf(ix^d,1)=0.d0
970 end if
971 if(blocal(2)/=0.d0) then
972 mf(ix^d,2)=sign(1.d0,blocal(2))/dsqrt(1.d0+(blocal(1)/blocal(2))**2)
973 else
974 mf(ix^d,2)=0.d0
975 end if
976 }
977 {^ifthreed
978 if(blocal(1)/=0.d0) then
979 mf(ix^d,1)=sign(1.d0,blocal(1))/dsqrt(1.d0+(blocal(2)/blocal(1))**2+(blocal(3)/blocal(1))**2)
980 else
981 mf(ix^d,1)=0.d0
982 end if
983 if(blocal(2)/=0.d0) then
984 mf(ix^d,2)=sign(1.d0,blocal(2))/dsqrt(1.d0+(blocal(1)/blocal(2))**2+(blocal(3)/blocal(2))**2)
985 else
986 mf(ix^d,2)=0.d0
987 end if
988 if(blocal(3)/=0.d0) then
989 mf(ix^d,3)=sign(1.d0,blocal(3))/dsqrt(1.d0+(blocal(1)/blocal(3))**2+(blocal(2)/blocal(3))**2)
990 else
991 mf(ix^d,3)=0.d0
992 end if
993 }
994 {end do\}
995 else
996 {do ix^db=ixmin^db,ixmax^db\}
997 {^iftwod
998 if(w(ix^d,iw_mag(1))/=0.d0) then
999 mf(ix^d,1)=sign(1.d0,w(ix^d,iw_mag(1)))/dsqrt(1.d0+(w(ix^d,iw_mag(2))/w(ix^d,iw_mag(1)))**2)
1000 else
1001 mf(ix^d,1)=0.d0
1002 end if
1003 if(w(ix^d,iw_mag(2))/=0.d0) then
1004 mf(ix^d,2)=sign(1.d0,w(ix^d,iw_mag(2)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(2)))**2)
1005 else
1006 mf(ix^d,2)=0.d0
1007 end if
1008 }
1009 {^ifthreed
1010 if(w(ix^d,iw_mag(1))/=0.d0) then
1011 mf(ix^d,1)=sign(1.d0,w(ix^d,iw_mag(1)))/dsqrt(1.d0+(w(ix^d,iw_mag(2))/w(ix^d,iw_mag(1)))**2+&
1012 (w(ix^d,iw_mag(3))/w(ix^d,iw_mag(1)))**2)
1013 else
1014 mf(ix^d,1)=0.d0
1015 end if
1016 if(w(ix^d,iw_mag(2))/=0.d0) then
1017 mf(ix^d,2)=sign(1.d0,w(ix^d,iw_mag(2)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(2)))**2+&
1018 (w(ix^d,iw_mag(3))/w(ix^d,iw_mag(2)))**2)
1019 else
1020 mf(ix^d,2)=0.d0
1021 end if
1022 if(w(ix^d,iw_mag(3))/=0.d0) then
1023 mf(ix^d,3)=sign(1.d0,w(ix^d,iw_mag(3)))/dsqrt(1.d0+(w(ix^d,iw_mag(1))/w(ix^d,iw_mag(3)))**2+&
1024 (w(ix^d,iw_mag(2))/w(ix^d,iw_mag(3)))**2)
1025 else
1026 mf(ix^d,3)=0.d0
1027 end if
1028 }
1029 {end do\}
1030 end if
1031 else
1032 mf(ix^s,1:ndim)=block%B0(ix^s,1:ndim,0)
1033 endif
1034 ! ixC is cell-corner index
1035 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
1036 ! b unit vector at cell corner
1037 {^ifthreed
1038 do idims=1,3
1039 {do ix^db=ixcmin^db,ixcmax^db\}
1040 bc(ix^d,idims)=(mf(ix1,ix2,ix3,idims)*block%dvolume(ix1,ix2,ix3)+mf(ix1+1,ix2,ix3,idims)*block%dvolume(ix1+1,ix2,ix3)&
1041 +mf(ix1,ix2+1,ix3,idims)*block%dvolume(ix1,ix2+1,ix3)+mf(ix1+1,ix2+1,ix3,idims)*block%dvolume(ix1+1,ix2+1,ix3)&
1042 +mf(ix1,ix2,ix3+1,idims)*block%dvolume(ix1,ix2,ix3+1)+mf(ix1+1,ix2,ix3+1,idims)*block%dvolume(ix1+1,ix2,ix3+1)&
1043 +mf(ix1,ix2+1,ix3+1,idims)*block%dvolume(ix1,ix2+1,ix3+1)+mf(ix1+1,ix2+1,ix3+1,idims)*block%dvolume(ix1+1,ix2+1,ix3+1))&
1044 /(block%dvolume(ix1,ix2,ix3)+block%dvolume(ix1+1,ix2,ix3)+block%dvolume(ix1,ix2+1,ix3)+block%dvolume(ix1+1,ix2+1,ix3)&
1045 +block%dvolume(ix1,ix2,ix3+1)+block%dvolume(ix1+1,ix2,ix3+1)+block%dvolume(ix1,ix2+1,ix3+1)+block%dvolume(ix1+1,ix2+1,ix3+1))
1046 {end do\}
1047 end do
1048 }
1049 {^iftwod
1050 do idims=1,2
1051 {do ix^db=ixcmin^db,ixcmax^db\}
1052 bc(ix^d,idims)=(mf(ix1,ix2,idims)*block%dvolume(ix1,ix2)+mf(ix1+1,ix2,idims)*block%dvolume(ix1+1,ix2)&
1053 +mf(ix1,ix2+1,idims)*block%dvolume(ix1,ix2+1)+mf(ix1+1,ix2+1,idims)*block%dvolume(ix1+1,ix2+1))&
1054 /(block%dvolume(ix1,ix2)+block%dvolume(ix1+1,ix2)+block%dvolume(ix1,ix2+1)+block%dvolume(ix1+1,ix2+1))
1055 {end do\}
1056 end do
1057 }
1058 ! T gradient at cell faces
1059 do idims=1,ndim
1060 ixbmin^d=ixmin^d;
1061 ixbmax^d=ixmax^d-kr(idims,^d);
1062 call gradientf(te,x,ixi^l,ixb^l,idims,gradt(ixi^s,idims))
1063 end do
1064 if(fl%tc_constant) then
1065 if(fl%tc_perpendicular) then
1066 ka(ixc^s)=fl%tc_k_para-fl%tc_k_perp
1067 ke(ixc^s)=fl%tc_k_perp
1068 else
1069 ka(ixc^s)=fl%tc_k_para
1070 end if
1071 else
1072 ! conductivity at cell center
1073 if(phys_trac) then
1074 {do ix^db=ixmin^db,ixmax^db\}
1075 if(te(ix^d) < block%wextra(ix^d,fl%Tcoff_)) then
1076 qdd(ix^d)=fl%tc_k_para*dsqrt(block%wextra(ix^d,fl%Tcoff_)**5)
1077 else
1078 qdd(ix^d)=fl%tc_k_para*dsqrt(te(ix^d)**5)
1079 end if
1080 {end do\}
1081 else
1082 qdd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s)**5)
1083 end if
1084 ! cell corner parallel conductivity in ka
1085 {^ifthreed
1086 {do ix^db=ixcmin^db,ixcmax^db\}
1087 ka(ix^d)=(qdd(ix1,ix2,ix3)*block%dvolume(ix1,ix2,ix3)+qdd(ix1+1,ix2,ix3)*block%dvolume(ix1+1,ix2,ix3)&
1088 +qdd(ix1,ix2+1,ix3)*block%dvolume(ix1,ix2+1,ix3)+qdd(ix1+1,ix2+1,ix3)*block%dvolume(ix1+1,ix2+1,ix3)&
1089 +qdd(ix1,ix2,ix3+1)*block%dvolume(ix1,ix2,ix3+1)+qdd(ix1+1,ix2,ix3+1)*block%dvolume(ix1+1,ix2,ix3+1)&
1090 +qdd(ix1,ix2+1,ix3+1)*block%dvolume(ix1,ix2+1,ix3+1)+qdd(ix1+1,ix2+1,ix3+1)*block%dvolume(ix1+1,ix2+1,ix3+1))&
1091 /(block%dvolume(ix1,ix2,ix3)+block%dvolume(ix1+1,ix2,ix3)+block%dvolume(ix1,ix2+1,ix3)+block%dvolume(ix1+1,ix2+1,ix3)&
1092 +block%dvolume(ix1,ix2,ix3+1)+block%dvolume(ix1+1,ix2,ix3+1)+block%dvolume(ix1,ix2+1,ix3+1)+block%dvolume(ix1+1,ix2+1,ix3+1))
1093 {end do\}
1094 }
1095 {^iftwod
1096 {do ix^db=ixcmin^db,ixcmax^db\}
1097 ka(ix^d)=(qdd(ix1,ix2)*block%dvolume(ix1,ix2)+qdd(ix1+1,ix2)*block%dvolume(ix1+1,ix2)&
1098 +qdd(ix1,ix2+1)*block%dvolume(ix1,ix2+1)+qdd(ix1+1,ix2+1)*block%dvolume(ix1+1,ix2+1))&
1099 /(block%dvolume(ix1,ix2)+block%dvolume(ix1+1,ix2)+block%dvolume(ix1,ix2+1)+block%dvolume(ix1+1,ix2+1))
1100 {end do\}
1101 }
1102 ! compensate with perpendicular conductivity
1103 if(fl%tc_perpendicular) then
1104 if(b0field) then
1105 qdd(ix^s)=fl%tc_k_perp*rho(ix^s)**2/((^d&(w(ix^s,iw_mag(^d))+block%B0(ix^s,^d,0))**2+)*dsqrt(te(ix^s))+smalldouble)
1106 else
1107 qdd(ix^s)=fl%tc_k_perp*rho(ix^s)**2/((^d&w(ix^s,iw_mag(^d))**2+)*dsqrt(te(ix^s))+smalldouble)
1108 end if
1109 {^ifthreed
1110 {do ix^db=ixcmin^db,ixcmax^db\}
1111 ke(ix^d)=(qdd(ix1,ix2,ix3)*block%dvolume(ix1,ix2,ix3)+qdd(ix1+1,ix2,ix3)*block%dvolume(ix1+1,ix2,ix3)&
1112 +qdd(ix1,ix2+1,ix3)*block%dvolume(ix1,ix2+1,ix3)+qdd(ix1+1,ix2+1,ix3)*block%dvolume(ix1+1,ix2+1,ix3)&
1113 +qdd(ix1,ix2,ix3+1)*block%dvolume(ix1,ix2,ix3+1)+qdd(ix1+1,ix2,ix3+1)*block%dvolume(ix1+1,ix2,ix3+1)&
1114 +qdd(ix1,ix2+1,ix3+1)*block%dvolume(ix1,ix2+1,ix3+1)+qdd(ix1+1,ix2+1,ix3+1)*block%dvolume(ix1+1,ix2+1,ix3+1))&
1115 /(block%dvolume(ix1,ix2,ix3)+block%dvolume(ix1+1,ix2,ix3)+block%dvolume(ix1,ix2+1,ix3)+block%dvolume(ix1+1,ix2+1,ix3)&
1116 +block%dvolume(ix1,ix2,ix3+1)+block%dvolume(ix1+1,ix2,ix3+1)+block%dvolume(ix1,ix2+1,ix3+1)+block%dvolume(ix1+1,ix2+1,ix3+1))
1117 if(ke(ix^d)<ka(ix^d)) then
1118 ka(ix^d)=ka(ix^d)-ke(ix^d)
1119 else
1120 ke(ix^d)=ka(ix^d)
1121 ka(ix^d)=0.d0
1122 end if
1123 {end do\}
1124 }
1125 {^iftwod
1126 {do ix^db=ixcmin^db,ixcmax^db\}
1127 ke(ix^d)=(qdd(ix1,ix2)*block%dvolume(ix1,ix2)+qdd(ix1+1,ix2)*block%dvolume(ix1+1,ix2)&
1128 +qdd(ix1,ix2+1)*block%dvolume(ix1,ix2+1)+qdd(ix1+1,ix2+1)*block%dvolume(ix1+1,ix2+1))&
1129 /(block%dvolume(ix1,ix2)+block%dvolume(ix1+1,ix2)+block%dvolume(ix1,ix2+1)+block%dvolume(ix1+1,ix2+1))
1130 if(ke(ix^d)<ka(ix^d)) then
1131 ka(ix^d)=ka(ix^d)-ke(ix^d)
1132 else
1133 ke(ix^d)=ka(ix^d)
1134 ka(ix^d)=0.d0
1135 end if
1136 {end do\}
1137 }
1138 end if
1139 end if
1140 if(fl%tc_slope_limiter==0) then
1141 ! calculate thermal conduction flux with symmetric scheme
1142 do idims=1,ndim
1143 !qdd corner values
1144 qdd=0.d0
1145 {do ix^db=0,1 \}
1146 if({ ix^d==0 .and. ^d==idims | .or.}) then
1147 ixbmin^d=ixcmin^d+ix^d;
1148 ixbmax^d=ixcmax^d+ix^d;
1149 qdd(ixc^s)=qdd(ixc^s)+gradt(ixb^s,idims)
1150 end if
1151 {end do\}
1152 ! temperature gradient at cell corner
1153 qvec(ixc^s,idims)=qdd(ixc^s)*0.5d0**(ndim-1)
1154 end do
1155 ! b grad T at cell corner
1156 qdd(ixc^s)=sum(qvec(ixc^s,1:ndim)*bc(ixc^s,1:ndim),dim=ndim+1)
1157 do idims=1,ndim
1158 ! TC flux at cell corner
1159 gradt(ixc^s,idims)=ka(ixc^s)*bc(ixc^s,idims)*qdd(ixc^s)
1160 if(fl%tc_perpendicular) gradt(ixc^s,idims)=gradt(ixc^s,idims)+ke(ixc^s)*qvec(ixc^s,idims)
1161 end do
1162 ! TC flux at cell face
1163 qvec=0.d0
1164 do idims=1,ndim
1165 ixb^l=ixo^l-kr(idims,^d);
1166 ixamax^d=ixomax^d; ixamin^d=ixbmin^d;
1167 {do ix^db=0,1 \}
1168 if({ ix^d==0 .and. ^d==idims | .or.}) then
1169 ixbmin^d=ixamin^d-ix^d;
1170 ixbmax^d=ixamax^d-ix^d;
1171 qvec(ixa^s,idims)=qvec(ixa^s,idims)+gradt(ixb^s,idims)
1172 end if
1173 {end do\}
1174 qvec(ixa^s,idims)=qvec(ixa^s,idims)*0.5d0**(ndim-1)
1175 if(fl%tc_saturate) then
1176 ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
1177 ! unsigned saturated TC flux = 5 phi rho c**3, c is isothermal sound speed
1178 bcf=0.d0
1179 {do ix^db=0,1 \}
1180 if({ ix^d==0 .and. ^d==idims | .or.}) then
1181 ixbmin^d=ixamin^d-ix^d;
1182 ixbmax^d=ixamax^d-ix^d;
1183 bcf(ixa^s,idims)=bcf(ixa^s,idims)+bc(ixb^s,idims)
1184 end if
1185 {end do\}
1186 ! averaged b at face centers
1187 bcf(ixa^s,idims)=bcf(ixa^s,idims)*0.5d0**(ndim-1)
1188 ixb^l=ixa^l+kr(idims,^d);
1189 qdd(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))
1190 {do ix^db=ixamin^db,ixamax^db\}
1191 if(dabs(qvec(ix^d,idims))>qdd(ix^d)) then
1192 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qdd(ix^d)
1193 end if
1194 {end do\}
1195 end if
1196 end do
1197 else
1198 ! calculate thermal conduction flux with slope-limited symmetric scheme
1199 do idims=1,ndim
1200 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
1201 {^ifthreed
1202 if(idims==1) then
1203 {do ix^db=ixamin^db,ixamax^db\}
1204 ! averaged b at face centers
1205 ^d&bcf({ix^d},^d)=0.25d0*(bc({ix^d},^d)+bc(ix1,ix2-1,ix3,^d)&
1206 +bc(ix1,ix2,ix3-1,^d)+bc(ix1,ix2-1,ix3-1,^d))\
1207 kaf(ix^d)=0.25d0*(ka(ix1,ix2,ix3)+ka(ix1,ix2-1,ix3)&
1208 +ka(ix1,ix2,ix3-1)+ka(ix1,ix2-1,ix3-1))
1209 ! averaged thermal conductivity at face centers
1210 if(fl%tc_perpendicular) &
1211 kef(ix^d)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1,ix2-1,ix3)&
1212 +ke(ix1,ix2,ix3-1)+ke(ix1,ix2-1,ix3-1))
1213 {end do\}
1214 else if(idims==2) then
1215 {do ix^db=ixamin^db,ixamax^db\}
1216 ^d&bcf({ix^d},^d)=0.25d0*(bc({ix^d},^d)+bc(ix1-1,ix2,ix3,^d)&
1217 +bc(ix1,ix2,ix3-1,^d)+bc(ix1-1,ix2,ix3-1,^d))\
1218 kaf(ix^d)=0.25d0*(ka(ix1,ix2,ix3)+ka(ix1-1,ix2,ix3)&
1219 +ka(ix1,ix2,ix3-1)+ka(ix1-1,ix2,ix3-1))
1220 if(fl%tc_perpendicular) &
1221 kef(ix^d)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1-1,ix2,ix3)&
1222 +ke(ix1,ix2,ix3-1)+ke(ix1-1,ix2,ix3-1))
1223 {end do\}
1224 else
1225 {do ix^db=ixamin^db,ixamax^db\}
1226 ^d&bcf({ix^d},^d)=0.25d0*(bc({ix^d},^d)+bc(ix1,ix2-1,ix3,^d)&
1227 +bc(ix1-1,ix2,ix3,^d)+bc(ix1-1,ix2-1,ix3,^d))\
1228 kaf(ix^d)=0.25d0*(ka(ix1,ix2,ix3)+ka(ix1,ix2-1,ix3)&
1229 +ka(ix1-1,ix2,ix3)+ka(ix1-1,ix2-1,ix3))
1230 if(fl%tc_perpendicular) &
1231 kef(ix^d)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1,ix2-1,ix3)&
1232 +ke(ix1-1,ix2,ix3)+ke(ix1-1,ix2-1,ix3))
1233 {end do\}
1234 end if
1235 }
1236 {^iftwod
1237 if(idims==1) then
1238 {do ix^db=ixamin^db,ixamax^db\}
1239 ^d&bcf({ix^d},^d)=0.5d0*(bc(ix1,ix2,^d)+bc(ix1,ix2-1,^d))\
1240 kaf(ix^d)=0.5d0*(ka(ix1,ix2)+ka(ix1,ix2-1))
1241 if(fl%tc_perpendicular) &
1242 kef(ix^d)=0.5d0*(ke(ix1,ix2)+ke(ix1,ix2-1))
1243 {end do\}
1244 else
1245 {do ix^db=ixamin^db,ixamax^db\}
1246 ^d&bcf({ix^d},^d)=0.5d0*(bc(ix1,ix2,^d)+bc(ix1-1,ix2,^d))\
1247 kaf(ix^d)=0.5d0*(ka(ix1,ix2)+ka(ix1-1,ix2))
1248 if(fl%tc_perpendicular) &
1249 kef(ix^d)=0.5d0*(ke(ix1,ix2)+ke(ix1-1,ix2))
1250 {end do\}
1251 end if
1252 }
1253 ! eq (19)
1254 ! temperature gradient at cell corner
1255 {^ifthreed
1256 if(idims==1) then
1257 {do ix^db=ixcmin^db,ixcmax^db\}
1258 qdd(ix^d)=(gradt(ix1,ix2,ix3,idims)*block%surfaceC(ix1,ix2,ix3,1)&
1259 +gradt(ix1,ix2+1,ix3,idims)*block%surfaceC(ix1,ix2+1,ix3,1)&
1260 +gradt(ix1,ix2,ix3+1,idims)*block%surfaceC(ix1,ix2,ix3+1,1)&
1261 +gradt(ix1,ix2+1,ix3+1,idims)*block%surfaceC(ix1,ix2+1,ix3+1,1))/&
1262 (block%surfaceC(ix1,ix2,ix3,1)+block%surfaceC(ix1,ix2+1,ix3,1)&
1263 +block%surfaceC(ix1,ix2,ix3+1,1)+block%surfaceC(ix1,ix2+1,ix3+1,1))
1264 {end do\}
1265 else if(idims==2) then
1266 {do ix^db=ixcmin^db,ixcmax^db\}
1267 qdd(ix^d)=(gradt(ix1,ix2,ix3,idims)*block%surfaceC(ix1,ix2,ix3,2)&
1268 +gradt(ix1+1,ix2,ix3,idims)*block%surfaceC(ix1+1,ix2,ix3,2)&
1269 +gradt(ix1,ix2,ix3+1,idims)*block%surfaceC(ix1,ix2,ix3+1,2)&
1270 +gradt(ix1+1,ix2,ix3+1,idims)*block%surfaceC(ix1+1,ix2,ix3+1,2))/&
1271 (block%surfaceC(ix1,ix2,ix3,2)+block%surfaceC(ix1+1,ix2,ix3,2)&
1272 +block%surfaceC(ix1,ix2,ix3+1,2)+block%surfaceC(ix1+1,ix2,ix3+1,2)+1.d-300)
1273 {end do\}
1274 else
1275 {do ix^db=ixcmin^db,ixcmax^db\}
1276 qdd(ix^d)=(gradt(ix1,ix2,ix3,idims)*block%surfaceC(ix1,ix2,ix3,3)&
1277 +gradt(ix1+1,ix2,ix3,idims)*block%surfaceC(ix1+1,ix2,ix3,3)&
1278 +gradt(ix1,ix2+1,ix3,idims)*block%surfaceC(ix1,ix2+1,ix3,3)&
1279 +gradt(ix1+1,ix2+1,ix3,idims)*block%surfaceC(ix1+1,ix2+1,ix3,3))/&
1280 (block%surfaceC(ix1,ix2,ix3,3)+block%surfaceC(ix1+1,ix2,ix3,3)&
1281 +block%surfaceC(ix1,ix2+1,ix3,3)+block%surfaceC(ix1+1,ix2+1,ix3,3))
1282 {end do\}
1283 end if
1284 }
1285 {^iftwod
1286 if(idims==1) then
1287 {do ix^db=ixcmin^db,ixcmax^db\}
1288 qdd(ix^d)=(gradt(ix1,ix2,idims)*block%surfaceC(ix1,ix2,1)&
1289 +gradt(ix1,ix2+1,idims)*block%surfaceC(ix1,ix2+1,1))/&
1290 (block%surfaceC(ix1,ix2,1)+block%surfaceC(ix1,ix2+1,1))
1291 {end do\}
1292 else
1293 {do ix^db=ixcmin^db,ixcmax^db\}
1294 qdd(ix^d)=(gradt(ix1,ix2,idims)*block%surfaceC(ix1,ix2,2)&
1295 +gradt(ix1+1,ix2,idims)*block%surfaceC(ix1+1,ix2,2))/&
1296 (block%surfaceC(ix1,ix2,2)+block%surfaceC(ix1+1,ix2,2))
1297 {end do\}
1298 end if
1299 }
1300 ! eq (21)
1301 {^ifthreed
1302 if(idims==1) then
1303 {do ix^db=ixamin^db,ixamax^db\}
1304 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1305 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1306 if(qdd(ix^d)<minq) then
1307 qd(ix^d,1)=minq
1308 else if(qdd(ix^d)>maxq) then
1309 qd(ix^d,1)=maxq
1310 else
1311 qd(ix^d,1)=qdd(ix^d)
1312 end if
1313 if(qdd(ix1,ix2-1,ix3)<minq) then
1314 qd(ix^d,2)=minq
1315 else if(qdd(ix1,ix2-1,ix3)>maxq) then
1316 qd(ix^d,2)=maxq
1317 else
1318 qd(ix^d,2)=qdd(ix1,ix2-1,ix3)
1319 end if
1320 if(qdd(ix1,ix2,ix3-1)<minq) then
1321 qd(ix^d,3)=minq
1322 else if(qdd(ix1,ix2,ix3-1)>maxq) then
1323 qd(ix^d,3)=maxq
1324 else
1325 qd(ix^d,3)=qdd(ix1,ix2,ix3-1)
1326 end if
1327 if(qdd(ix1,ix2-1,ix3-1)<minq) then
1328 qd(ix^d,4)=minq
1329 else if(qdd(ix1,ix2-1,ix3-1)>maxq) then
1330 qd(ix^d,4)=maxq
1331 else
1332 qd(ix^d,4)=qdd(ix1,ix2-1,ix3-1)
1333 end if
1334 qvec(ix^d,idims)=kaf(ix^d)*0.25d0*(bc(ix^d,idims)**2*qd(ix^d,1)+bc(ix1,ix2-1,ix3,idims)**2*qd(ix^d,2)&
1335 +bc(ix1,ix2,ix3-1,idims)**2*qd(ix^d,3)+bc(ix1,ix2-1,ix3-1,idims)**2*qd(ix^d,4))
1336 if(fl%tc_perpendicular) &
1337 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.25d0*(qd(ix^d,1)+qd(ix^d,2)+qd(ix^d,3)+qd(ix^d,4))
1338 {end do\}
1339 else if(idims==2) then
1340 {do ix^db=ixamin^db,ixamax^db\}
1341 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1342 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1343 if(qdd(ix^d)<minq) then
1344 qd(ix^d,1)=minq
1345 else if(qdd(ix^d)>maxq) then
1346 qd(ix^d,1)=maxq
1347 else
1348 qd(ix^d,1)=qdd(ix^d)
1349 end if
1350 if(qdd(ix1-1,ix2,ix3)<minq) then
1351 qd(ix^d,2)=minq
1352 else if(qdd(ix1-1,ix2,ix3)>maxq) then
1353 qd(ix^d,2)=maxq
1354 else
1355 qd(ix^d,2)=qdd(ix1-1,ix2,ix3)
1356 end if
1357 if(qdd(ix1,ix2,ix3-1)<minq) then
1358 qd(ix^d,3)=minq
1359 else if(qdd(ix1,ix2,ix3-1)>maxq) then
1360 qd(ix^d,3)=maxq
1361 else
1362 qd(ix^d,3)=qdd(ix1,ix2,ix3-1)
1363 end if
1364 if(qdd(ix1-1,ix2,ix3-1)<minq) then
1365 qd(ix^d,4)=minq
1366 else if(qdd(ix1-1,ix2,ix3-1)>maxq) then
1367 qd(ix^d,4)=maxq
1368 else
1369 qd(ix^d,4)=qdd(ix1-1,ix2,ix3-1)
1370 end if
1371 qvec(ix^d,idims)=kaf(ix^d)*0.25d0*(bc(ix^d,idims)**2*qd(ix^d,1)+bc(ix1-1,ix2,ix3,idims)**2*qd(ix^d,2)&
1372 +bc(ix1,ix2,ix3-1,idims)**2*qd(ix^d,3)+bc(ix1-1,ix2,ix3-1,idims)**2*qd(ix^d,4))
1373 if(fl%tc_perpendicular) &
1374 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.25d0*(qd(ix^d,1)+qd(ix^d,2)+qd(ix^d,3)+qd(ix^d,4))
1375 {end do\}
1376 else
1377 {do ix^db=ixamin^db,ixamax^db\}
1378 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1379 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1380 if(qdd(ix^d)<minq) then
1381 qd(ix^d,1)=minq
1382 else if(qdd(ix^d)>maxq) then
1383 qd(ix^d,1)=maxq
1384 else
1385 qd(ix^d,1)=qdd(ix^d)
1386 end if
1387 if(qdd(ix1-1,ix2,ix3)<minq) then
1388 qd(ix^d,2)=minq
1389 else if(qdd(ix1-1,ix2,ix3)>maxq) then
1390 qd(ix^d,2)=maxq
1391 else
1392 qd(ix^d,2)=qdd(ix1-1,ix2,ix3)
1393 end if
1394 if(qdd(ix1,ix2-1,ix3)<minq) then
1395 qd(ix^d,3)=minq
1396 else if(qdd(ix1,ix2-1,ix3)>maxq) then
1397 qd(ix^d,3)=maxq
1398 else
1399 qd(ix^d,3)=qdd(ix1,ix2-1,ix3)
1400 end if
1401 if(qdd(ix1-1,ix2-1,ix3)<minq) then
1402 qd(ix^d,4)=minq
1403 else if(qdd(ix1-1,ix2-1,ix3)>maxq) then
1404 qd(ix^d,4)=maxq
1405 else
1406 qd(ix^d,4)=qdd(ix1-1,ix2-1,ix3)
1407 end if
1408 qvec(ix^d,idims)=kaf(ix^d)*0.25d0*(bc(ix^d,idims)**2*qd(ix^d,1)+bc(ix1-1,ix2,ix3,idims)**2*qd(ix^d,2)&
1409 +bc(ix1,ix2-1,ix3,idims)**2*qd(ix^d,3)+bc(ix1-1,ix2-1,ix3,idims)**2*qd(ix^d,4))
1410 if(fl%tc_perpendicular) &
1411 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.25d0*(qd(ix^d,1)+qd(ix^d,2)+qd(ix^d,3)+qd(ix^d,4))
1412 {end do\}
1413 end if
1414 }
1415 {^iftwod
1416 if(idims==1) then
1417 {do ix^db=ixamin^db,ixamax^db\}
1418 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1419 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1420 if(qdd(ix^d)<minq) then
1421 qd(ix^d,1)=minq
1422 else if(qdd(ix^d)>maxq) then
1423 qd(ix^d,1)=maxq
1424 else
1425 qd(ix^d,1)=qdd(ix^d)
1426 end if
1427 if(qdd(ix1,ix2-1)<minq) then
1428 qd(ix^d,2)=minq
1429 else if(qdd(ix1,ix2-1)>maxq) then
1430 qd(ix^d,2)=maxq
1431 else
1432 qd(ix^d,2)=qdd(ix1,ix2-1)
1433 end if
1434 qvec(ix^d,idims)=kaf(ix^d)*0.5d0*(bc(ix1,ix2,idims)**2*qd(ix^d,1)+bc(ix1,ix2-1,idims)**2*qd(ix^d,2))
1435 if(fl%tc_perpendicular) &
1436 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.5d0*(qd(ix^d,1)+qd(ix^d,2))
1437 {end do\}
1438 else
1439 {do ix^db=ixamin^db,ixamax^db\}
1440 minq=min(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1441 maxq=max(alpha*gradt(ix^d,idims),gradt(ix^d,idims)/alpha)
1442 if(qdd(ix^d)<minq) then
1443 qd(ix^d,1)=minq
1444 else if(qdd(ix^d)>maxq) then
1445 qd(ix^d,1)=maxq
1446 else
1447 qd(ix^d,1)=qdd(ix^d)
1448 end if
1449 if(qdd(ix1-1,ix2)<minq) then
1450 qd(ix^d,2)=minq
1451 else if(qdd(ix1-1,ix2)>maxq) then
1452 qd(ix^d,2)=maxq
1453 else
1454 qd(ix^d,2)=qdd(ix1-1,ix2)
1455 end if
1456 qvec(ix^d,idims)=kaf(ix^d)*0.5d0*(bc(ix1,ix2,idims)**2*qd(ix^d,1)+bc(ix1-1,ix2,idims)**2*qd(ix^d,2))
1457 if(fl%tc_perpendicular) &
1458 qvec(ix^d,idims)=qvec(ix^d,idims)+kef(ix^d)*0.5d0*(qd(ix^d,1)+qd(ix^d,2))
1459 {end do\}
1460 end if
1461 }
1462 ! calculate normal of magnetic field
1463 ixb^l=ixa^l+kr(idims,^d);
1464 bnorm(ixa^s)=0.5d0*(mf(ixa^s,idims)+mf(ixb^s,idims))
1465 ! limited transverse component, eq (17)
1466 ixbmin^d=ixamin^d;
1467 ixbmax^d=ixamax^d+kr(idims,^d);
1468 do idir=1,ndim
1469 if(idir==idims) cycle
1470 qdd(ixi^s)=slope_limiter(gradt(ixi^s,idir),ixi^l,ixb^l,idir,-1,fl%tc_slope_limiter)
1471 qdd(ixi^s)=slope_limiter(qdd,ixi^l,ixa^l,idims,1,fl%tc_slope_limiter)
1472 qvec(ixa^s,idims)=qvec(ixa^s,idims)+kaf(ixa^s)*bnorm(ixa^s)*bcf(ixa^s,idir)*qdd(ixa^s)
1473 end do
1474 if(fl%tc_saturate) then
1475 ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
1476 ! unsigned saturated TC flux = 5 phi rho c**3, c=sqrt(p/rho) is isothermal sound speed, phi=1.1
1477 ixb^l=ixa^l+kr(idims,^d);
1478 qdd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3*dabs(bnorm(ixa^s))
1479 {do ix^db=ixamin^db,ixamax^db\}
1480 if(dabs(qvec(ix^d,idims))>qdd(ix^d)) then
1481 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qdd(ix^d)
1482 end if
1483 {end do\}
1484 end if
1485 end do
1486 end if
1487 end subroutine set_source_tc_mhd_geo
1488
1489 function slope_limiter(f,ixI^L,ixO^L,idims,pm,tc_slope_limiter) result(lf)
1491 integer, intent(in) :: ixi^l, ixo^l, idims, pm
1492 double precision, intent(in) :: f(ixi^s)
1493 double precision :: lf(ixi^s)
1494 integer, intent(in) :: tc_slope_limiter
1495
1496 double precision :: signf(ixi^s)
1497 integer :: ixb^l
1498
1499 ixb^l=ixo^l+pm*kr(idims,^d);
1500 signf(ixo^s)=sign(1.d0,f(ixo^s))
1501 select case(tc_slope_limiter)
1502 case(1)
1503 ! 'MC' monotonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
1504 lf(ixo^s)=two*signf(ixo^s)* &
1505 max(zero,min(dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s),&
1506 signf(ixo^s)*quarter*(f(ixb^s)+f(ixo^s))))
1507 case(2)
1508 ! 'minmod' limiter
1509 lf(ixo^s)=signf(ixo^s)*max(0.d0,min(abs(f(ixo^s)),signf(ixo^s)*f(ixb^s)))
1510 case(3)
1511 ! 'superbee' Roe superbee limiter (eq.3.51i)
1512 lf(ixo^s)=signf(ixo^s)* &
1513 max(zero,min(two*dabs(f(ixo^s)),signf(ixo^s)*f(ixb^s)),&
1514 min(dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s)))
1515 case(4)
1516 ! 'koren' Barry Koren Right variant
1517 lf(ixo^s)=signf(ixo^s)* &
1518 max(zero,min(two*dabs(f(ixo^s)),two*signf(ixo^s)*f(ixb^s),&
1519 (two*f(ixb^s)*signf(ixo^s)+dabs(f(ixo^s)))*third))
1520 case default
1521 call mpistop("Unknown slope limiter for thermal conduction")
1522 end select
1523 end function slope_limiter
1524
1525 !> Get the explicit timestep for the TC (hd implementation)
1526 !> Note: also used in 1D MHD (or for neutrals in twofl)
1527 function get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x,fl) result(dtnew)
1528 ! Check diffusion time limit dt < dx_i**2 / ((gamma-1)*tc_k_para/rho)
1530
1531 integer, intent(in) :: ixi^l, ixo^l
1532 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
1533 double precision, intent(in) :: w(ixi^s,1:nw)
1534 type(tc_fluid), intent(in) :: fl
1535 double precision :: dtnew
1536
1537 double precision :: tmp(ixo^s),tmp2(ixo^s),te(ixi^s),rho(ixi^s),hfs(ixo^s),gradt(ixi^s)
1538 double precision :: dtdiff_tcond,maxtmp2
1539 integer :: idim
1540
1541 call fl%get_temperature_from_conserved(w,x,ixi^l,ixi^l,te)
1542 call fl%get_rho(w,x,ixi^l,ixo^l,rho)
1543
1544 if(fl%tc_constant) then
1545 tmp(ixo^s)=fl%tc_k_para/rho(ixo^s)
1546 else
1547 if(fl%tc_saturate) then
1548 ! Kannan 2016 MN 458, 410
1549 ! 3^1.5*kB^2/(4*sqrt(pi)*e^4)
1550 ! l_mfpe=3.d0**1.5d0*kB_cgs**2/(4.d0*sqrt(dpi)*e_cgs**4*37.d0)=7093.9239487765044d0
1551 tmp2(ixo^s)=te(ixo^s)**2/rho(ixo^s)*7093.9239487765044d0*unit_temperature**2/(unit_numberdensity*unit_length)
1552 hfs=0.d0
1553 do idim=1,ndim
1554 call gradient(te,ixi^l,ixo^l,idim,gradt)
1555 hfs(ixo^s)=hfs(ixo^s)+gradt(ixo^s)**2
1556 end do
1557 ! kappa=kappa_Spitzer/(1+4.2*l_mfpe/(T/|gradT|))
1558 tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/(rho(ixo^s)*(1.d0+4.2d0*tmp2(ixo^s)*dsqrt(hfs(ixo^s))/te(ixo^s)))
1559 else
1560 tmp(ixo^s)=fl%tc_k_para*dsqrt((te(ixo^s))**5)/rho(ixo^s)
1561 end if
1562 end if
1563
1564 dtnew = bigdouble
1565 do idim=1,ndim
1566 ! approximate thermal conduction flux: tc_k_para/rho/dx**2
1567 maxtmp2=maxval(tmp(ixo^s)/(block%ds(ixo^s,idim)**2))
1568 ! dt< dx_idim**2/((gamma-1)*tc_k_para/rho)
1569 dtdiff_tcond=1.d0/(tc_gamma_1*maxtmp2+smalldouble)
1570 ! limit the time step
1571 dtnew=min(dtnew,dtdiff_tcond)
1572 end do
1573 dtnew=dtnew/dble(ndim)
1574 end function get_tc_dt_hd
1575
1576 subroutine sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,fl)
1579
1580 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
1581 double precision, intent(in) :: x(ixi^s,1:ndim)
1582 double precision, intent(in) :: w(ixi^s,1:nw)
1583 double precision, intent(inout) :: wres(ixi^s,1:nw)
1584 double precision, intent(in) :: my_dt
1585 logical, intent(in) :: fix_conserve_at_step
1586 type(tc_fluid), intent(in) :: fl
1587
1588 double precision :: te(ixi^s),rho(ixi^s)
1589 double precision :: qvec(ixi^s,1:ndim),qd(ixi^s)
1590 double precision, allocatable, dimension(:^D&,:) :: qvec_equi
1591 double precision :: fluxall(ixi^s,1,1:ndim)
1592
1593 double precision :: dxinv(ndim)
1594 integer :: idims,ixa^l
1595
1596 dxinv=1.d0/dxlevel
1597
1598 !calculate Te in whole domain (+ghosts)
1599 call fl%get_temperature_from_eint(w, x, ixi^l, ixi^l, te)
1600 call fl%get_rho(w, x, ixi^l, ixi^l, rho)
1601 if(slab_uniform) then
1602 call set_source_tc_hd(ixi^l,ixo^l,w,x,fl,qvec,rho,te)
1603 else
1604 call set_source_tc_hd_geo(ixi^l,ixo^l,w,x,fl,qvec,rho,te)
1605 end if
1606 if(fl%has_equi) then
1607 allocate(qvec_equi(ixi^s,1:ndim))
1608 call fl%get_temperature_equi(w, x, ixi^l, ixi^l, te) !calculate Te in whole domain (+ghosts)
1609 call fl%get_rho_equi(w, x, ixi^l, ixi^l, rho) !calculate rho in whole domain (+ghosts)
1610 if(slab_uniform) then
1611 call set_source_tc_hd(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te)
1612 else
1613 call set_source_tc_hd_geo(ixi^l,ixo^l,w,x,fl,qvec_equi,rho,te)
1614 end if
1615 do idims=1,ndim
1616 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
1617 qvec(ixa^s,idims)=qvec(ixa^s,idims)-qvec_equi(ixa^s,idims)
1618 end do
1619 deallocate(qvec_equi)
1620 endif
1621
1622 if(slab_uniform) then
1623 do idims=1,ndim
1624 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
1625 qvec(ixa^s,idims)=dxinv(idims)*qvec(ixa^s,idims)
1626 ixa^l=ixo^l-kr(idims,^d);
1627 if(idims==1) then
1628 qd(ixo^s)=qvec(ixo^s,idims)-qvec(ixa^s,idims)
1629 else
1630 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixa^s,idims)
1631 end if
1632 end do
1633 else
1634 do idims=1,ndim
1635 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
1636 qvec(ixa^s,idims)=qvec(ixa^s,idims)*block%surfaceC(ixa^s,idims)
1637 ixa^l=ixo^l-kr(idims,^d);
1638 if(idims==1) then
1639 qd(ixo^s)=qvec(ixo^s,idims)-qvec(ixa^s,idims)
1640 else
1641 qd(ixo^s)=qd(ixo^s)+qvec(ixo^s,idims)-qvec(ixa^s,idims)
1642 end if
1643 end do
1644 qd(ixo^s)=qd(ixo^s)/block%dvolume(ixo^s)
1645 end if
1646
1647 if(fix_conserve_at_step) then
1648 do idims=1,ndim
1649 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
1650 fluxall(ixa^s,1,idims)=my_dt*qvec(ixa^s,idims)
1651 end do
1652 call store_flux(igrid,fluxall,1,ndim,nflux)
1653 end if
1654
1655 wres(ixo^s,fl%e_)=qd(ixo^s)
1656 end subroutine sts_set_source_tc_hd
1657
1658 !> get HD thermal conduction source and flux in uniform Cartesian grid
1659 subroutine set_source_tc_hd(ixI^L,ixO^L,w,x,fl,qvec,rho,Te)
1661 integer, intent(in) :: ixI^L, ixO^L
1662 double precision, intent(in) :: x(ixI^S,1:ndim)
1663 double precision, intent(in) :: w(ixI^S,1:nw)
1664 type(tc_fluid), intent(in) :: fl
1665 double precision, intent(in) :: Te(ixI^S),rho(ixI^S)
1666 double precision, intent(out) :: qvec(ixI^S,1:ndim)
1667 double precision :: gradT(ixI^S,1:ndim),ke(ixI^S),qd(ixI^S)
1668 integer :: idims,ix^D,ix^L,ixC^L,ixA^L,ixB^L
1669
1670 ix^l=ixo^l^ladd1;
1671 ! ixC is cell-corner index
1672 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
1673
1674 ! calculate thermal conduction flux with symmetric scheme
1675 ! T gradient (central difference) at cell corners
1676 do idims=1,ndim
1677 ixbmin^d=ixmin^d;
1678 ixbmax^d=ixmax^d-kr(idims,^d);
1679 call gradientf(te,x,ixi^l,ixb^l,idims,ke)
1680 {^ifthreed
1681 if(idims==1) then
1682 {do ix^db=ixcmin^db,ixcmax^db\}
1683 qvec(ix^d,idims)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1,ix2+1,ix3)&
1684 +ke(ix1,ix2,ix3+1)+ke(ix1,ix2+1,ix3+1))
1685 {end do\}
1686 else if(idims==2) then
1687 {do ix^db=ixcmin^db,ixcmax^db\}
1688 qvec(ix^d,idims)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1+1,ix2,ix3)&
1689 +ke(ix1,ix2,ix3+1)+ke(ix1+1,ix2,ix3+1))
1690 {end do\}
1691 else
1692 {do ix^db=ixcmin^db,ixcmax^db\}
1693 qvec(ix^d,idims)=0.25d0*(ke(ix1,ix2,ix3)+ke(ix1+1,ix2,ix3)&
1694 +ke(ix1,ix2+1,ix3)+ke(ix1+1,ix2+1,ix3))
1695 {end do\}
1696 end if
1697 }
1698 {^iftwod
1699 if(idims==1) then
1700 {do ix^db=ixcmin^db,ixcmax^db\}
1701 qvec(ix^d,idims)=0.5d0*(ke(ix1,ix2)+ke(ix1,ix2+1))
1702 {end do\}
1703 else
1704 {do ix^db=ixcmin^db,ixcmax^db\}
1705 qvec(ix^d,idims)=0.5d0*(ke(ix1,ix2)+ke(ix1+1,ix2))
1706 {end do\}
1707 end if
1708 }
1709 {^ifoned
1710 do ix1=ixcmin1,ixcmax1
1711 qvec(ix1,idims)=ke(ix1)
1712 end do
1713 }
1714 end do
1715 ! conductivity at cell center
1716 if(fl%tc_constant) then
1717 qd(ix^s)=fl%tc_k_para
1718 else
1719 if(phys_trac) then
1720 {do ix^db=ixmin^db,ixmax^db\}
1721 if(te(ix^d) < block%wextra(ix^d,fl%Tcoff_)) then
1722 qd(ix^d)=fl%tc_k_para*dsqrt(block%wextra(ix^d,fl%Tcoff_)**5)
1723 else
1724 qd(ix^d)=fl%tc_k_para*dsqrt(te(ix^d)**5)
1725 end if
1726 {end do\}
1727 else
1728 qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s)**5)
1729 end if
1730 end if
1731 ! conductivity Ke at cell corner
1732 ! cell corner conduction flux gradT
1733 {^ifthreed
1734 {do ix^db=ixcmin^db,ixcmax^db\}
1735 ke(ix^d)=0.125d0*(qd(ix1,ix2,ix3)+qd(ix1+1,ix2,ix3)&
1736 +qd(ix1,ix2+1,ix3)+qd(ix1+1,ix2+1,ix3)&
1737 +qd(ix1,ix2,ix3+1)+qd(ix1+1,ix2,ix3+1)&
1738 +qd(ix1,ix2+1,ix3+1)+qd(ix1+1,ix2+1,ix3+1))
1739 gradt(ix^d,1)=ke(ix^d)*qvec(ix^d,1)
1740 gradt(ix^d,2)=ke(ix^d)*qvec(ix^d,2)
1741 gradt(ix^d,3)=ke(ix^d)*qvec(ix^d,3)
1742 {end do\}
1743 }
1744 {^iftwod
1745 {do ix^db=ixcmin^db,ixcmax^db\}
1746 ke(ix^d)=0.25d0*(qd(ix1,ix2)+qd(ix1+1,ix2)+qd(ix1,ix2+1)+qd(ix1+1,ix2+1))
1747 gradt(ix^d,1)=ke(ix^d)*qvec(ix^d,1)
1748 gradt(ix^d,2)=ke(ix^d)*qvec(ix^d,2)
1749 {end do\}
1750 }
1751 {^ifoned
1752 do ix1=ixcmin1,ixcmax1
1753 gradt(ix^d,1)=0.5d0*(qd(ix1)+qd(ix1+1))*qvec(ix^d,1)
1754 end do
1755 }
1756
1757 ! conduction flux qvec at cell face
1758 do idims=1,ndim
1759 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
1760 {^ifthreed
1761 if(idims==1) then
1762 {do ix^db=ixamin^db,ixamax^db\}
1763 qvec(ix^d,idims)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1,ix2-1,ix3,idims)&
1764 +gradt(ix1,ix2,ix3-1,idims)+gradt(ix1,ix2-1,ix3-1,idims))
1765 {end do\}
1766 else if(idims==2) then
1767 {do ix^db=ixamin^db,ixamax^db\}
1768 qvec(ix^d,idims)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1-1,ix2,ix3,idims)&
1769 +gradt(ix1,ix2,ix3-1,idims)+gradt(ix1-1,ix2,ix3-1,idims))
1770 {end do\}
1771 else
1772 {do ix^db=ixamin^db,ixamax^db\}
1773 qvec(ix^d,idims)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1,ix2-1,ix3,idims)&
1774 +gradt(ix1-1,ix2,ix3,idims)+gradt(ix1-1,ix2-1,ix3,idims))
1775 {end do\}
1776 end if
1777 }
1778 {^iftwod
1779 if(idims==1) then
1780 {do ix^db=ixamin^db,ixamax^db\}
1781 qvec(ix^d,idims)=0.5d0*(gradt(ix1,ix2,idims)+gradt(ix1,ix2-1,idims))
1782 {end do\}
1783 else
1784 {do ix^db=ixamin^db,ixamax^db\}
1785 qvec(ix^d,idims)=0.5d0*(gradt(ix1,ix2,idims)+gradt(ix1-1,ix2,idims))
1786 {end do\}
1787 end if
1788 }
1789 {^ifoned
1790 do ix1=ixamin1,ixamax1
1791 qvec(ix1,idims)=gradt(ix1,idims)
1792 end do
1793 }
1794 if(fl%tc_saturate) then
1795 ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
1796 ! unsigned saturated TC flux = 5 phi rho c**3, c=sqrt(p/rho) is isothermal sound speed, phi=1.1
1797 ixb^l=ixa^l+kr(idims,^d);
1798 qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3
1799 {do ix^db=ixamin^db,ixamax^db\}
1800 if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
1801 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
1802 end if
1803 {end do\}
1804 end if
1805 end do
1806 end subroutine set_source_tc_hd
1807
1808 !> get HD thermal conduction source and flux in non-uniform geometry grid
1809 subroutine set_source_tc_hd_geo(ixI^L,ixO^L,w,x,fl,qvec,rho,Te)
1811 integer, intent(in) :: ixI^L, ixO^L
1812 double precision, intent(in) :: x(ixI^S,1:ndim)
1813 double precision, intent(in) :: w(ixI^S,1:nw)
1814 type(tc_fluid), intent(in) :: fl
1815 double precision, intent(in) :: Te(ixI^S),rho(ixI^S)
1816 double precision, intent(out) :: qvec(ixI^S,1:ndim)
1817 double precision :: gradT(ixI^S,1:ndim),ke(ixI^S),qd(ixI^S)
1818 integer :: idims,ix^D,ix^L,ixC^L,ixA^L,ixB^L
1819
1820 ix^l=ixo^l^ladd1;
1821 ! ixC is cell-corner index
1822 ixcmax^d=ixomax^d; ixcmin^d=ixomin^d-1;
1823
1824 ! calculate thermal conduction flux with symmetric scheme
1825 ! T gradient from face centers to cell corners: surface weighted average
1826 do idims=1,ndim
1827 ixbmin^d=ixmin^d;
1828 ixbmax^d=ixmax^d-kr(idims,^d);
1829 call gradientf(te,x,ixi^l,ixb^l,idims,ke)
1830 {^ifthreed
1831 if(idims==1) then
1832 {do ix^db=ixcmin^db,ixcmax^db\}
1833 qvec(ix^d,1)=(ke(ix1,ix2,ix3)*block%surfaceC(ix1,ix2,ix3,1)&
1834 +ke(ix1,ix2+1,ix3)*block%surfaceC(ix1,ix2+1,ix3,1)&
1835 +ke(ix1,ix2,ix3+1)*block%surfaceC(ix1,ix2,ix3+1,1)&
1836 +ke(ix1,ix2+1,ix3+1)*block%surfaceC(ix1,ix2+1,ix3+1,1))/&
1837 (block%surfaceC(ix1,ix2,ix3,1)+block%surfaceC(ix1,ix2+1,ix3,1)&
1838 +block%surfaceC(ix1,ix2,ix3+1,1)+block%surfaceC(ix1,ix2+1,ix3+1,1))
1839 {end do\}
1840 else if(idims==2) then
1841 {do ix^db=ixcmin^db,ixcmax^db\}
1842 qvec(ix^d,2)=(ke(ix1,ix2,ix3)*block%surfaceC(ix1,ix2,ix3,2)&
1843 +ke(ix1+1,ix2,ix3)*block%surfaceC(ix1+1,ix2,ix3,2)&
1844 +ke(ix1,ix2,ix3+1)*block%surfaceC(ix1,ix2,ix3+1,2)&
1845 +ke(ix1+1,ix2,ix3+1)*block%surfaceC(ix1+1,ix2,ix3+1,2))/&
1846 (block%surfaceC(ix1,ix2,ix3,2)+block%surfaceC(ix1+1,ix2,ix3,2)&
1847 +block%surfaceC(ix1,ix2,ix3+1,2)+block%surfaceC(ix1+1,ix2,ix3+1,2)+1.d-300)
1848 ! zero theta-normal surface area at pole axis
1849 {end do\}
1850 else
1851 {do ix^db=ixcmin^db,ixcmax^db\}
1852 qvec(ix^d,3)=(ke(ix1,ix2,ix3)*block%surfaceC(ix1,ix2,ix3,3)&
1853 +ke(ix1+1,ix2,ix3)*block%surfaceC(ix1+1,ix2,ix3,3)&
1854 +ke(ix1,ix2+1,ix3)*block%surfaceC(ix1,ix2+1,ix3,3)&
1855 +ke(ix1+1,ix2+1,ix3)*block%surfaceC(ix1+1,ix2+1,ix3,3))/&
1856 (block%surfaceC(ix1,ix2,ix3,3)+block%surfaceC(ix1+1,ix2,ix3,3)&
1857 +block%surfaceC(ix1,ix2+1,ix3,3)+block%surfaceC(ix1+1,ix2+1,ix3,3))
1858 {end do\}
1859 end if
1860 }
1861 {^iftwod
1862 if(idims==1) then
1863 {do ix^db=ixcmin^db,ixcmax^db\}
1864 qvec(ix^d,1)=(ke(ix1,ix2)*block%surfaceC(ix1,ix2,1)+ke(ix1,ix2+1)*block%surfaceC(ix1,ix2+1,1))&
1865 /(block%surfaceC(ix1,ix2,1)+block%surfaceC(ix1,ix2+1,1))
1866 {end do\}
1867 else
1868 {do ix^db=ixcmin^db,ixcmax^db\}
1869 qvec(ix^d,2)=(ke(ix1,ix2)*block%surfaceC(ix1,ix2,2)+ke(ix1+1,ix2)*block%surfaceC(ix1+1,ix2,2))&
1870 /(block%surfaceC(ix1,ix2,2)+block%surfaceC(ix1+1,ix2,2))
1871 {end do\}
1872 end if
1873 }
1874 {^ifoned
1875 do ix1=ixcmin1,ixcmax1
1876 qvec(ix1,idims)=ke(ix1)
1877 end do
1878 }
1879 end do
1880 ! conductivity at cell center
1881 if(fl%tc_constant) then
1882 qd(ix^s)=fl%tc_k_para
1883 else
1884 if(phys_trac) then
1885 {do ix^db=ixmin^db,ixmax^db\}
1886 if(te(ix^d) < block%wextra(ix^d,fl%Tcoff_)) then
1887 qd(ix^d)=fl%tc_k_para*dsqrt(block%wextra(ix^d,fl%Tcoff_)**5)
1888 else
1889 qd(ix^d)=fl%tc_k_para*dsqrt(te(ix^d)**5)
1890 end if
1891 {end do\}
1892 else
1893 qd(ix^s)=fl%tc_k_para*dsqrt(te(ix^s)**5)
1894 end if
1895 end if
1896 ! conductivity Ke at cell corner
1897 ! cell corner conduction flux gradT
1898 {^ifthreed
1899 {do ix^db=ixcmin^db,ixcmax^db\}
1900 ke(ix^d)=(qd(ix1,ix2,ix3)*block%dvolume(ix1,ix2,ix3)+qd(ix1+1,ix2,ix3)*block%dvolume(ix1+1,ix2,ix3)&
1901 +qd(ix1,ix2+1,ix3)*block%dvolume(ix1,ix2+1,ix3)+qd(ix1+1,ix2+1,ix3)*block%dvolume(ix1+1,ix2+1,ix3)&
1902 +qd(ix1,ix2,ix3+1)*block%dvolume(ix1,ix2,ix3+1)+qd(ix1+1,ix2,ix3+1)*block%dvolume(ix1+1,ix2,ix3+1)&
1903 +qd(ix1,ix2+1,ix3+1)*block%dvolume(ix1,ix2+1,ix3+1)+qd(ix1+1,ix2+1,ix3+1)*block%dvolume(ix1+1,ix2+1,ix3+1))&
1904 /(block%dvolume(ix1,ix2,ix3)+block%dvolume(ix1+1,ix2,ix3)+block%dvolume(ix1,ix2+1,ix3)+block%dvolume(ix1+1,ix2+1,ix3)&
1905 +block%dvolume(ix1,ix2,ix3+1)+block%dvolume(ix1+1,ix2,ix3+1)+block%dvolume(ix1,ix2+1,ix3+1)+block%dvolume(ix1+1,ix2+1,ix3+1))
1906 gradt(ix^d,1)=ke(ix^d)*qvec(ix^d,1)
1907 gradt(ix^d,2)=ke(ix^d)*qvec(ix^d,2)
1908 gradt(ix^d,3)=ke(ix^d)*qvec(ix^d,3)
1909 {end do\}
1910 }
1911 {^iftwod
1912 {do ix^db=ixcmin^db,ixcmax^db\}
1913 ke(ix^d)=(qd(ix1,ix2)*block%dvolume(ix1,ix2)+qd(ix1+1,ix2)*block%dvolume(ix1+1,ix2)&
1914 +qd(ix1,ix2+1)*block%dvolume(ix1,ix2+1)+qd(ix1+1,ix2+1)*block%dvolume(ix1+1,ix2+1))&
1915 /(block%dvolume(ix1,ix2)+block%dvolume(ix1+1,ix2)+block%dvolume(ix1,ix2+1)+block%dvolume(ix1+1,ix2+1))
1916 gradt(ix^d,1)=ke(ix^d)*qvec(ix^d,1)
1917 gradt(ix^d,2)=ke(ix^d)*qvec(ix^d,2)
1918 {end do\}
1919 }
1920 {^ifoned
1921 do ix1=ixcmin1,ixcmax1
1922 gradt(ix^d,1)=(qd(ix1)*block%dvolume(ix1)+qd(ix1+1)*block%dvolume(ix1+1))/(block%dvolume(ix1)+block%dvolume(ix1+1))*qvec(ix^d,1)
1923 end do
1924 }
1925
1926 ! conduction flux qvec at cell face
1927 do idims=1,ndim
1928 ixamax^d=ixomax^d; ixamin^d=ixomin^d-kr(idims,^d);
1929 {^ifthreed
1930 if(idims==1) then
1931 {do ix^db=ixamin^db,ixamax^db\}
1932 qvec(ix^d,idims)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1,ix2-1,ix3,idims)&
1933 +gradt(ix1,ix2,ix3-1,idims)+gradt(ix1,ix2-1,ix3-1,idims))
1934 {end do\}
1935 else if(idims==2) then
1936 {do ix^db=ixamin^db,ixamax^db\}
1937 qvec(ix^d,idims)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1-1,ix2,ix3,idims)&
1938 +gradt(ix1,ix2,ix3-1,idims)+gradt(ix1-1,ix2,ix3-1,idims))
1939 {end do\}
1940 else
1941 {do ix^db=ixamin^db,ixamax^db\}
1942 qvec(ix^d,idims)=0.25d0*(gradt(ix1,ix2,ix3,idims)+gradt(ix1,ix2-1,ix3,idims)&
1943 +gradt(ix1-1,ix2,ix3,idims)+gradt(ix1-1,ix2-1,ix3,idims))
1944 {end do\}
1945 end if
1946 }
1947 {^iftwod
1948 if(idims==1) then
1949 {do ix^db=ixamin^db,ixamax^db\}
1950 qvec(ix^d,idims)=0.5d0*(gradt(ix1,ix2,idims)+gradt(ix1,ix2-1,idims))
1951 {end do\}
1952 else
1953 {do ix^db=ixamin^db,ixamax^db\}
1954 qvec(ix^d,idims)=0.5d0*(gradt(ix1,ix2,idims)+gradt(ix1-1,ix2,idims))
1955 {end do\}
1956 end if
1957 }
1958 {^ifoned
1959 do ix1=ixamin1,ixamax1
1960 qvec(ix1,idims)=gradt(ix1,idims)
1961 end do
1962 }
1963 if(fl%tc_saturate) then
1964 ! consider saturation (Cowie and Mckee 1977 ApJ, 211, 135)
1965 ! unsigned saturated TC flux = 5 phi rho c**3, c=sqrt(p/rho) is isothermal sound speed, phi=1.1
1966 ixb^l=ixa^l+kr(idims,^d);
1967 qd(ixa^s)=2.75d0*(rho(ixa^s)+rho(ixb^s))*dsqrt(0.5d0*(te(ixa^s)+te(ixb^s)))**3
1968 {do ix^db=ixamin^db,ixamax^db\}
1969 if(dabs(qvec(ix^d,idims))>qd(ix^d)) then
1970 qvec(ix^d,idims)=sign(1.d0,qvec(ix^d,idims))*qd(ix^d)
1971 end if
1972 {end do\}
1973 end if
1974 end do
1975 end subroutine set_source_tc_hd_geo
1976
1977end module mod_thermal_conduction
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
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
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.
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.
double precision, dimension(:), allocatable, parameter d
logical b0field
split magnetic field as background B0 field
double precision unit_temperature
Physical scaling factor for temperature.
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
Thermal conduction for HD and MHD or RHD and RMHD or twofl (plasma-neutral) module Adaptation of mod_...
subroutine, public tc_get_hd_params(fl, read_hd_params)
Init TC coefficients: HD case.
double precision function, public get_tc_dt_mhd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (mhd implementation) Note: for multi-D MHD (1D MHD will use HD f...
double precision function, public get_tc_dt_hd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (hd implementation) Note: also used in 1D MHD (or for neutrals i...
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_hd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
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...
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_hd_geo(ixil, ixol, w, x, fl, qvec, rho, te)
get HD thermal conduction source and flux in non-uniform geometry grid
double precision function, dimension(ixi^s) slope_limiter(f, ixil, ixol, idims, pm, tc_slope_limiter)
subroutine set_source_tc_hd(ixil, ixol, w, x, fl, qvec, rho, te)
get HD thermal conduction source and flux in uniform Cartesian grid
subroutine set_source_tc_mhd_geo(ixil, ixol, w, x, fl, qvec, rho, te, alpha)
subroutine set_source_tc_mhd(ixil, ixol, w, x, fl, qvec, rho, te, alpha)