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