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_a2max => ffhd_get_a2max
314 phys_get_tcutoff => ffhd_get_tcutoff
315 phys_get_cbounds => ffhd_get_cbounds
316 phys_to_primitive => ffhd_to_primitive_origin
317 ffhd_to_primitive => ffhd_to_primitive_origin
318 phys_to_conserved => ffhd_to_conserved_origin
319 ffhd_to_conserved => ffhd_to_conserved_origin
320 phys_get_flux => ffhd_get_flux
321 phys_get_v => ffhd_get_v_origin
322 ffhd_get_v => ffhd_get_v_origin
323 phys_get_rho => ffhd_get_rho
324 ffhd_kin_en => ffhd_kin_en_origin
325 phys_add_source_geom => ffhd_add_source_geom
326 phys_add_source => ffhd_add_source
327 phys_check_params => ffhd_check_params
328 phys_write_info => ffhd_write_info
329 phys_handle_small_values => ffhd_handle_small_values_origin
330 ffhd_handle_small_values => ffhd_handle_small_values_origin
331 phys_check_w => ffhd_check_w_origin
332
333 if(.not.ffhd_energy) then
334 phys_get_pthermal => ffhd_get_pthermal_iso
335 ffhd_get_pthermal => ffhd_get_pthermal_iso
336 else
337 phys_get_pthermal => ffhd_get_pthermal_origin
338 ffhd_get_pthermal => ffhd_get_pthermal_origin
339 end if
340
341 ! choose Rfactor in ideal gas law
343 ffhd_get_rfactor=>rfactor_from_temperature_ionization
344 phys_update_temperature => ffhd_update_temperature
345 else if(associated(usr_rfactor)) then
346 ffhd_get_rfactor=>usr_rfactor
347 else
348 ffhd_get_rfactor=>rfactor_from_constant_ionization
349 end if
350
352 ffhd_get_temperature => ffhd_get_temperature_from_te
353 else
354 ffhd_get_temperature => ffhd_get_temperature_from_etot
355 end if
356
357 ! derive units from basic units
358 call ffhd_physical_units()
359
361 if(si_unit)then
362 ! parallel conduction Spitzer
364 else
365 ! in cgs
367 endif
368 end if
369 if(.not. ffhd_energy .and. ffhd_thermal_conduction) then
370 call mpistop("thermal conduction needs ffhd_energy=T")
371 end if
373 call mpistop("hyperbolic thermal conduction needs ffhd_energy=T")
374 end if
375 if(.not. ffhd_energy .and. ffhd_radiative_cooling) then
376 call mpistop("radiative cooling needs ffhd_energy=T")
377 end if
378
379 ! initialize thermal conduction module
381 call sts_init()
383
384 allocate(tc_fl)
385 call tc_get_hd_params(tc_fl,tc_params_read_ffhd)
386 call add_sts_method(ffhd_get_tc_dt_ffhd,ffhd_sts_set_source_tc_ffhd,e_,1,e_,1,.false.)
387 tc_fl%get_temperature_from_conserved => ffhd_get_temperature_from_etot
388 tc_fl%get_temperature_from_eint => ffhd_get_temperature_from_eint
390 call set_error_handling_to_head(ffhd_tc_handle_small_e)
391 tc_fl%get_rho => ffhd_get_rho
392 tc_fl%e_ = e_
393 tc_fl%Tcoff_ = tcoff_
394 end if
395
396 ! Initialize radiative cooling module
399 allocate(rc_fl)
400 call radiative_cooling_init(rc_fl,rc_params_read)
401 rc_fl%get_rho => ffhd_get_rho
402 rc_fl%get_pthermal => ffhd_get_pthermal
403 rc_fl%get_var_Rfactor => ffhd_get_rfactor
404 rc_fl%e_ = e_
405 rc_fl%Tcoff_ = tcoff_
406 rc_fl%has_equi = .false.
407 end if
408 allocate(te_fl_ffhd)
409 te_fl_ffhd%get_rho=> ffhd_get_rho
410 te_fl_ffhd%get_pthermal=> ffhd_get_pthermal
411 te_fl_ffhd%get_var_Rfactor => ffhd_get_rfactor
412{^ifthreed
413 phys_te_images => ffhd_te_images
414}
415
416 ! Initialize gravity module
417 if(ffhd_gravity) then
418 call gravity_init()
419 end if
420
421 ! initialize ionization degree table
423 end subroutine ffhd_phys_init
424
425{^ifthreed
426 subroutine ffhd_te_images
429
430 select case(convert_type)
431 case('EIvtiCCmpi','EIvtuCCmpi')
433 case('ESvtiCCmpi','ESvtuCCmpi')
435 case('SIvtiCCmpi','SIvtuCCmpi')
437 case('WIvtiCCmpi','WIvtuCCmpi')
439 case default
440 call mpistop("Error in synthesize emission: Unknown convert_type")
441 end select
442 end subroutine ffhd_te_images
443}
444
445 subroutine ffhd_sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
449 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
450 double precision, intent(in) :: x(ixi^s,1:ndim)
451 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
452 double precision, intent(in) :: my_dt
453 logical, intent(in) :: fix_conserve_at_step
454 call sts_set_source_tc_mhd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
455 end subroutine ffhd_sts_set_source_tc_ffhd
456
457 function ffhd_get_tc_dt_ffhd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
458 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
459 !where tc_k_para_i=tc_k_para*B_i**2/B**2
460 !and T=p/rho
463
464 integer, intent(in) :: ixi^l, ixo^l
465 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
466 double precision, intent(in) :: w(ixi^s,1:nw)
467 double precision :: dtnew
468
469 dtnew=get_tc_dt_mhd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
470 end function ffhd_get_tc_dt_ffhd
471
472 subroutine ffhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
474
475 integer, intent(in) :: ixi^l,ixo^l
476 double precision, intent(inout) :: w(ixi^s,1:nw)
477 double precision, intent(in) :: x(ixi^s,1:ndim)
478 integer, intent(in) :: step
479 character(len=140) :: error_msg
480
481 write(error_msg,"(a,i3)") "Thermal conduction step ", step
482 call ffhd_handle_small_ei(w,x,ixi^l,ixo^l,e_,error_msg)
483 end subroutine ffhd_tc_handle_small_e
484
485 subroutine tc_params_read_ffhd(fl)
487 type(tc_fluid), intent(inout) :: fl
488 integer :: n
489 ! list parameters
490 logical :: tc_saturate=.false.
491 double precision :: tc_k_para=0d0
492 character(len=std_len) :: tc_slope_limiter="MC"
493
494 namelist /tc_list/ tc_saturate, tc_slope_limiter, tc_k_para
495
496 do n = 1, size(par_files)
497 open(unitpar, file=trim(par_files(n)), status="old")
498 read(unitpar, tc_list, end=111)
499111 close(unitpar)
500 end do
501
502 fl%tc_saturate = tc_saturate
503 fl%tc_k_para = tc_k_para
504 select case(tc_slope_limiter)
505 case ('no','none')
506 fl%tc_slope_limiter = 0
507 case ('MC')
508 ! montonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
509 fl%tc_slope_limiter = 1
510 case('minmod')
511 ! minmod limiter
512 fl%tc_slope_limiter = 2
513 case ('superbee')
514 ! Roes superbee limiter (eq.3.51i)
515 fl%tc_slope_limiter = 3
516 case ('koren')
517 ! Barry Koren Right variant
518 fl%tc_slope_limiter = 4
519 case default
520 call mpistop("Unknown tc_slope_limiter, choose MC, minmod")
521 end select
522 end subroutine tc_params_read_ffhd
523
524 subroutine rc_params_read(fl)
526 use mod_constants, only: bigdouble
527 type(rc_fluid), intent(inout) :: fl
528 integer :: n
529 integer :: ncool = 4000
530 double precision :: cfrac=0.1d0
531
532 !> Name of cooling curve
533 character(len=std_len) :: coolcurve='JCcorona'
534
535 !> Name of cooling method
536 character(len=std_len) :: coolmethod='exact'
537
538 !> Fixed temperature not lower than tlow
539 logical :: tfix=.false.
540
541 !> Lower limit of temperature
542 double precision :: tlow=bigdouble
543
544 !> Add cooling source in a split way (.true.) or un-split way (.false.)
545 logical :: rc_split=.false.
546 logical :: rad_cut=.false.
547 double precision :: rad_cut_hgt=0.5d0
548 double precision :: rad_cut_dey=0.15d0
549
550 namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split, rad_cut, rad_cut_hgt, rad_cut_dey
551
552 do n = 1, size(par_files)
553 open(unitpar, file=trim(par_files(n)), status="old")
554 read(unitpar, rc_list, end=111)
555111 close(unitpar)
556 end do
557
558 fl%ncool=ncool
559 fl%coolcurve=coolcurve
560 fl%coolmethod=coolmethod
561 fl%tlow=tlow
562 fl%Tfix=tfix
563 fl%rc_split=rc_split
564 fl%cfrac=cfrac
565 fl%rad_cut=rad_cut
566 fl%rad_cut_hgt=rad_cut_hgt
567 fl%rad_cut_dey=rad_cut_dey
568 end subroutine rc_params_read
569
570 subroutine ffhd_check_params
574
575 gamma_1=ffhd_gamma-1.d0
576 if (.not. ffhd_energy) then
577 if (ffhd_gamma <= 0.0d0) call mpistop ("Error: ffhd_gamma <= 0")
578 if (ffhd_adiab < 0.0d0) call mpistop ("Error: ffhd_adiab < 0")
580 else
581 if (ffhd_gamma <= 0.0d0 .or. ffhd_gamma == 1.0d0) &
582 call mpistop ("Error: ffhd_gamma <= 0 or ffhd_gamma == 1")
583 inv_gamma_1=1.d0/gamma_1
584 small_e = small_pressure * inv_gamma_1
585 end if
586
587 if (number_equi_vars > 0 .and. .not. associated(usr_set_equi_vars)) then
588 call mpistop("usr_set_equi_vars has to be implemented in the user file")
589 end if
590 end subroutine ffhd_check_params
591
592 subroutine ffhd_physical_units()
594 double precision :: mp,kb
595 double precision :: a,b
596
597 if(si_unit) then
598 mp=mp_si
599 kb=kb_si
600 else
601 mp=mp_cgs
602 kb=kb_cgs
603 end if
604 if(eq_state_units) then
605 a=1d0+4d0*he_abundance
607 b=1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0)
608 else
609 b=2d0+3d0*he_abundance
610 end if
611 rr=1d0
612 else
613 a=1d0
614 b=1d0
615 rr=(1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0))/(1d0+4d0*he_abundance)
616 end if
617 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
618 if(unit_density/=1.d0) then
620 else if(unit_numberdensity/=1.d0) then
622 end if
623 if(unit_temperature/=1.d0) then
626 if(unit_length/=1.d0) then
628 else if(unit_time/=1.d0) then
630 end if
631 else if(unit_pressure/=1.d0) then
634 if(unit_length/=1.d0) then
636 else if(unit_time/=1.d0) then
638 end if
639 else if(unit_velocity/=1.d0) then
642 if(unit_length/=1.d0) then
644 else if(unit_time/=1.d0) then
646 end if
647 else if(unit_time/=1.d0) then
651 end if
652 else if(unit_temperature/=1.d0) then
653 ! units of temperature and velocity are dependent
654 if(unit_pressure/=1.d0) then
658 if(unit_length/=1.d0) then
660 else if(unit_time/=1.d0) then
662 end if
663 end if
664 else if(unit_pressure/=1.d0) then
665 if(unit_velocity/=1.d0) then
669 if(unit_length/=1.d0) then
671 else if(unit_time/=1.d0) then
673 end if
674 else if(unit_time/=0.d0) then
679 end if
680 end if
682 end subroutine ffhd_physical_units
683
684 subroutine ffhd_check_w_origin(primitive,ixI^L,ixO^L,w,flag)
686 logical, intent(in) :: primitive
687 integer, intent(in) :: ixi^l, ixo^l
688 double precision, intent(in) :: w(ixi^s,nw)
689 double precision :: tmp(ixi^s)
690 logical, intent(inout) :: flag(ixi^s,1:nw)
691
692 flag=.false.
693 where(w(ixo^s,rho_) < small_density) flag(ixo^s,rho_) = .true.
694
695 if(ffhd_energy) then
696 if(primitive) then
697 where(w(ixo^s,e_) < small_pressure) flag(ixo^s,e_) = .true.
698 else
699 tmp(ixo^s)=w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l)
700 where(tmp(ixo^s) < small_e) flag(ixo^s,e_) = .true.
701 end if
702 end if
703 end subroutine ffhd_check_w_origin
704
705 subroutine ffhd_to_conserved_origin(ixI^L,ixO^L,w,x)
707 integer, intent(in) :: ixi^l, ixo^l
708 double precision, intent(inout) :: w(ixi^s, nw)
709 double precision, intent(in) :: x(ixi^s, 1:ndim)
710
711 if(ffhd_energy) then
712 w(ixo^s,e_)=w(ixo^s,p_)*inv_gamma_1+half*w(ixo^s,mom(1))**2*w(ixo^s,rho_)
713 end if
714 w(ixo^s,mom(1))=w(ixo^s,rho_)*w(ixo^s,mom(1))
715 end subroutine ffhd_to_conserved_origin
716
717 subroutine ffhd_to_primitive_origin(ixI^L,ixO^L,w,x)
719 integer, intent(in) :: ixi^l, ixo^l
720 double precision, intent(inout) :: w(ixi^s, nw)
721 double precision, intent(in) :: x(ixi^s, 1:ndim)
722
723 if(fix_small_values) then
724 call ffhd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'ffhd_to_primitive_origin')
725 end if
726
727 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/w(ixo^s,rho_)
728 if(ffhd_energy) then
729 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-half*w(ixo^s,rho_)*w(ixo^s,mom(1))**2)
730 end if
731 end subroutine ffhd_to_primitive_origin
732
733 subroutine ffhd_ei_to_e(ixI^L,ixO^L,w,x)
735 integer, intent(in) :: ixi^l, ixo^l
736 double precision, intent(inout) :: w(ixi^s, nw)
737 double precision, intent(in) :: x(ixi^s, 1:ndim)
738
739 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
740 end subroutine ffhd_ei_to_e
741
742 subroutine ffhd_e_to_ei(ixI^L,ixO^L,w,x)
744 integer, intent(in) :: ixi^l, ixo^l
745 double precision, intent(inout) :: w(ixi^s, nw)
746 double precision, intent(in) :: x(ixi^s, 1:ndim)
747
748 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
749 if(fix_small_values) then
750 call ffhd_handle_small_ei(w,x,ixi^l,ixi^l,e_,'ffhd_e_to_ei')
751 end if
752 end subroutine ffhd_e_to_ei
753
754 subroutine ffhd_handle_small_values_origin(primitive, w, x, ixI^L, ixO^L, subname)
757 logical, intent(in) :: primitive
758 integer, intent(in) :: ixi^l,ixo^l
759 double precision, intent(inout) :: w(ixi^s,1:nw)
760 double precision, intent(in) :: x(ixi^s,1:ndim)
761 character(len=*), intent(in) :: subname
762
763 logical :: flag(ixi^s,1:nw)
764 double precision :: tmp2(ixi^s)
765
766 call phys_check_w(primitive, ixi^l, ixi^l, w, flag)
767
768 if(any(flag)) then
769 select case (small_values_method)
770 case ("replace")
771 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
772 if(small_values_fix_iw(mom(1))) then
773 where(flag(ixo^s,rho_)) w(ixo^s, mom(1)) = 0.0d0
774 end if
775 if(ffhd_energy) then
776 if(primitive) then
777 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
778 else
779 where(flag(ixo^s,e_))
780 w(ixo^s,e_) = small_e+ffhd_kin_en(w,ixi^l,ixo^l)
781 end where
782 end if
783 end if
784 case ("average")
785 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
786 if(ffhd_energy) then
787 if(primitive) then
788 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
789 else
790 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
791 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
792 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
793 end if
794 end if
795 case default
796 if(.not.primitive) then
797 if(ffhd_energy) then
798 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
799 end if
800 w(ixo^s,mom(1))=w(ixo^s,mom(1))/w(ixo^s,rho_)
801 end if
802 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
803 end select
804 end if
805 end subroutine ffhd_handle_small_values_origin
806
807 subroutine ffhd_get_v_origin(w,x,ixI^L,ixO^L,v)
809 integer, intent(in) :: ixi^l, ixo^l
810 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
811 double precision, intent(out) :: v(ixi^s,ndir)
812 double precision :: rho(ixi^s)
813 integer :: idir
814
815 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
816 rho(ixo^s)=1.d0/rho(ixo^s)
817 do idir=1,ndir
818 v(ixo^s,ndir) = w(ixo^s,mom(1))*block%B0(ixo^s,idir,0)*rho(ixo^s)
819 end do
820 end subroutine ffhd_get_v_origin
821
822 subroutine ffhd_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
824 integer, intent(in) :: ixi^l, ixo^l, idim
825 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
826 double precision, intent(out) :: v(ixi^s)
827 double precision :: rho(ixi^s)
828
829 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
830 v(ixo^s) = (w(ixo^s, mom(1))*block%B0(ixo^s,idim,0)) / rho(ixo^s)
831 end subroutine ffhd_get_v_idim
832
833 subroutine ffhd_get_cmax_origin(wprim,x,ixI^L,ixO^L,idim,cmax)
835 integer, intent(in) :: ixi^l, ixo^l, idim
836 ! w in primitive form
837 double precision, intent(in) :: wprim(ixi^s, nw), x(ixi^s,1:ndim)
838 double precision, intent(inout) :: cmax(ixi^s)
839
840 if(ffhd_energy) then
841 cmax(ixo^s)=dsqrt(ffhd_gamma*wprim(ixo^s,p_)/wprim(ixo^s,rho_))
842 else
843 cmax(ixo^s)=dsqrt(ffhd_gamma*ffhd_adiab*wprim(ixo^s,rho_)**gamma_1)
844 end if
845 cmax(ixo^s)=dabs(wprim(ixo^s,mom(1))*block%B0(ixo^s,idim,0))+cmax(ixo^s)
846
847 end subroutine ffhd_get_cmax_origin
848
849 subroutine ffhd_get_a2max(w,x,ixI^L,ixO^L,a2max)
851 use mod_geometry
852 integer, intent(in) :: ixi^l, ixo^l
853 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
854 double precision, intent(inout) :: a2max(ndim)
855 double precision :: a2(ixi^s,ndim,nw)
856 integer :: gxo^l,hxo^l,jxo^l,kxo^l,i,j
857
858 if(.not.slab_uniform)then
859 call mpistop("subroutine get_a2max in mod_ffhd_phys adopts cartesian setting")
860 endif
861 a2=zero
862 do i = 1,ndim
863 !> 4th order
864 hxo^l=ixo^l-kr(i,^d);
865 gxo^l=hxo^l-kr(i,^d);
866 jxo^l=ixo^l+kr(i,^d);
867 kxo^l=jxo^l+kr(i,^d);
868 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
869 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
870 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/dxlevel(i)**2
871 end do
872 end subroutine ffhd_get_a2max
873
874 subroutine ffhd_get_tcutoff(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
876 use mod_geometry
877 integer, intent(in) :: ixi^l,ixo^l
878 double precision, intent(in) :: x(ixi^s,1:ndim)
879 double precision, intent(inout) :: w(ixi^s,1:nw)
880 double precision, intent(out) :: tco_local,tmax_local
881 double precision, parameter :: trac_delta=0.25d0
882 double precision :: te(ixi^s),lts(ixi^s)
883 double precision, dimension(ixI^S,1:ndim) :: gradt
884 double precision :: bdir(ndim)
885 double precision :: ltrc,ltrp,altr
886 integer :: idims,ix^d,jxo^l,hxo^l,ixa^d,ixb^d
887 integer :: jxp^l,hxp^l,ixp^l,ixq^l
888 logical :: lrlt(ixi^s)
889
890 call ffhd_get_temperature(w,x,ixi^l,ixi^l,te)
891 tco_local=zero
892 tmax_local=maxval(te(ixo^s))
893
894 {^ifoned
895 select case(ffhd_trac_type)
896 case(0)
897 !> test case, fixed cutoff temperature
898 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
899 case(1)
900 hxo^l=ixo^l-1;
901 jxo^l=ixo^l+1;
902 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
903 lrlt=.false.
904 where(lts(ixo^s) > trac_delta)
905 lrlt(ixo^s)=.true.
906 end where
907 if(any(lrlt(ixo^s))) then
908 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
909 end if
910 case(2)
911 !> iijima et al. 2021, LTRAC method
912 ltrc=1.5d0
913 ltrp=4.d0
914 ixp^l=ixo^l^ladd1;
915 hxo^l=ixo^l-1;
916 jxo^l=ixo^l+1;
917 hxp^l=ixp^l-1;
918 jxp^l=ixp^l+1;
919 lts(ixp^s)=0.5d0*dabs(te(jxp^s)-te(hxp^s))/te(ixp^s)
920 lts(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
921 lts(ixo^s)=0.25d0*(lts(jxo^s)+two*lts(ixo^s)+lts(hxo^s))
922 block%wextra(ixo^s,tcoff_)=te(ixo^s)*lts(ixo^s)**0.4d0
923 case default
924 call mpistop("ffhd_trac_type not allowed for 1D simulation")
925 end select
926 }
927 {^nooned
928 select case(ffhd_trac_type)
929 case(0)
930 !> test case, fixed cutoff temperature
931 if(slab_uniform) then
932 !> assume cgs units
933 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)
934 else
935 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
936 end if
937 case(1,4,6)
938 do idims=1,ndim
939 call gradient(te,ixi^l,ixo^l,idims,gradt(ixi^s,idims))
940 end do
941 if(ffhd_trac_type .gt. 1) then
942 ! B direction at cell center
943 bdir=zero
944 {do ixa^d=0,1\}
945 ixb^d=(ixomin^d+ixomax^d-1)/2+ixa^d;
946 bdir(1:ndim)=bdir(1:ndim)+block%B0(ixb^d,1:ndim,0)
947 {end do\}
948 if(sum(bdir(:)**2) .gt. zero) then
949 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
950 end if
951 block%special_values(3:ndim+2)=bdir(1:ndim)
952 end if
953 {do ix^db=ixomin^db,ixomax^db\}
954 {^iftwod
955 ! temperature length scale inversed
956 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2))*&
957 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
958 }
959 {^ifthreed
960 ! temperature length scale inversed
961 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2),block%ds(ix^d,3))*&
962 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
963 }
964 if(lts(ix^d)>trac_delta) then
965 block%special_values(1)=max(block%special_values(1),te(ix^d))
966 end if
967 {end do\}
968 block%special_values(2)=tmax_local
969 case(2)
970 !> iijima et al. 2021, LTRAC method
971 ltrc=1.5d0
972 ltrp=4.d0
973 ixp^l=ixo^l^ladd2;
974 do idims=1,ndim
975 ixq^l=ixp^l;
976 hxp^l=ixp^l;
977 jxp^l=ixp^l;
978 select case(idims)
979 {case(^d)
980 ixqmin^d=ixqmin^d+1
981 ixqmax^d=ixqmax^d-1
982 hxpmax^d=ixpmin^d
983 jxpmin^d=ixpmax^d
984 \}
985 end select
986 call gradient(te,ixi^l,ixq^l,idims,gradt(ixi^s,idims))
987 call gradientf(te,x,ixi^l,hxp^l,idims,gradt(ixi^s,idims),nghostcells,.true.)
988 call gradientf(te,x,ixi^l,jxp^l,idims,gradt(ixi^s,idims),nghostcells,.false.)
989 end do
990 {do ix^db=ixpmin^db,ixpmax^db\}
991 ! temperature length scale inversed
992 lts(ix^d)=abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
993 ! fraction of cells size to temperature length scale
994 lts(ix^d)=min(^d&block%ds({ix^d},^d))*lts(ix^d)
995 lts(ix^d)=max(one,(exp(lts(ix^d))/ltrc)**ltrp)
996 {end do\}
997
998 ! need one ghost layer for thermal conductivity
999 ixp^l=ixo^l^ladd1;
1000 {do ix^db=ixpmin^db,ixpmax^db\}
1001 {^iftwod
1002 altr=0.25d0*((lts(ix1-1,ix2)+two*lts(ix^d)+lts(ix1+1,ix2))*block%B0(ix^d,1,0)**2+&
1003 (lts(ix1,ix2-1)+two*lts(ix^d)+lts(ix1,ix2+1))*block%B0(ix^d,2,0)**2)
1004 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
1005 }
1006 {^ifthreed
1007 altr=0.25d0*((lts(ix1-1,ix2,ix3)+two*lts(ix^d)+lts(ix1+1,ix2,ix3))*block%B0(ix^d,1,0)**2+&
1008 (lts(ix1,ix2-1,ix3)+two*lts(ix^d)+lts(ix1,ix2+1,ix3))*block%B0(ix^d,2,0)**2+&
1009 (lts(ix1,ix2,ix3-1)+two*lts(ix^d)+lts(ix1,ix2,ix3+1))*block%B0(ix^d,3,0)**2)
1010 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
1011 }
1012 {end do\}
1013 case(3,5)
1014 !> do nothing here
1015 case default
1016 call mpistop("unknown ffhd_trac_type")
1017 end select
1018 }
1019 end subroutine ffhd_get_tcutoff
1020
1021 subroutine ffhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
1023 integer, intent(in) :: ixi^l, ixo^l, idim
1024 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1025 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1026 double precision, intent(in) :: x(ixi^s,1:ndim)
1027 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1028 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1029 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
1030 double precision :: wmean(ixi^s,nw)
1031 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1032
1033 select case (boundspeed)
1034 case (1)
1035 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1036 ! Methods for Fluid Dynamics" by Toro.
1037 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1038 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1039 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1040 umean(ixo^s)=(wlp(ixo^s,mom(1))*tmp1(ixo^s)&
1041 +wrp(ixo^s,mom(1))*tmp2(ixo^s))*tmp3(ixo^s)
1042 umean(ixo^s)=umean(ixo^s)*block%B0(ixo^s,idim,idim)
1043 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1044 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1045 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1046 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1047 ((wrp(ixo^s,mom(1))-wlp(ixo^s,mom(1)))*block%B0(ixo^s,idim,idim))**2
1048 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1049 if(present(cmin)) then
1050 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1051 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1052 else
1053 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1054 end if
1055 case (2)
1056 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1057 tmp1(ixo^s)=wmean(ixo^s,mom(1))*block%B0(ixo^s,idim,idim)/wmean(ixo^s,rho_)
1058 call ffhd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1059 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1060 if(present(cmin)) then
1061 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1062 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1063 else
1064 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1065 end if
1066 case (3)
1067 ! Miyoshi 2005 JCP 208, 315 equation (67)
1068 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1069 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1070 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1071 if(present(cmin)) then
1072 cmin(ixo^s,1)=min(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1073 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))-csoundl(ixo^s)
1074 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1075 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1076 else
1077 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1078 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1079 end if
1080 end select
1081 end subroutine ffhd_get_cbounds
1082
1083 subroutine ffhd_get_pthermal_iso(w,x,ixI^L,ixO^L,pth)
1085
1086 integer, intent(in) :: ixi^l, ixo^l
1087 double precision, intent(in) :: w(ixi^s,nw)
1088 double precision, intent(in) :: x(ixi^s,1:ndim)
1089 double precision, intent(out):: pth(ixi^s)
1090
1091 call ffhd_get_rho(w,x,ixi^l,ixo^l,pth)
1092 pth(ixo^s)=ffhd_adiab*pth(ixo^s)**ffhd_gamma
1093 end subroutine ffhd_get_pthermal_iso
1094
1095 subroutine ffhd_get_pthermal_origin(w,x,ixI^L,ixO^L,pth)
1098 integer, intent(in) :: ixi^l, ixo^l
1099 double precision, intent(in) :: w(ixi^s,nw)
1100 double precision, intent(in) :: x(ixi^s,1:ndim)
1101 double precision, intent(out):: pth(ixi^s)
1102 integer :: iw, ix^d
1103
1104 pth(ixo^s)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
1105 if (fix_small_values) then
1106 {do ix^db= ixo^lim^db\}
1107 if(pth(ix^d)<small_pressure) then
1108 pth(ix^d)=small_pressure
1109 end if
1110 {end do^D&\}
1111 elseif(check_small_values) then
1112 {do ix^db= ixo^lim^db\}
1113 if(pth(ix^d)<small_pressure) then
1114 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1115 " encountered when call ffhd_get_pthermal"
1116 write(*,*) "Iteration: ", it, " Time: ", global_time
1117 write(*,*) "Location: ", x(ix^d,:)
1118 write(*,*) "Cell number: ", ix^d
1119 do iw=1,nw
1120 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1121 end do
1122 ! use erroneous arithmetic operation to crash the run
1123 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1124 write(*,*) "Saving status at the previous time step"
1125 crash=.true.
1126 end if
1127 {end do^D&\}
1128 end if
1129 end subroutine ffhd_get_pthermal_origin
1130
1131 subroutine ffhd_get_temperature_from_te(w, x, ixI^L, ixO^L, res)
1133 integer, intent(in) :: ixi^l, ixo^l
1134 double precision, intent(in) :: w(ixi^s, 1:nw)
1135 double precision, intent(in) :: x(ixi^s, 1:ndim)
1136 double precision, intent(out):: res(ixi^s)
1137
1138 res(ixo^s) = w(ixo^s, te_)
1139 end subroutine ffhd_get_temperature_from_te
1140
1141 subroutine ffhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1143 integer, intent(in) :: ixi^l, ixo^l
1144 double precision, intent(in) :: w(ixi^s, 1:nw)
1145 double precision, intent(in) :: x(ixi^s, 1:ndim)
1146 double precision, intent(out):: res(ixi^s)
1147 double precision :: r(ixi^s)
1148
1149 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1150 res(ixo^s) = gamma_1 * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1151 end subroutine ffhd_get_temperature_from_eint
1152
1153 subroutine ffhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1155 integer, intent(in) :: ixi^l, ixo^l
1156 double precision, intent(in) :: w(ixi^s, 1:nw)
1157 double precision, intent(in) :: x(ixi^s, 1:ndim)
1158 double precision, intent(out):: res(ixi^s)
1159
1160 double precision :: r(ixi^s)
1161
1162 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1163 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,res)
1164 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1165 end subroutine ffhd_get_temperature_from_etot
1166
1167 subroutine ffhd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1169 integer, intent(in) :: ixi^l, ixo^l
1170 double precision, intent(in) :: w(ixi^s,nw)
1171 double precision, intent(in) :: x(ixi^s,1:ndim)
1172 double precision, intent(out) :: csound2(ixi^s)
1173 double precision :: rho(ixi^s)
1174
1175 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1176 if(ffhd_energy) then
1177 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1178 csound2(ixo^s)=ffhd_gamma*csound2(ixo^s)/rho(ixo^s)
1179 else
1180 csound2(ixo^s)=ffhd_gamma*ffhd_adiab*rho(ixo^s)**gamma_1
1181 end if
1182 end subroutine ffhd_get_csound2
1183
1184 subroutine ffhd_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
1186 use mod_geometry
1187 integer, intent(in) :: ixi^l, ixo^l, idim
1188 ! conservative w
1189 double precision, intent(in) :: wc(ixi^s,nw)
1190 ! primitive w
1191 double precision, intent(in) :: w(ixi^s,nw)
1192 double precision, intent(in) :: x(ixi^s,1:ndim)
1193 double precision,intent(out) :: f(ixi^s,nwflux)
1194 double precision :: ptotal(ixo^s)
1195
1196 f(ixo^s,rho_)=w(ixo^s,mom(1))*w(ixo^s,rho_)*block%B0(ixo^s,idim,idim)
1197
1198 if(ffhd_energy) then
1199 ptotal(ixo^s)=w(ixo^s,p_)
1200 else
1201 ptotal(ixo^s)=ffhd_adiab*w(ixo^s,rho_)**ffhd_gamma
1202 end if
1203
1204 ! Get flux of momentum
1205 f(ixo^s,mom(1))=(wc(ixo^s,mom(1))*w(ixo^s,mom(1))+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1206
1207 ! Get flux of energy
1208 if(ffhd_energy) then
1209 f(ixo^s,e_)=w(ixo^s,mom(1))*(wc(ixo^s,e_)+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1211 f(ixo^s,e_)=f(ixo^s,e_)+w(ixo^s,q_)*block%B0(ixo^s,idim,idim)
1212 f(ixo^s,q_)=zero
1213 end if
1214 end if
1215 end subroutine ffhd_get_flux
1216
1217 subroutine ffhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1221 integer, intent(in) :: ixi^l, ixo^l
1222 double precision, intent(in) :: qdt,dtfactor
1223 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:ndim)
1224 double precision, intent(inout) :: w(ixi^s,1:nw)
1225 logical, intent(in) :: qsourcesplit
1226 logical, intent(inout) :: active
1227
1228 if (.not. qsourcesplit) then
1229 active = .true.
1230 call add_punitb(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1232 call add_hypertc_source(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1233 end if
1234 end if
1235
1236 if(ffhd_radiative_cooling) then
1237 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1238 w,x,qsourcesplit,active, rc_fl)
1239 end if
1240
1241 if(ffhd_gravity) then
1242 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1243 w,x,ffhd_energy,qsourcesplit,active)
1244 end if
1245
1246 ! update temperature from new pressure, density, and old ionization degree
1248 if(.not.qsourcesplit) then
1249 active = .true.
1250 call ffhd_update_temperature(ixi^l,ixo^l,wct,w,x)
1251 end if
1252 end if
1253 end subroutine ffhd_add_source
1254
1255 subroutine add_punitb(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1257 use mod_geometry
1258 integer, intent(in) :: ixi^l,ixo^l
1259 double precision, intent(in) :: qdt
1260 double precision, intent(in) :: wct(ixi^s,1:nw),x(ixi^s,1:ndim)
1261 double precision, intent(in) :: wctprim(ixi^s,1:nw)
1262 double precision, intent(inout) :: w(ixi^s,1:nw)
1263
1264 integer :: idims,hxo^l
1265 double precision :: divb(ixi^s)
1266 double precision :: rhovpar(ixi^s),gradrhov(ixi^s)
1267
1268 divb=zero
1269 if(slab_uniform) then
1270 do idims=1,ndim
1271 hxo^l=ixo^l-kr(idims,^d);
1272 divb(ixo^s)=divb(ixo^s)+(block%B0(ixo^s,idims,idims)-block%B0(hxo^s,idims,idims))/dxlevel(idims)
1273 end do
1274 else
1275 call divvector(block%B0(ixi^s,1:ndir,0),ixi^l,ixo^l,divb)
1276 end if
1277 w(ixo^s,mom(1))=w(ixo^s,mom(1))+qdt*wctprim(ixo^s,p_)*divb(ixo^s)
1278 end subroutine add_punitb
1279
1280 subroutine ffhd_get_rho(w,x,ixI^L,ixO^L,rho)
1282 integer, intent(in) :: ixi^l, ixo^l
1283 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
1284 double precision, intent(out) :: rho(ixi^s)
1285
1286 rho(ixo^s) = w(ixo^s,rho_)
1287 end subroutine ffhd_get_rho
1288
1289 subroutine ffhd_handle_small_ei(w, x, ixI^L, ixO^L, ie, subname)
1292 integer, intent(in) :: ixi^l,ixo^l, ie
1293 double precision, intent(inout) :: w(ixi^s,1:nw)
1294 double precision, intent(in) :: x(ixi^s,1:ndim)
1295 character(len=*), intent(in) :: subname
1296 integer :: idir
1297 logical :: flag(ixi^s,1:nw)
1298 double precision :: rho(ixi^s)
1299
1300 flag=.false.
1301 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
1302 if(any(flag(ixo^s,ie))) then
1303 select case (small_values_method)
1304 case ("replace")
1305 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
1306 case ("average")
1307 call small_values_average(ixi^l, ixo^l, w, x, flag, ie)
1308 case default
1309 w(ixo^s,e_)=w(ixo^s,e_)*gamma_1
1310 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1311 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/rho(ixo^s)
1312 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1313 end select
1314 end if
1315 end subroutine ffhd_handle_small_ei
1316
1317 subroutine ffhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1320 integer, intent(in) :: ixi^l, ixo^l
1321 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1322 double precision, intent(inout) :: w(ixi^s,1:nw)
1323 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1324
1325 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1326 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1327 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1328 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1329 end subroutine ffhd_update_temperature
1330
1331 subroutine ffhd_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
1333 use mod_usr_methods
1335 use mod_gravity, only: gravity_get_dt
1336 integer, intent(in) :: ixi^l, ixo^l
1337 double precision, intent(inout) :: dtnew
1338 double precision, intent(in) :: dx^d
1339 double precision, intent(in) :: w(ixi^s,1:nw)
1340 double precision, intent(in) :: x(ixi^s,1:ndim)
1341
1342 dtnew = bigdouble
1343
1344 if(ffhd_radiative_cooling) then
1345 call cooling_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x,rc_fl)
1346 end if
1347
1348 if(ffhd_gravity) then
1349 call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1350 end if
1351 end subroutine ffhd_get_dt
1352
1353 subroutine ffhd_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
1355 integer, intent(in) :: ixi^l, ixo^l
1356 double precision, intent(in) :: qdt, dtfactor,x(ixi^s,1:ndim)
1357 double precision, intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw),w(ixi^s,1:nw)
1358
1359 ! no geometric source terms needed for ffhd, no divergences of tensors
1360 end subroutine ffhd_add_source_geom
1361
1362 function ffhd_kin_en_origin(w, ixI^L, ixO^L, inv_rho) result(ke)
1363 use mod_global_parameters, only: nw, ndim,block
1364 integer, intent(in) :: ixi^l, ixo^l
1365 double precision, intent(in) :: w(ixi^s, nw)
1366 double precision :: ke(ixo^s)
1367 double precision, intent(in), optional :: inv_rho(ixo^s)
1368
1369 if(present(inv_rho)) then
1370 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2*inv_rho(ixo^s)
1371 else
1372 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2/w(ixo^s,rho_)
1373 end if
1374 end function ffhd_kin_en_origin
1375
1376 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1379 integer, intent(in) :: ixi^l, ixo^l
1380 double precision, intent(in) :: w(ixi^s,1:nw)
1381 double precision, intent(in) :: x(ixi^s,1:ndim)
1382 double precision, intent(out):: rfactor(ixi^s)
1383 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1384
1385 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1386 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)
1387 end subroutine rfactor_from_temperature_ionization
1388
1389 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1391 integer, intent(in) :: ixi^l, ixo^l
1392 double precision, intent(in) :: w(ixi^s,1:nw)
1393 double precision, intent(in) :: x(ixi^s,1:ndim)
1394 double precision, intent(out):: rfactor(ixi^s)
1395
1396 rfactor(ixo^s)=rr
1397 end subroutine rfactor_from_constant_ionization
1398
1399 subroutine add_hypertc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1401 use mod_geometry
1402 integer, intent(in) :: ixi^l,ixo^l
1403 double precision, intent(in) :: qdt
1404 double precision, dimension(ixI^S,1:ndim), intent(in) :: x
1405 double precision, dimension(ixI^S,1:nw), intent(in) :: wct,wctprim
1406 double precision, dimension(ixI^S,1:nw), intent(inout) :: w
1407
1408 double precision, dimension(ixI^S) :: te,r,bgradt,gradt
1409 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1410 integer :: ix^d,idims
1411
1412 call ffhd_get_rfactor(wct,x,ixi^l,ixi^l,r)
1413 te(ixi^s)=wctprim(ixi^s,p_)/(r(ixi^s)*wct(ixi^s,rho_))
1414 {^ifoned
1415 call gradient(te,ixi^l,ixo^l,1,bgradt,2)
1416 }
1417 {^nooned
1418 bgradt(ixo^s)=zero
1419 do idims=1,ndim
1420 ! compute gradient conform the geometry, 4th order CD for uniform cartesian by setting 2
1421 call gradient(te,ixi^l,ixo^l,idims,gradt,2)
1422 bgradt(ixo^s)=bgradt(ixo^s)+(block%B0(ixo^s,idims,0))*gradt(ixo^s)
1423 enddo
1424 }
1425 {do ix^db=ixomin^db,ixomax^db\}
1426 if(ffhd_trac) then
1427 r(ix^d)=max(te(ix^d),block%wextra(ix^d,tcoff_))
1428 else
1429 r(ix^d)=te(ix^d)
1430 endif
1431 sigma_t5=hypertc_kappa*dsqrt(r(ix^d)**5)
1432 sigma_t7=sigma_t5*r(ix^d)
1433 sigmat5_bgradt=sigma_t5*bgradt(ix^d)
1434 if(ffhd_htc_sat) then
1435 f_sat=one/(one+dabs(sigmat5_bgradt)/(1.5d0*wct(ix^d,rho_)*(ffhd_gamma*wctprim(ix^d,p_)/wct(ix^d,rho_))**1.5d0))
1436 tau=max(4.d0*dt, f_sat*sigma_t7/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1437 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1438 else
1439 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1440 max(4.d0*dt, sigma_t7/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1441 end if
1442 {end do\}
1443
1444
1445 end subroutine add_hypertc_source
1446
1447end 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 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 for HD and MHD
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine cooling_get_dt(w, ixil, ixol, dtnew, dxd, x, fl)
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 1) in amrvac.par in sts_list set the following parameters which have...
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