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