MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_ffhd_phys.t
Go to the documentation of this file.
1!> Frozen-field hydrodynamics module
3
4#include "amrvac.h"
5
6 use mod_global_parameters, only: std_len, const_c
10 use mod_physics
11 use mod_comm_lib, only: mpistop
12
13 implicit none
14 private
15
16 !> Whether an energy equation is used
17 logical, public, protected :: ffhd_energy = .true.
18
19 !> Whether thermal conduction is used
20 logical, public, protected :: ffhd_thermal_conduction = .false.
21 !> Whether hyperbolic type thermal conduction is used
22 logical, public, protected :: ffhd_hyperbolic_thermal_conduction = .false.
23 !> Wheterh saturation is considered for hyperbolic TC
24 logical, public, protected :: ffhd_htc_sat = .false.
25 !> type of fluid for thermal conduction
26 type(tc_fluid), public, allocatable :: tc_fl
27 !> type of fluid for thermal emission synthesis
28 type(te_fluid), public, allocatable :: te_fl_ffhd
29
30 !> Whether radiative cooling is added
31 logical, public, protected :: ffhd_radiative_cooling = .false.
32 !> type of fluid for radiative cooling
33 type(rc_fluid), public, allocatable :: rc_fl
34
35 !> Whether gravity is added
36 logical, public, protected :: ffhd_gravity = .false.
37
38 !> Whether TRAC method is used
39 logical, public, protected :: ffhd_trac = .false.
40
41 !> Which TRAC method is used
42 integer, public, protected :: ffhd_trac_type=1
43
44 !> Height of the mask used in the TRAC method
45 double precision, public, protected :: ffhd_trac_mask = 0.d0
46
47 !> Distance between two adjacent traced magnetic field lines (in finest cell size)
48 integer, public, protected :: ffhd_trac_finegrid=4
49
50 !> Whether plasma is partially ionized
51 logical, public, protected :: ffhd_partial_ionization = .false.
52
53 !> Index of the density (in the w array)
54 integer, public, protected :: rho_
55
56 !> Indices of the momentum density
57 integer, allocatable, public, protected :: mom(:)
58
59 !> Index of the energy density (-1 if not present)
60 integer, public, protected :: e_
61
62 !> Index of the gas pressure (-1 if not present) should equal e_
63 integer, public, protected :: p_
64
65 !> Indices of temperature
66 integer, public, protected :: te_
67
68 !> Index of the cutoff temperature for the TRAC method
69 integer, public, protected :: tcoff_
70 integer, public, protected :: tweight_
71 integer, public, protected :: q_
72
73 !> The adiabatic index
74 double precision, public :: ffhd_gamma = 5.d0/3.0d0
75
76 !> The adiabatic constant
77 double precision, public :: ffhd_adiab = 1.0d0
78
79 !> The small_est allowed energy
80 double precision, protected :: small_e
81
82 !> The thermal conductivity kappa in hyperbolic thermal conduction
83 double precision, public :: hypertc_kappa
84
85 !> Helium abundance over Hydrogen
86 double precision, public, protected :: he_abundance=0.1d0
87 !> Ionization fraction of H
88 !> H_ion_fr = H+/(H+ + H)
89 double precision, public, protected :: h_ion_fr=1d0
90 !> Ionization fraction of He
91 !> He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
92 double precision, public, protected :: he_ion_fr=1d0
93 !> Ratio of number He2+ / number He+ + He2+
94 !> He_ion_fr2 = He2+/(He2+ + He+)
95 double precision, public, protected :: he_ion_fr2=1d0
96 ! used for eq of state when it is not defined by units,
97 ! the units do not contain terms related to ionization fraction
98 ! and it is p = RR * rho * T
99 double precision, public, protected :: rr=1d0
100 ! remove the below flag and assume default value = .false.
101 ! when eq state properly implemented everywhere
102 ! and not anymore through units
103 logical, public, protected :: eq_state_units = .true.
104
105 !> gamma minus one and its inverse
106 double precision :: gamma_1, inv_gamma_1
107
108 !define the function interface for the kinetic energy
109 abstract interface
110
111 function fun_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
112 use mod_global_parameters, only: nw, ndim,block
113 integer, intent(in) :: ixi^l, ixo^l
114 double precision, intent(in) :: w(ixi^s, nw)
115 double precision :: ke(ixo^s)
116 double precision, intent(in), optional :: inv_rho(ixo^s)
117 end function fun_kin_en
118
119 end interface
120
121 procedure(sub_convert), pointer :: ffhd_to_primitive => null()
122 procedure(sub_convert), pointer :: ffhd_to_conserved => null()
123 procedure(sub_small_values), pointer :: ffhd_handle_small_values => null()
124 procedure(sub_get_pthermal), pointer :: ffhd_get_pthermal => null()
125 procedure(sub_get_pthermal), pointer :: ffhd_get_rfactor => null()
126 procedure(sub_get_pthermal), pointer :: ffhd_get_temperature => null()
127 procedure(sub_get_v), pointer :: ffhd_get_v => null()
128 procedure(fun_kin_en), pointer :: ffhd_kin_en => null()
129 ! Public methods
130 public :: ffhd_phys_init
131 public :: ffhd_kin_en
132 public :: ffhd_get_pthermal
133 public :: ffhd_get_temperature
134 public :: ffhd_get_v
135 public :: ffhd_get_rho
136 public :: ffhd_get_v_idim
137 public :: ffhd_to_conserved
138 public :: ffhd_to_primitive
139 public :: ffhd_get_csound2
140 public :: ffhd_e_to_ei
141 public :: ffhd_ei_to_e
142
143contains
144
145 subroutine ffhd_read_params(files)
147 use mod_particles, only: particles_eta, particles_etah
148 character(len=*), intent(in) :: files(:)
149 integer :: n
150
151 namelist /ffhd_list/ ffhd_energy, ffhd_gamma, ffhd_adiab, &
155
156 do n = 1, size(files)
157 open(unitpar, file=trim(files(n)), status="old")
158 read(unitpar, ffhd_list, end=111)
159111 close(unitpar)
160 end do
161 end subroutine ffhd_read_params
162
163 !> Write this module's parameters to a snapsoht
164 subroutine ffhd_write_info(fh)
166 integer, intent(in) :: fh
167 integer, parameter :: n_par = 1
168 double precision :: values(n_par)
169 character(len=name_len) :: names(n_par)
170 integer, dimension(MPI_STATUS_SIZE) :: st
171 integer :: er
172
173 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
174
175 names(1) = "gamma"
176 values(1) = ffhd_gamma
177 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
178 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
179 end subroutine ffhd_write_info
180
181 subroutine ffhd_phys_init()
185 use mod_gravity, only: gravity_init
190 integer :: itr, idir
191
192 call ffhd_read_params(par_files)
193
194 if(.not. ffhd_energy) then
197 if(mype==0) write(*,*) 'WARNING: set ffhd_thermal_conduction=F when ffhd_energy=F'
198 end if
201 if(mype==0) write(*,*) 'WARNING: set ffhd_hyperbolic_thermal_conduction=F when ffhd_energy=F'
202 end if
205 if(mype==0) write(*,*) 'WARNING: set ffhd_radiative_cooling=F when ffhd_energy=F'
206 end if
207 if(ffhd_trac) then
208 ffhd_trac=.false.
209 if(mype==0) write(*,*) 'WARNING: set ffhd_trac=F when ffhd_energy=F'
210 end if
213 if(mype==0) write(*,*) 'WARNING: set ffhd_partial_ionization=F when ffhd_energy=F'
214 end if
215 end if
216 if(.not.eq_state_units) then
219 if(mype==0) write(*,*) 'WARNING: set ffhd_partial_ionization=F when eq_state_units=F'
220 end if
221 end if
222
225 if(mype==0) write(*,*) 'WARNING: turn off parabolic TC when using hyperbolic TC'
226 end if
227
228 physics_type = "ffhd"
229 phys_energy=ffhd_energy
230 phys_internal_e=.false.
233 phys_partial_ionization=ffhd_partial_ionization
234
235 phys_gamma = ffhd_gamma
236 phys_total_energy=ffhd_energy
238
239 {^ifoned
240 if(ffhd_trac .and. ffhd_trac_type .gt. 2) then
242 if(mype==0) write(*,*) 'WARNING: reset ffhd_trac_type=1 for 1D simulation'
243 end if
244 }
245 if(ffhd_trac .and. ffhd_trac_type .le. 4) then
246 ffhd_trac_mask=bigdouble
247 if(mype==0) write(*,*) 'WARNING: set ffhd_trac_mask==bigdouble for global TRAC method'
248 end if
250
251 allocate(start_indices(number_species),stop_indices(number_species))
252 start_indices(1)=1
253 ! Determine flux variables
254 rho_ = var_set_rho()
255
256 allocate(mom(1))
257 mom(:) = var_set_momentum(1)
258
259 ! Set index of energy variable
260 if(ffhd_energy) then
261 e_ = var_set_energy() ! energy density
262 p_ = e_ ! gas pressure
263 else
264 e_ = -1
265 p_ = -1
266 end if
267
269 q_ = var_set_q()
270 need_global_cmax=.true.
271 else
272 q_=-1
273 end if
274
275 ! set number of variables which need update ghostcells
276 nwgc=nwflux
277
278 ! set the index of the last flux variable for species 1
279 stop_indices(1)=nwflux
280
281 ! set temperature as an auxiliary variable to get ionization degree
283 te_ = var_set_auxvar('Te','Te')
284 else
285 te_ = -1
286 end if
287
288 ! set cutoff temperature when using the TRAC method, as well as an auxiliary weight
289 tweight_ = -1
290 if(ffhd_trac) then
291 tcoff_ = var_set_wextra()
292 iw_tcoff=tcoff_
293 if(ffhd_trac_type .ge. 3) then
294 tweight_ = var_set_wextra()
295 iw_tweight=tweight_
296 end if
297 else
298 tcoff_ = -1
299 end if
300
301 nvector = 0 ! No. vector vars
302
303 ! Check whether custom flux types have been defined
304 if(.not. allocated(flux_type)) then
305 allocate(flux_type(ndir, nwflux))
306 flux_type = flux_default
307 else if(any(shape(flux_type) /= [ndir, nwflux])) then
308 call mpistop("phys_check error: flux_type has wrong shape")
309 end if
310
311 phys_get_dt => ffhd_get_dt
312 phys_get_cmax => ffhd_get_cmax_origin
313 phys_get_tcutoff => ffhd_get_tcutoff
314 phys_get_cbounds => ffhd_get_cbounds
315 phys_to_primitive => ffhd_to_primitive_origin
316 ffhd_to_primitive => ffhd_to_primitive_origin
317 phys_to_conserved => ffhd_to_conserved_origin
318 ffhd_to_conserved => ffhd_to_conserved_origin
319 phys_get_flux => ffhd_get_flux
320 phys_get_v => ffhd_get_v_origin
321 ffhd_get_v => ffhd_get_v_origin
322 phys_get_rho => ffhd_get_rho
323 ffhd_kin_en => ffhd_kin_en_origin
324 phys_add_source_geom => ffhd_add_source_geom
325 phys_add_source => ffhd_add_source
326 phys_check_params => ffhd_check_params
327 phys_write_info => ffhd_write_info
328 phys_handle_small_values => ffhd_handle_small_values_origin
329 ffhd_handle_small_values => ffhd_handle_small_values_origin
330 phys_check_w => ffhd_check_w_origin
331
332 if(.not.ffhd_energy) then
333 phys_get_pthermal => ffhd_get_pthermal_iso
334 ffhd_get_pthermal => ffhd_get_pthermal_iso
335 else
336 phys_get_pthermal => ffhd_get_pthermal_origin
337 ffhd_get_pthermal => ffhd_get_pthermal_origin
338 end if
339
340 ! choose Rfactor in ideal gas law
342 ffhd_get_rfactor=>rfactor_from_temperature_ionization
343 phys_update_temperature => ffhd_update_temperature
344 else if(associated(usr_rfactor)) then
345 ffhd_get_rfactor=>usr_rfactor
346 else
347 ffhd_get_rfactor=>rfactor_from_constant_ionization
348 end if
349
351 ffhd_get_temperature => ffhd_get_temperature_from_te
352 else
353 ffhd_get_temperature => ffhd_get_temperature_from_etot
354 end if
355
356 ! derive units from basic units
357 call ffhd_physical_units()
358
360 if(si_unit)then
361 ! parallel conduction Spitzer
363 else
364 ! in cgs
366 endif
367 end if
368 if(.not. ffhd_energy .and. ffhd_thermal_conduction) then
369 call mpistop("thermal conduction needs ffhd_energy=T")
370 end if
372 call mpistop("hyperbolic thermal conduction needs ffhd_energy=T")
373 end if
374 if(.not. ffhd_energy .and. ffhd_radiative_cooling) then
375 call mpistop("radiative cooling needs ffhd_energy=T")
376 end if
377
378 ! initialize thermal conduction module
380 call sts_init()
382
383 allocate(tc_fl)
384 call tc_get_hd_params(tc_fl,tc_params_read_ffhd)
385 call add_sts_method(ffhd_get_tc_dt_ffhd,ffhd_sts_set_source_tc_ffhd,e_,1,e_,1,.false.)
386 tc_fl%get_temperature_from_conserved => ffhd_get_temperature_from_etot
387 tc_fl%get_temperature_from_eint => ffhd_get_temperature_from_eint
389 call set_error_handling_to_head(ffhd_tc_handle_small_e)
390 tc_fl%get_rho => ffhd_get_rho
391 tc_fl%e_ = e_
392 tc_fl%Tcoff_ = tcoff_
393 end if
394
395 ! Initialize radiative cooling module
398 allocate(rc_fl)
399 call radiative_cooling_init(rc_fl,rc_params_read)
400 rc_fl%get_rho => ffhd_get_rho
401 rc_fl%get_pthermal => ffhd_get_pthermal
402 rc_fl%get_var_Rfactor => ffhd_get_rfactor
403 rc_fl%e_ = e_
404 rc_fl%Tcoff_ = tcoff_
405 rc_fl%has_equi = .false.
406 rc_fl%subtract_equi = .false.
407 end if
408
409{^ifthreed
410 allocate(te_fl_ffhd)
411 te_fl_ffhd%get_rho=> ffhd_get_rho
412 te_fl_ffhd%get_pthermal=> ffhd_get_pthermal
413 te_fl_ffhd%get_var_Rfactor => ffhd_get_rfactor
414 phys_te_images => ffhd_te_images
415}
416
417 ! Initialize gravity module
418 if(ffhd_gravity) then
419 call gravity_init()
420 end if
421
422 ! initialize ionization degree table
424 end subroutine ffhd_phys_init
425
426{^ifthreed
427 subroutine ffhd_te_images
430
431 select case(convert_type)
432 case('EIvtiCCmpi','EIvtuCCmpi')
434 case('ESvtiCCmpi','ESvtuCCmpi')
436 case('SIvtiCCmpi','SIvtuCCmpi')
438 case('WIvtiCCmpi','WIvtuCCmpi')
440 case default
441 call mpistop("Error in synthesize emission: Unknown convert_type")
442 end select
443 end subroutine ffhd_te_images
444}
445
446 subroutine ffhd_sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
450 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
451 double precision, intent(in) :: x(ixi^s,1:ndim)
452 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
453 double precision, intent(in) :: my_dt
454 logical, intent(in) :: fix_conserve_at_step
455 call sts_set_source_tc_mhd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
456 end subroutine ffhd_sts_set_source_tc_ffhd
457
458 function ffhd_get_tc_dt_ffhd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
459 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
460 !where tc_k_para_i=tc_k_para*B_i**2/B**2
461 !and T=p/rho
464
465 integer, intent(in) :: ixi^l, ixo^l
466 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
467 double precision, intent(in) :: w(ixi^s,1:nw)
468 double precision :: dtnew
469
470 dtnew=get_tc_dt_mhd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
471 end function ffhd_get_tc_dt_ffhd
472
473 subroutine ffhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
475
476 integer, intent(in) :: ixi^l,ixo^l
477 double precision, intent(inout) :: w(ixi^s,1:nw)
478 double precision, intent(in) :: x(ixi^s,1:ndim)
479 integer, intent(in) :: step
480 character(len=140) :: error_msg
481
482 write(error_msg,"(a,i3)") "Thermal conduction step ", step
483 call ffhd_handle_small_ei(w,x,ixi^l,ixo^l,e_,error_msg)
484 end subroutine ffhd_tc_handle_small_e
485
486 subroutine tc_params_read_ffhd(fl)
488 type(tc_fluid), intent(inout) :: fl
489 integer :: n
490 ! list parameters
491 logical :: tc_saturate=.false.
492 double precision :: tc_k_para=0d0
493 character(len=std_len) :: tc_slope_limiter="MC"
494
495 namelist /tc_list/ tc_saturate, tc_slope_limiter, tc_k_para
496
497 do n = 1, size(par_files)
498 open(unitpar, file=trim(par_files(n)), status="old")
499 read(unitpar, tc_list, end=111)
500111 close(unitpar)
501 end do
502
503 fl%tc_saturate = tc_saturate
504 fl%tc_k_para = tc_k_para
505 select case(tc_slope_limiter)
506 case ('no','none')
507 fl%tc_slope_limiter = 0
508 case ('MC')
509 ! montonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
510 fl%tc_slope_limiter = 1
511 case('minmod')
512 ! minmod limiter
513 fl%tc_slope_limiter = 2
514 case ('superbee')
515 ! Roes superbee limiter (eq.3.51i)
516 fl%tc_slope_limiter = 3
517 case ('koren')
518 ! Barry Koren Right variant
519 fl%tc_slope_limiter = 4
520 case default
521 call mpistop("Unknown tc_slope_limiter, choose MC, minmod")
522 end select
523 end subroutine tc_params_read_ffhd
524
525 subroutine rc_params_read(fl)
527 use mod_constants, only: bigdouble
528 type(rc_fluid), intent(inout) :: fl
529 integer :: n
530 integer :: ncool = 4000
531
532 !> Name of cooling curve
533 character(len=std_len) :: coolcurve='JCcorona'
534
535 !> Fixed temperature not lower than tlow
536 logical :: tfix=.false.
537
538 !> Lower limit of temperature
539 double precision :: tlow=bigdouble
540
541 !> Add cooling source in a split way (.true.) or un-split way (.false.)
542 logical :: rc_split=.false.
543 logical :: rad_cut=.false.
544 double precision :: rad_cut_hgt=0.5d0
545 double precision :: rad_cut_dey=0.15d0
546
547 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split, rad_cut, rad_cut_hgt, rad_cut_dey
548
549 do n = 1, size(par_files)
550 open(unitpar, file=trim(par_files(n)), status="old")
551 read(unitpar, rc_list, end=111)
552111 close(unitpar)
553 end do
554
555 fl%ncool=ncool
556 fl%coolcurve=coolcurve
557 fl%tlow=tlow
558 fl%Tfix=tfix
559 fl%rc_split=rc_split
560 fl%rad_cut=rad_cut
561 fl%rad_cut_hgt=rad_cut_hgt
562 fl%rad_cut_dey=rad_cut_dey
563 end subroutine rc_params_read
564
565 subroutine ffhd_check_params
569
570 gamma_1=ffhd_gamma-1.d0
571 if (.not. ffhd_energy) then
572 if (ffhd_gamma <= 0.0d0) call mpistop ("Error: ffhd_gamma <= 0")
573 if (ffhd_adiab < 0.0d0) call mpistop ("Error: ffhd_adiab < 0")
575 else
576 if (ffhd_gamma <= 0.0d0 .or. ffhd_gamma == 1.0d0) &
577 call mpistop ("Error: ffhd_gamma <= 0 or ffhd_gamma == 1")
578 inv_gamma_1=1.d0/gamma_1
579 small_e = small_pressure * inv_gamma_1
580 end if
581
582 if (number_equi_vars > 0 .and. .not. associated(usr_set_equi_vars)) then
583 call mpistop("usr_set_equi_vars has to be implemented in the user file")
584 end if
585 end subroutine ffhd_check_params
586
587 subroutine ffhd_physical_units()
589 double precision :: mp,kb
590 double precision :: a,b
591
592 if(si_unit) then
593 mp=mp_si
594 kb=kb_si
595 else
596 mp=mp_cgs
597 kb=kb_cgs
598 end if
599 if(eq_state_units) then
600 a=1d0+4d0*he_abundance
602 b=1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0)
603 else
604 b=2d0+3d0*he_abundance
605 end if
606 rr=1d0
607 else
608 a=1d0
609 b=1d0
610 rr=(1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0))/(1d0+4d0*he_abundance)
611 end if
612 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
613 if(unit_density/=1.d0) then
615 else if(unit_numberdensity/=1.d0) then
617 end if
618 if(unit_temperature/=1.d0) then
621 if(unit_length/=1.d0) then
623 else if(unit_time/=1.d0) then
625 end if
626 else if(unit_pressure/=1.d0) then
629 if(unit_length/=1.d0) then
631 else if(unit_time/=1.d0) then
633 end if
634 else if(unit_velocity/=1.d0) then
637 if(unit_length/=1.d0) then
639 else if(unit_time/=1.d0) then
641 end if
642 else if(unit_time/=1.d0) then
646 end if
647 else if(unit_temperature/=1.d0) then
648 ! units of temperature and velocity are dependent
649 if(unit_pressure/=1.d0) then
653 if(unit_length/=1.d0) then
655 else if(unit_time/=1.d0) then
657 end if
658 end if
659 else if(unit_pressure/=1.d0) then
660 if(unit_velocity/=1.d0) then
664 if(unit_length/=1.d0) then
666 else if(unit_time/=1.d0) then
668 end if
669 else if(unit_time/=0.d0) then
674 end if
675 end if
677 end subroutine ffhd_physical_units
678
679 subroutine ffhd_check_w_origin(primitive,ixI^L,ixO^L,w,flag)
681 logical, intent(in) :: primitive
682 integer, intent(in) :: ixi^l, ixo^l
683 double precision, intent(in) :: w(ixi^s,nw)
684 double precision :: tmp(ixi^s)
685 logical, intent(inout) :: flag(ixi^s,1:nw)
686
687 flag=.false.
688 where(w(ixo^s,rho_) < small_density) flag(ixo^s,rho_) = .true.
689
690 if(ffhd_energy) then
691 if(primitive) then
692 where(w(ixo^s,e_) < small_pressure) flag(ixo^s,e_) = .true.
693 else
694 tmp(ixo^s)=w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l)
695 where(tmp(ixo^s) < small_e) flag(ixo^s,e_) = .true.
696 end if
697 end if
698 end subroutine ffhd_check_w_origin
699
700 subroutine ffhd_to_conserved_origin(ixI^L,ixO^L,w,x)
702 integer, intent(in) :: ixi^l, ixo^l
703 double precision, intent(inout) :: w(ixi^s, nw)
704 double precision, intent(in) :: x(ixi^s, 1:ndim)
705
706 if(ffhd_energy) then
707 w(ixo^s,e_)=w(ixo^s,p_)*inv_gamma_1+half*w(ixo^s,mom(1))**2*w(ixo^s,rho_)
708 end if
709 w(ixo^s,mom(1))=w(ixo^s,rho_)*w(ixo^s,mom(1))
710 end subroutine ffhd_to_conserved_origin
711
712 subroutine ffhd_to_primitive_origin(ixI^L,ixO^L,w,x)
714 integer, intent(in) :: ixi^l, ixo^l
715 double precision, intent(inout) :: w(ixi^s, nw)
716 double precision, intent(in) :: x(ixi^s, 1:ndim)
717
718 if(fix_small_values) then
719 call ffhd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'ffhd_to_primitive_origin')
720 end if
721
722 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/w(ixo^s,rho_)
723 if(ffhd_energy) then
724 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-half*w(ixo^s,rho_)*w(ixo^s,mom(1))**2)
725 end if
726 end subroutine ffhd_to_primitive_origin
727
728 subroutine ffhd_ei_to_e(ixI^L,ixO^L,w,x)
730 integer, intent(in) :: ixi^l, ixo^l
731 double precision, intent(inout) :: w(ixi^s, nw)
732 double precision, intent(in) :: x(ixi^s, 1:ndim)
733
734 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
735 end subroutine ffhd_ei_to_e
736
737 subroutine ffhd_e_to_ei(ixI^L,ixO^L,w,x)
739 integer, intent(in) :: ixi^l, ixo^l
740 double precision, intent(inout) :: w(ixi^s, nw)
741 double precision, intent(in) :: x(ixi^s, 1:ndim)
742
743 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
744 if(fix_small_values) then
745 call ffhd_handle_small_ei(w,x,ixi^l,ixi^l,e_,'ffhd_e_to_ei')
746 end if
747 end subroutine ffhd_e_to_ei
748
749 subroutine ffhd_handle_small_values_origin(primitive, w, x, ixI^L, ixO^L, subname)
752 logical, intent(in) :: primitive
753 integer, intent(in) :: ixi^l,ixo^l
754 double precision, intent(inout) :: w(ixi^s,1:nw)
755 double precision, intent(in) :: x(ixi^s,1:ndim)
756 character(len=*), intent(in) :: subname
757
758 logical :: flag(ixi^s,1:nw)
759 double precision :: tmp2(ixi^s)
760
761 call phys_check_w(primitive, ixi^l, ixi^l, w, flag)
762
763 if(any(flag)) then
764 select case (small_values_method)
765 case ("replace")
766 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
767 if(small_values_fix_iw(mom(1))) then
768 where(flag(ixo^s,rho_)) w(ixo^s, mom(1)) = 0.0d0
769 end if
770 if(ffhd_energy) then
771 if(primitive) then
772 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
773 else
774 where(flag(ixo^s,e_))
775 w(ixo^s,e_) = small_e+ffhd_kin_en(w,ixi^l,ixo^l)
776 end where
777 end if
778 end if
779 case ("average")
780 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
781 if(ffhd_energy) then
782 if(primitive) then
783 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
784 else
785 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
786 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
787 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
788 end if
789 end if
790 case default
791 if(.not.primitive) then
792 if(ffhd_energy) then
793 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
794 end if
795 w(ixo^s,mom(1))=w(ixo^s,mom(1))/w(ixo^s,rho_)
796 end if
797 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
798 end select
799 end if
800 end subroutine ffhd_handle_small_values_origin
801
802 subroutine ffhd_get_v_origin(w,x,ixI^L,ixO^L,v)
804 integer, intent(in) :: ixi^l, ixo^l
805 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
806 double precision, intent(out) :: v(ixi^s,ndir)
807 double precision :: rho(ixi^s)
808 integer :: idir
809
810 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
811 rho(ixo^s)=1.d0/rho(ixo^s)
812 do idir=1,ndir
813 v(ixo^s,ndir) = w(ixo^s,mom(1))*block%B0(ixo^s,idir,0)*rho(ixo^s)
814 end do
815 end subroutine ffhd_get_v_origin
816
817 subroutine ffhd_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
819 integer, intent(in) :: ixi^l, ixo^l, idim
820 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
821 double precision, intent(out) :: v(ixi^s)
822 double precision :: rho(ixi^s)
823
824 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
825 v(ixo^s) = (w(ixo^s, mom(1))*block%B0(ixo^s,idim,0)) / rho(ixo^s)
826 end subroutine ffhd_get_v_idim
827
828 subroutine ffhd_get_cmax_origin(wprim,x,ixI^L,ixO^L,idim,cmax)
830 integer, intent(in) :: ixi^l, ixo^l, idim
831 ! w in primitive form
832 double precision, intent(in) :: wprim(ixi^s, nw), x(ixi^s,1:ndim)
833 double precision, intent(inout) :: cmax(ixi^s)
834
835 if(ffhd_energy) then
836 cmax(ixo^s)=dsqrt(ffhd_gamma*wprim(ixo^s,p_)/wprim(ixo^s,rho_))
837 else
838 cmax(ixo^s)=dsqrt(ffhd_gamma*ffhd_adiab*wprim(ixo^s,rho_)**gamma_1)
839 end if
840 cmax(ixo^s)=dabs(wprim(ixo^s,mom(1))*block%B0(ixo^s,idim,0))+cmax(ixo^s)
841
842 end subroutine ffhd_get_cmax_origin
843
844 subroutine ffhd_get_tcutoff(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
846 use mod_geometry
847 integer, intent(in) :: ixi^l,ixo^l
848 double precision, intent(in) :: x(ixi^s,1:ndim)
849 double precision, intent(inout) :: w(ixi^s,1:nw)
850 double precision, intent(out) :: tco_local,tmax_local
851 double precision, parameter :: trac_delta=0.25d0
852 double precision :: te(ixi^s),lts(ixi^s)
853 double precision, dimension(ixI^S,1:ndim) :: gradt
854 double precision :: bdir(ndim)
855 double precision :: ltrc,ltrp,altr
856 integer :: idims,ix^d,jxo^l,hxo^l,ixa^d,ixb^d
857 integer :: jxp^l,hxp^l,ixp^l,ixq^l
858 logical :: lrlt(ixi^s)
859
860 call ffhd_get_temperature(w,x,ixi^l,ixi^l,te)
861 tco_local=zero
862 tmax_local=maxval(te(ixo^s))
863
864 {^ifoned
865 select case(ffhd_trac_type)
866 case(0)
867 !> test case, fixed cutoff temperature
868 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
869 case(1)
870 hxo^l=ixo^l-1;
871 jxo^l=ixo^l+1;
872 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
873 lrlt=.false.
874 where(lts(ixo^s) > trac_delta)
875 lrlt(ixo^s)=.true.
876 end where
877 if(any(lrlt(ixo^s))) then
878 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
879 end if
880 case(2)
881 !> iijima et al. 2021, LTRAC method
882 ltrc=1.5d0
883 ltrp=4.d0
884 ixp^l=ixo^l^ladd1;
885 hxo^l=ixo^l-1;
886 jxo^l=ixo^l+1;
887 hxp^l=ixp^l-1;
888 jxp^l=ixp^l+1;
889 lts(ixp^s)=0.5d0*dabs(te(jxp^s)-te(hxp^s))/te(ixp^s)
890 lts(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
891 lts(ixo^s)=0.25d0*(lts(jxo^s)+two*lts(ixo^s)+lts(hxo^s))
892 block%wextra(ixo^s,tcoff_)=te(ixo^s)*lts(ixo^s)**0.4d0
893 case default
894 call mpistop("ffhd_trac_type not allowed for 1D simulation")
895 end select
896 }
897 {^nooned
898 select case(ffhd_trac_type)
899 case(0)
900 !> test case, fixed cutoff temperature
901 if(slab_uniform) then
902 !> assume cgs units
903 block%wextra(ixi^s,tcoff_)=max(min(3.d5/unit_temperature,6.d5/unit_temperature-3.d-4/unit_temperature*unit_length*x(ixi^s,ndim)),zero)
904 else
905 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
906 end if
907 case(1,4,6)
908 do idims=1,ndim
909 call gradient(te,ixi^l,ixo^l,idims,gradt(ixi^s,idims))
910 end do
911 if(ffhd_trac_type .gt. 1) then
912 ! B direction at cell center
913 bdir=zero
914 {do ixa^d=0,1\}
915 ixb^d=(ixomin^d+ixomax^d-1)/2+ixa^d;
916 bdir(1:ndim)=bdir(1:ndim)+block%B0(ixb^d,1:ndim,0)
917 {end do\}
918 if(sum(bdir(:)**2) .gt. zero) then
919 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
920 end if
921 block%special_values(3:ndim+2)=bdir(1:ndim)
922 end if
923 {do ix^db=ixomin^db,ixomax^db\}
924 {^iftwod
925 ! temperature length scale inversed
926 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2))*&
927 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
928 }
929 {^ifthreed
930 ! temperature length scale inversed
931 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2),block%ds(ix^d,3))*&
932 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
933 }
934 if(lts(ix^d)>trac_delta) then
935 block%special_values(1)=max(block%special_values(1),te(ix^d))
936 end if
937 {end do\}
938 block%special_values(2)=tmax_local
939 case(2)
940 !> iijima et al. 2021, LTRAC method
941 ltrc=1.5d0
942 ltrp=4.d0
943 ixp^l=ixo^l^ladd2;
944 do idims=1,ndim
945 ixq^l=ixp^l;
946 hxp^l=ixp^l;
947 jxp^l=ixp^l;
948 select case(idims)
949 {case(^d)
950 ixqmin^d=ixqmin^d+1
951 ixqmax^d=ixqmax^d-1
952 hxpmax^d=ixpmin^d
953 jxpmin^d=ixpmax^d
954 \}
955 end select
956 call gradient(te,ixi^l,ixq^l,idims,gradt(ixi^s,idims))
957 call gradientf(te,x,ixi^l,hxp^l,idims,gradt(ixi^s,idims),nghostcells,.true.)
958 call gradientf(te,x,ixi^l,jxp^l,idims,gradt(ixi^s,idims),nghostcells,.false.)
959 end do
960 {do ix^db=ixpmin^db,ixpmax^db\}
961 ! temperature length scale inversed
962 lts(ix^d)=abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
963 ! fraction of cells size to temperature length scale
964 lts(ix^d)=min(^d&block%ds({ix^d},^d))*lts(ix^d)
965 lts(ix^d)=max(one,(exp(lts(ix^d))/ltrc)**ltrp)
966 {end do\}
967
968 ! need one ghost layer for thermal conductivity
969 ixp^l=ixo^l^ladd1;
970 {do ix^db=ixpmin^db,ixpmax^db\}
971 {^iftwod
972 altr=0.25d0*((lts(ix1-1,ix2)+two*lts(ix^d)+lts(ix1+1,ix2))*block%B0(ix^d,1,0)**2+&
973 (lts(ix1,ix2-1)+two*lts(ix^d)+lts(ix1,ix2+1))*block%B0(ix^d,2,0)**2)
974 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
975 }
976 {^ifthreed
977 altr=0.25d0*((lts(ix1-1,ix2,ix3)+two*lts(ix^d)+lts(ix1+1,ix2,ix3))*block%B0(ix^d,1,0)**2+&
978 (lts(ix1,ix2-1,ix3)+two*lts(ix^d)+lts(ix1,ix2+1,ix3))*block%B0(ix^d,2,0)**2+&
979 (lts(ix1,ix2,ix3-1)+two*lts(ix^d)+lts(ix1,ix2,ix3+1))*block%B0(ix^d,3,0)**2)
980 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
981 }
982 {end do\}
983 case(3,5)
984 !> do nothing here
985 case default
986 call mpistop("unknown ffhd_trac_type")
987 end select
988 }
989 end subroutine ffhd_get_tcutoff
990
991 subroutine ffhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
993 integer, intent(in) :: ixi^l, ixo^l, idim
994 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
995 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
996 double precision, intent(in) :: x(ixi^s,1:ndim)
997 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
998 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
999 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
1000 double precision :: wmean(ixi^s,nw)
1001 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1002
1003 select case (boundspeed)
1004 case (1)
1005 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1006 ! Methods for Fluid Dynamics" by Toro.
1007 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1008 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1009 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1010 umean(ixo^s)=(wlp(ixo^s,mom(1))*tmp1(ixo^s)&
1011 +wrp(ixo^s,mom(1))*tmp2(ixo^s))*tmp3(ixo^s)
1012 umean(ixo^s)=umean(ixo^s)*block%B0(ixo^s,idim,idim)
1013 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1014 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1015 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1016 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1017 ((wrp(ixo^s,mom(1))-wlp(ixo^s,mom(1)))*block%B0(ixo^s,idim,idim))**2
1018 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1019 if(present(cmin)) then
1020 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1021 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1022 else
1023 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1024 end if
1025 case (2)
1026 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1027 tmp1(ixo^s)=wmean(ixo^s,mom(1))*block%B0(ixo^s,idim,idim)/wmean(ixo^s,rho_)
1028 call ffhd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1029 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1030 if(present(cmin)) then
1031 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1032 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1033 else
1034 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1035 end if
1036 case (3)
1037 ! Miyoshi 2005 JCP 208, 315 equation (67)
1038 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1039 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1040 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1041 if(present(cmin)) then
1042 cmin(ixo^s,1)=min(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1043 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))-csoundl(ixo^s)
1044 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1045 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1046 else
1047 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1048 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1049 end if
1050 end select
1051 end subroutine ffhd_get_cbounds
1052
1053 subroutine ffhd_get_pthermal_iso(w,x,ixI^L,ixO^L,pth)
1055
1056 integer, intent(in) :: ixi^l, ixo^l
1057 double precision, intent(in) :: w(ixi^s,nw)
1058 double precision, intent(in) :: x(ixi^s,1:ndim)
1059 double precision, intent(out):: pth(ixi^s)
1060
1061 call ffhd_get_rho(w,x,ixi^l,ixo^l,pth)
1062 pth(ixo^s)=ffhd_adiab*pth(ixo^s)**ffhd_gamma
1063 end subroutine ffhd_get_pthermal_iso
1064
1065 subroutine ffhd_get_pthermal_origin(w,x,ixI^L,ixO^L,pth)
1068 integer, intent(in) :: ixi^l, ixo^l
1069 double precision, intent(in) :: w(ixi^s,nw)
1070 double precision, intent(in) :: x(ixi^s,1:ndim)
1071 double precision, intent(out):: pth(ixi^s)
1072 integer :: iw, ix^d
1073
1074 pth(ixo^s)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
1075 if (fix_small_values) then
1076 {do ix^db= ixo^lim^db\}
1077 if(pth(ix^d)<small_pressure) then
1078 pth(ix^d)=small_pressure
1079 end if
1080 {end do^D&\}
1081 elseif(check_small_values) then
1082 {do ix^db= ixo^lim^db\}
1083 if(pth(ix^d)<small_pressure) then
1084 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1085 " encountered when call ffhd_get_pthermal"
1086 write(*,*) "Iteration: ", it, " Time: ", global_time
1087 write(*,*) "Location: ", x(ix^d,:)
1088 write(*,*) "Cell number: ", ix^d
1089 do iw=1,nw
1090 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1091 end do
1092 ! use erroneous arithmetic operation to crash the run
1093 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1094 write(*,*) "Saving status at the previous time step"
1095 crash=.true.
1096 end if
1097 {end do^D&\}
1098 end if
1099 end subroutine ffhd_get_pthermal_origin
1100
1101 subroutine ffhd_get_temperature_from_te(w, x, ixI^L, ixO^L, res)
1103 integer, intent(in) :: ixi^l, ixo^l
1104 double precision, intent(in) :: w(ixi^s, 1:nw)
1105 double precision, intent(in) :: x(ixi^s, 1:ndim)
1106 double precision, intent(out):: res(ixi^s)
1107
1108 res(ixo^s) = w(ixo^s, te_)
1109 end subroutine ffhd_get_temperature_from_te
1110
1111 subroutine ffhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1113 integer, intent(in) :: ixi^l, ixo^l
1114 double precision, intent(in) :: w(ixi^s, 1:nw)
1115 double precision, intent(in) :: x(ixi^s, 1:ndim)
1116 double precision, intent(out):: res(ixi^s)
1117 double precision :: r(ixi^s)
1118
1119 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1120 res(ixo^s) = gamma_1 * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1121 end subroutine ffhd_get_temperature_from_eint
1122
1123 subroutine ffhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1125 integer, intent(in) :: ixi^l, ixo^l
1126 double precision, intent(in) :: w(ixi^s, 1:nw)
1127 double precision, intent(in) :: x(ixi^s, 1:ndim)
1128 double precision, intent(out):: res(ixi^s)
1129
1130 double precision :: r(ixi^s)
1131
1132 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1133 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,res)
1134 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1135 end subroutine ffhd_get_temperature_from_etot
1136
1137 subroutine ffhd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1139 integer, intent(in) :: ixi^l, ixo^l
1140 double precision, intent(in) :: w(ixi^s,nw)
1141 double precision, intent(in) :: x(ixi^s,1:ndim)
1142 double precision, intent(out) :: csound2(ixi^s)
1143 double precision :: rho(ixi^s)
1144
1145 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1146 if(ffhd_energy) then
1147 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1148 csound2(ixo^s)=ffhd_gamma*csound2(ixo^s)/rho(ixo^s)
1149 else
1150 csound2(ixo^s)=ffhd_gamma*ffhd_adiab*rho(ixo^s)**gamma_1
1151 end if
1152 end subroutine ffhd_get_csound2
1153
1154 subroutine ffhd_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
1156 use mod_geometry
1157 integer, intent(in) :: ixi^l, ixo^l, idim
1158 ! conservative w
1159 double precision, intent(in) :: wc(ixi^s,nw)
1160 ! primitive w
1161 double precision, intent(in) :: w(ixi^s,nw)
1162 double precision, intent(in) :: x(ixi^s,1:ndim)
1163 double precision,intent(out) :: f(ixi^s,nwflux)
1164 double precision :: ptotal(ixo^s)
1165
1166 f(ixo^s,rho_)=w(ixo^s,mom(1))*w(ixo^s,rho_)*block%B0(ixo^s,idim,idim)
1167
1168 if(ffhd_energy) then
1169 ptotal(ixo^s)=w(ixo^s,p_)
1170 else
1171 ptotal(ixo^s)=ffhd_adiab*w(ixo^s,rho_)**ffhd_gamma
1172 end if
1173
1174 ! Get flux of momentum
1175 f(ixo^s,mom(1))=(wc(ixo^s,mom(1))*w(ixo^s,mom(1))+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1176
1177 ! Get flux of energy
1178 if(ffhd_energy) then
1179 f(ixo^s,e_)=w(ixo^s,mom(1))*(wc(ixo^s,e_)+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1181 f(ixo^s,e_)=f(ixo^s,e_)+w(ixo^s,q_)*block%B0(ixo^s,idim,idim)
1182 f(ixo^s,q_)=zero
1183 end if
1184 end if
1185 end subroutine ffhd_get_flux
1186
1187 subroutine ffhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1191 integer, intent(in) :: ixi^l, ixo^l
1192 double precision, intent(in) :: qdt,dtfactor
1193 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:ndim)
1194 double precision, intent(inout) :: w(ixi^s,1:nw)
1195 logical, intent(in) :: qsourcesplit
1196 logical, intent(inout) :: active
1197
1198 if (.not. qsourcesplit) then
1199 active = .true.
1200 call add_punitb(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1202 call add_hypertc_source(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1203 end if
1204 end if
1205
1206 if(ffhd_radiative_cooling) then
1207 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1208 w,x,qsourcesplit,active, rc_fl)
1209 end if
1210
1211 if(ffhd_gravity) then
1212 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1213 w,x,ffhd_energy,qsourcesplit,active)
1214 end if
1215
1216 ! update temperature from new pressure, density, and old ionization degree
1218 if(.not.qsourcesplit) then
1219 active = .true.
1220 call ffhd_update_temperature(ixi^l,ixo^l,wct,w,x)
1221 end if
1222 end if
1223 end subroutine ffhd_add_source
1224
1225 subroutine add_punitb(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1227 use mod_geometry
1228 integer, intent(in) :: ixi^l,ixo^l
1229 double precision, intent(in) :: qdt
1230 double precision, intent(in) :: wct(ixi^s,1:nw),x(ixi^s,1:ndim)
1231 double precision, intent(in) :: wctprim(ixi^s,1:nw)
1232 double precision, intent(inout) :: w(ixi^s,1:nw)
1233
1234 integer :: idims,hxo^l
1235 double precision :: divb(ixi^s)
1236 double precision :: rhovpar(ixi^s),gradrhov(ixi^s)
1237
1238 divb=zero
1239 if(slab_uniform) then
1240 do idims=1,ndim
1241 hxo^l=ixo^l-kr(idims,^d);
1242 divb(ixo^s)=divb(ixo^s)+(block%B0(ixo^s,idims,idims)-block%B0(hxo^s,idims,idims))/dxlevel(idims)
1243 end do
1244 else
1245 call divvector(block%B0(ixi^s,1:ndir,0),ixi^l,ixo^l,divb)
1246 end if
1247 w(ixo^s,mom(1))=w(ixo^s,mom(1))+qdt*wctprim(ixo^s,p_)*divb(ixo^s)
1248 end subroutine add_punitb
1249
1250 subroutine ffhd_get_rho(w,x,ixI^L,ixO^L,rho)
1252 integer, intent(in) :: ixi^l, ixo^l
1253 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
1254 double precision, intent(out) :: rho(ixi^s)
1255
1256 rho(ixo^s) = w(ixo^s,rho_)
1257 end subroutine ffhd_get_rho
1258
1259 subroutine ffhd_handle_small_ei(w, x, ixI^L, ixO^L, ie, subname)
1262 integer, intent(in) :: ixi^l,ixo^l, ie
1263 double precision, intent(inout) :: w(ixi^s,1:nw)
1264 double precision, intent(in) :: x(ixi^s,1:ndim)
1265 character(len=*), intent(in) :: subname
1266 integer :: idir
1267 logical :: flag(ixi^s,1:nw)
1268 double precision :: rho(ixi^s)
1269
1270 flag=.false.
1271 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
1272 if(any(flag(ixo^s,ie))) then
1273 select case (small_values_method)
1274 case ("replace")
1275 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
1276 case ("average")
1277 call small_values_average(ixi^l, ixo^l, w, x, flag, ie)
1278 case default
1279 w(ixo^s,e_)=w(ixo^s,e_)*gamma_1
1280 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1281 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/rho(ixo^s)
1282 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1283 end select
1284 end if
1285 end subroutine ffhd_handle_small_ei
1286
1287 subroutine ffhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1290 integer, intent(in) :: ixi^l, ixo^l
1291 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1292 double precision, intent(inout) :: w(ixi^s,1:nw)
1293 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1294
1295 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1296 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1297 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1298 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1299 end subroutine ffhd_update_temperature
1300
1301 subroutine ffhd_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
1303 use mod_usr_methods
1304 use mod_gravity, only: gravity_get_dt
1305 integer, intent(in) :: ixi^l, ixo^l
1306 double precision, intent(inout) :: dtnew
1307 double precision, intent(in) :: dx^d
1308 double precision, intent(in) :: w(ixi^s,1:nw)
1309 double precision, intent(in) :: x(ixi^s,1:ndim)
1310
1311 dtnew = bigdouble
1312
1313 if(ffhd_gravity) then
1314 call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1315 end if
1316 end subroutine ffhd_get_dt
1317
1318 subroutine ffhd_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
1320 integer, intent(in) :: ixi^l, ixo^l
1321 double precision, intent(in) :: qdt, dtfactor,x(ixi^s,1:ndim)
1322 double precision, intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw),w(ixi^s,1:nw)
1323
1324 ! no geometric source terms needed for ffhd, no divergences of tensors
1325 end subroutine ffhd_add_source_geom
1326
1327 function ffhd_kin_en_origin(w, ixI^L, ixO^L, inv_rho) result(ke)
1328 use mod_global_parameters, only: nw, ndim,block
1329 integer, intent(in) :: ixi^l, ixo^l
1330 double precision, intent(in) :: w(ixi^s, nw)
1331 double precision :: ke(ixo^s)
1332 double precision, intent(in), optional :: inv_rho(ixo^s)
1333
1334 if(present(inv_rho)) then
1335 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2*inv_rho(ixo^s)
1336 else
1337 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2/w(ixo^s,rho_)
1338 end if
1339 end function ffhd_kin_en_origin
1340
1341 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1344 integer, intent(in) :: ixi^l, ixo^l
1345 double precision, intent(in) :: w(ixi^s,1:nw)
1346 double precision, intent(in) :: x(ixi^s,1:ndim)
1347 double precision, intent(out):: rfactor(ixi^s)
1348 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1349
1350 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1351 rfactor(ixo^s)=(1.d0+iz_h(ixo^s)+0.1d0*(1.d0+iz_he(ixo^s)*(1.d0+iz_he(ixo^s))))/(2.d0+3.d0*he_abundance)
1352 end subroutine rfactor_from_temperature_ionization
1353
1354 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1356 integer, intent(in) :: ixi^l, ixo^l
1357 double precision, intent(in) :: w(ixi^s,1:nw)
1358 double precision, intent(in) :: x(ixi^s,1:ndim)
1359 double precision, intent(out):: rfactor(ixi^s)
1360
1361 rfactor(ixo^s)=rr
1362 end subroutine rfactor_from_constant_ionization
1363
1364 subroutine add_hypertc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1366
1367 integer, intent(in) :: ixi^l,ixo^l
1368 double precision, intent(in) :: qdt
1369 double precision, dimension(ixI^S,1:ndim), intent(in) :: x
1370 double precision, dimension(ixI^S,1:nw), intent(in) :: wct,wctprim
1371 double precision, dimension(ixI^S,1:nw), intent(inout) :: w
1372
1373 double precision, dimension(ixI^S) :: te
1374 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1375 integer :: ix^d
1376
1377 te(ixi^s)=wctprim(ixi^s,p_)/wct(ixi^s,rho_)
1378 ! temperature on face T_(i+1/2)=(7(T_i+T_(i+1))-(T_(i-1)+T_(i+2)))/12
1379 ! T_(i+1/2)-T_(i-1/2)=(8(T_(i+1)-T_(i-1))-T_(i+2)+T_(i-2))/12
1380 {^iftwod
1381 do ix2=ixomin2,ixomax2
1382 do ix1=ixomin1,ixomax1
1383 if(ffhd_trac) then
1384 if(te(ix^d)<block%wextra(ix^d,tcoff_)) then
1385 sigma_t5=hypertc_kappa*sqrt(block%wextra(ix^d,tcoff_)**5)
1386 sigma_t7=sigma_t5*block%wextra(ix^d,tcoff_)
1387 else
1388 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1389 sigma_t7=sigma_t5*te(ix^d)
1390 end if
1391 else
1392 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1393 sigma_t7=sigma_t5*te(ix^d)
1394 end if
1395 sigmat5_bgradt=sigma_t5*(&
1396 block%B0(ix^d,1,0)*((8.d0*(te(ix1+1,ix2)-te(ix1-1,ix2))-te(ix1+2,ix2)+te(ix1-2,ix2))/12.d0)/block%ds(ix^d,1)&
1397 +block%B0(ix^d,2,0)*((8.d0*(te(ix1,ix2+1)-te(ix1,ix2-1))-te(ix1,ix2+2)+te(ix1,ix2-2))/12.d0)/block%ds(ix^d,2))
1398 if(ffhd_htc_sat) then
1399 ! 5 phi rho c^3, phi=0.3, c=sqrt(p/rho) isothermal sound speed
1400 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^d,rho_)*(wctprim(ix^d,p_)/wct(ix^d,rho_))**1.5d0))
1401 tau=max(4.d0*dt, f_sat*sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1402 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1403 else
1404 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1405 max(4.d0*dt, sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1406 end if
1407 end do
1408 end do
1409 }
1410 {^ifthreed
1411 do ix3=ixomin3,ixomax3
1412 do ix2=ixomin2,ixomax2
1413 do ix1=ixomin1,ixomax1
1414 if(ffhd_trac) then
1415 if(te(ix^d)<block%wextra(ix^d,tcoff_)) then
1416 sigma_t5=hypertc_kappa*sqrt(block%wextra(ix^d,tcoff_)**5)
1417 sigma_t7=sigma_t5*block%wextra(ix^d,tcoff_)
1418 else
1419 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1420 sigma_t7=sigma_t5*te(ix^d)
1421 end if
1422 else
1423 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1424 sigma_t7=sigma_t5*te(ix^d)
1425 end if
1426 sigmat5_bgradt=sigma_t5*(&
1427 block%B0(ix^d,1,0)*((8.d0*(te(ix1+1,ix2,ix3)-te(ix1-1,ix2,ix3))-te(ix1+2,ix2,ix3)+te(ix1-2,ix2,ix3))/12.d0)/block%ds(ix^d,1)&
1428 +block%B0(ix^d,2,0)*((8.d0*(te(ix1,ix2+1,ix3)-te(ix1,ix2-1,ix3))-te(ix1,ix2+2,ix3)+te(ix1,ix2-2,ix3))/12.d0)/block%ds(ix^d,2)&
1429 +block%B0(ix^d,3,0)*((8.d0*(te(ix1,ix2,ix3+1)-te(ix1,ix2,ix3-1))-te(ix1,ix2,ix3+2)+te(ix1,ix2,ix3-2))/12.d0)/block%ds(ix^d,3))
1430 if(ffhd_htc_sat) then
1431 ! 5 phi rho c^3, phi=0.3, c=sqrt(p/rho) isothermal sound speed
1432 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^d,rho_)*(wctprim(ix^d,p_)/wct(ix^d,rho_))**1.5d0))
1433 tau=max(4.d0*dt, f_sat*sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1434 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1435 else
1436 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1437 max(4.d0*dt, sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1438 end if
1439 end do
1440 end do
1441 end do
1442 }
1443 end subroutine add_hypertc_source
1444
1445end module mod_ffhd_phys
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter bigdouble
A very large real number.
subroutine add_convert_method(phys_convert_vars, nwc, dataset_names, file_suffix)
Definition mod_convert.t:59
Frozen-field hydrodynamics module.
integer, public, protected te_
Indices of temperature.
integer, public, protected ffhd_trac_type
Which TRAC method is used.
double precision, public hypertc_kappa
The thermal conductivity kappa in hyperbolic thermal conduction.
integer, public, protected e_
Index of the energy density (-1 if not present)
double precision, public ffhd_gamma
The adiabatic index.
logical, public, protected eq_state_units
double precision, public ffhd_adiab
The adiabatic constant.
procedure(sub_get_pthermal), pointer, public ffhd_get_temperature
logical, public, protected ffhd_hyperbolic_thermal_conduction
Whether hyperbolic type thermal conduction is used.
type(rc_fluid), allocatable, public rc_fl
type of fluid for radiative cooling
integer, public, protected rho_
Index of the density (in the w array)
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
procedure(sub_get_pthermal), pointer, public ffhd_get_pthermal
logical, public, protected ffhd_htc_sat
Wheterh saturation is considered for hyperbolic TC.
subroutine, public ffhd_get_rho(w, x, ixil, ixol, rho)
logical, public, protected ffhd_partial_ionization
Whether plasma is partially ionized.
procedure(sub_convert), pointer, public ffhd_to_conserved
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
procedure(fun_kin_en), pointer, public ffhd_kin_en
subroutine, public ffhd_phys_init()
procedure(sub_get_v), pointer, public ffhd_get_v
subroutine, public ffhd_get_csound2(w, x, ixil, ixol, csound2)
logical, public, protected ffhd_energy
Whether an energy equation is used.
double precision, public, protected rr
type(tc_fluid), allocatable, public tc_fl
type of fluid for thermal conduction
type(te_fluid), allocatable, public te_fl_ffhd
type of fluid for thermal emission synthesis
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
logical, public, protected ffhd_radiative_cooling
Whether radiative cooling is added.
subroutine, public ffhd_get_v_idim(w, x, ixil, ixol, idim, v)
integer, public, protected q_
procedure(sub_convert), pointer, public ffhd_to_primitive
integer, public, protected tweight_
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
logical, public, protected ffhd_trac
Whether TRAC method is used.
logical, public, protected ffhd_thermal_conduction
Whether thermal conduction is used.
subroutine, public ffhd_ei_to_e(ixil, ixol, w, x)
logical, public, protected ffhd_gravity
Whether gravity is added.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
double precision, public, protected ffhd_trac_mask
Height of the mask used in the TRAC method.
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
integer, public, protected ffhd_trac_finegrid
Distance between two adjacent traced magnetic field lines (in finest cell size)
subroutine, public ffhd_e_to_ei(ixil, ixol, w, x)
Module for flux conservation near refinement boundaries.
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
subroutine divvector(qvec, ixil, ixol, divq, nth_in)
subroutine gradient(q, ixil, ixol, idir, gradq, nth_in)
This module contains definitions of global parameters and variables and some generic functions/subrou...
type(state), pointer block
Block pointer for using one block and its previous state.
double precision unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
integer, parameter unitpar
file handle for IO
double precision unit_mass
Physical scaling factor for mass.
integer, dimension(3, 3) kr
Kronecker delta tensor.
double precision unit_numberdensity
Physical scaling factor for number density.
character(len=std_len) convert_type
Which format to use when converting.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
double precision cmax_global
global fastest wave speed needed in fd scheme and glm method
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
double precision dt
global time step
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision courantpar
The Courant (CFL) number used for the simulation.
double precision unit_velocity
Physical scaling factor for velocity.
logical b0field
split magnetic field as background B0 field
double precision unit_temperature
Physical scaling factor for temperature.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
spatial steps for all dimensions at all levels
logical phys_trac
Use TRAC for MHD or 1D HD.
logical need_global_cmax
need global maximal wave speed
logical fix_small_values
fix small values with average or replace methods
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical slab_uniform
uniform Cartesian geometry or not (stretched Cartesian)
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
integer number_equi_vars
number of equilibrium set variables, besides the mag field
Module for including gravity in (magneto)hydrodynamics simulations.
Definition mod_gravity.t:2
subroutine gravity_get_dt(w, ixil, ixol, dtnew, dxd, x)
Definition mod_gravity.t:81
subroutine gravity_init()
Initialize the module.
Definition mod_gravity.t:26
subroutine gravity_add_source(qdt, ixil, ixol, wct, wctprim, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
Definition mod_gravity.t:43
module ionization degree - get ionization degree for given temperature
subroutine ionization_degree_from_temperature(ixil, ixol, te, iz_h, iz_he)
Module containing all the particle routines.
This module defines the procedures of a physics module. It contains function pointers for the various...
Definition mod_physics.t:4
module radiative cooling – add optically thin radiative cooling
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine radiative_cooling_init(fl, read_params)
subroutine radiative_cooling_add_source(qdt, ixil, ixol, wct, wctprim, w, x, qsourcesplit, active, fl)
Module for handling problematic values in simulations, such as negative pressures.
subroutine, public small_values_average(ixil, ixol, w, x, w_flag, windex)
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
subroutine, public small_values_error(wprim, x, ixil, ixol, w_flag, subname)
logical, dimension(:), allocatable, public small_values_fix_iw
Whether to apply small value fixes to certain variables.
character(len=20), public small_values_method
How to handle small values.
Generic supertimestepping method which can be used for multiple source terms in the governing equatio...
subroutine, public add_sts_method(sts_getdt, sts_set_sources, startvar, nflux, startwbc, nwbc, evolve_b)
subroutine which added programatically a term to be calculated using STS Params: sts_getdt function c...
subroutine, public set_conversion_methods_to_head(sts_before_first_cycle, sts_after_last_cycle)
Set the hooks called before the first cycle and after the last cycle in the STS update This method sh...
subroutine, public set_error_handling_to_head(sts_error_handling)
Set the hook of error handling in the STS update. This method is called before updating the BC....
subroutine, public sts_init()
Initialize sts module.
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...
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_mhd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
anisotropic thermal conduction with slope limited symmetric scheme Sharma 2007 Journal of Computation...
subroutine get_euv_image(qunit, fl)
subroutine get_sxr_image(qunit, fl)
subroutine get_euv_spectrum(qunit, fl)
subroutine get_whitelight_image(qunit, fl)
Module with all the methods that users can customize in AMRVAC.
procedure(rfactor), pointer usr_rfactor
procedure(set_equi_vars), pointer usr_set_equi_vars