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 rc_fl%subtract_equi = .false.
408 end if
409
410{^ifthreed
411 allocate(te_fl_ffhd)
412 te_fl_ffhd%get_rho=> ffhd_get_rho
413 te_fl_ffhd%get_pthermal=> ffhd_get_pthermal
414 te_fl_ffhd%get_var_Rfactor => ffhd_get_rfactor
415 phys_te_images => ffhd_te_images
416}
417
418 ! Initialize gravity module
419 if(ffhd_gravity) then
420 call gravity_init()
421 end if
422
423 ! initialize ionization degree table
425 end subroutine ffhd_phys_init
426
427{^ifthreed
428 subroutine ffhd_te_images
431
432 select case(convert_type)
433 case('EIvtiCCmpi','EIvtuCCmpi')
435 case('ESvtiCCmpi','ESvtuCCmpi')
437 case('SIvtiCCmpi','SIvtuCCmpi')
439 case('WIvtiCCmpi','WIvtuCCmpi')
441 case default
442 call mpistop("Error in synthesize emission: Unknown convert_type")
443 end select
444 end subroutine ffhd_te_images
445}
446
447 subroutine ffhd_sts_set_source_tc_ffhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
451 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
452 double precision, intent(in) :: x(ixi^s,1:ndim)
453 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
454 double precision, intent(in) :: my_dt
455 logical, intent(in) :: fix_conserve_at_step
456 call sts_set_source_tc_mhd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
457 end subroutine ffhd_sts_set_source_tc_ffhd
458
459 function ffhd_get_tc_dt_ffhd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
460 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
461 !where tc_k_para_i=tc_k_para*B_i**2/B**2
462 !and T=p/rho
465
466 integer, intent(in) :: ixi^l, ixo^l
467 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
468 double precision, intent(in) :: w(ixi^s,1:nw)
469 double precision :: dtnew
470
471 dtnew=get_tc_dt_mhd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
472 end function ffhd_get_tc_dt_ffhd
473
474 subroutine ffhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
476
477 integer, intent(in) :: ixi^l,ixo^l
478 double precision, intent(inout) :: w(ixi^s,1:nw)
479 double precision, intent(in) :: x(ixi^s,1:ndim)
480 integer, intent(in) :: step
481 character(len=140) :: error_msg
482
483 write(error_msg,"(a,i3)") "Thermal conduction step ", step
484 call ffhd_handle_small_ei(w,x,ixi^l,ixo^l,e_,error_msg)
485 end subroutine ffhd_tc_handle_small_e
486
487 subroutine tc_params_read_ffhd(fl)
489 type(tc_fluid), intent(inout) :: fl
490 integer :: n
491 ! list parameters
492 logical :: tc_saturate=.false.
493 double precision :: tc_k_para=0d0
494 character(len=std_len) :: tc_slope_limiter="MC"
495
496 namelist /tc_list/ tc_saturate, tc_slope_limiter, tc_k_para
497
498 do n = 1, size(par_files)
499 open(unitpar, file=trim(par_files(n)), status="old")
500 read(unitpar, tc_list, end=111)
501111 close(unitpar)
502 end do
503
504 fl%tc_saturate = tc_saturate
505 fl%tc_k_para = tc_k_para
506 select case(tc_slope_limiter)
507 case ('no','none')
508 fl%tc_slope_limiter = 0
509 case ('MC')
510 ! montonized central limiter Woodward and Collela limiter (eq.3.51h), a factor of 2 is pulled out
511 fl%tc_slope_limiter = 1
512 case('minmod')
513 ! minmod limiter
514 fl%tc_slope_limiter = 2
515 case ('superbee')
516 ! Roes superbee limiter (eq.3.51i)
517 fl%tc_slope_limiter = 3
518 case ('koren')
519 ! Barry Koren Right variant
520 fl%tc_slope_limiter = 4
521 case default
522 call mpistop("Unknown tc_slope_limiter, choose MC, minmod")
523 end select
524 end subroutine tc_params_read_ffhd
525
526 subroutine rc_params_read(fl)
528 use mod_constants, only: bigdouble
529 type(rc_fluid), intent(inout) :: fl
530 integer :: n
531 integer :: ncool = 4000
532
533 !> Name of cooling curve
534 character(len=std_len) :: coolcurve='JCcorona'
535
536 !> Fixed temperature not lower than tlow
537 logical :: tfix=.false.
538
539 !> Lower limit of temperature
540 double precision :: tlow=bigdouble
541
542 !> Add cooling source in a split way (.true.) or un-split way (.false.)
543 logical :: rc_split=.false.
544 logical :: rad_cut=.false.
545 double precision :: rad_cut_hgt=0.5d0
546 double precision :: rad_cut_dey=0.15d0
547
548 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split, rad_cut, rad_cut_hgt, rad_cut_dey
549
550 do n = 1, size(par_files)
551 open(unitpar, file=trim(par_files(n)), status="old")
552 read(unitpar, rc_list, end=111)
553111 close(unitpar)
554 end do
555
556 fl%ncool=ncool
557 fl%coolcurve=coolcurve
558 fl%tlow=tlow
559 fl%Tfix=tfix
560 fl%rc_split=rc_split
561 fl%rad_cut=rad_cut
562 fl%rad_cut_hgt=rad_cut_hgt
563 fl%rad_cut_dey=rad_cut_dey
564 end subroutine rc_params_read
565
566 subroutine ffhd_check_params
570
571 gamma_1=ffhd_gamma-1.d0
572 if (.not. ffhd_energy) then
573 if (ffhd_gamma <= 0.0d0) call mpistop ("Error: ffhd_gamma <= 0")
574 if (ffhd_adiab < 0.0d0) call mpistop ("Error: ffhd_adiab < 0")
576 else
577 if (ffhd_gamma <= 0.0d0 .or. ffhd_gamma == 1.0d0) &
578 call mpistop ("Error: ffhd_gamma <= 0 or ffhd_gamma == 1")
579 inv_gamma_1=1.d0/gamma_1
580 small_e = small_pressure * inv_gamma_1
581 end if
582
583 if (number_equi_vars > 0 .and. .not. associated(usr_set_equi_vars)) then
584 call mpistop("usr_set_equi_vars has to be implemented in the user file")
585 end if
586 end subroutine ffhd_check_params
587
588 subroutine ffhd_physical_units()
590 double precision :: mp,kb
591 double precision :: a,b
592
593 if(si_unit) then
594 mp=mp_si
595 kb=kb_si
596 else
597 mp=mp_cgs
598 kb=kb_cgs
599 end if
600 if(eq_state_units) then
601 a=1d0+4d0*he_abundance
603 b=1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0)
604 else
605 b=2d0+3d0*he_abundance
606 end if
607 rr=1d0
608 else
609 a=1d0
610 b=1d0
611 rr=(1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0))/(1d0+4d0*he_abundance)
612 end if
613 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
614 if(unit_density/=1.d0) then
616 else if(unit_numberdensity/=1.d0) then
618 end if
619 if(unit_temperature/=1.d0) then
622 if(unit_length/=1.d0) then
624 else if(unit_time/=1.d0) then
626 end if
627 else if(unit_pressure/=1.d0) then
630 if(unit_length/=1.d0) then
632 else if(unit_time/=1.d0) then
634 end if
635 else if(unit_velocity/=1.d0) then
638 if(unit_length/=1.d0) then
640 else if(unit_time/=1.d0) then
642 end if
643 else if(unit_time/=1.d0) then
647 end if
648 else if(unit_temperature/=1.d0) then
649 ! units of temperature and velocity are dependent
650 if(unit_pressure/=1.d0) then
654 if(unit_length/=1.d0) then
656 else if(unit_time/=1.d0) then
658 end if
659 end if
660 else if(unit_pressure/=1.d0) then
661 if(unit_velocity/=1.d0) then
665 if(unit_length/=1.d0) then
667 else if(unit_time/=1.d0) then
669 end if
670 else if(unit_time/=0.d0) then
675 end if
676 end if
678 end subroutine ffhd_physical_units
679
680 subroutine ffhd_check_w_origin(primitive,ixI^L,ixO^L,w,flag)
682 logical, intent(in) :: primitive
683 integer, intent(in) :: ixi^l, ixo^l
684 double precision, intent(in) :: w(ixi^s,nw)
685 double precision :: tmp(ixi^s)
686 logical, intent(inout) :: flag(ixi^s,1:nw)
687
688 flag=.false.
689 where(w(ixo^s,rho_) < small_density) flag(ixo^s,rho_) = .true.
690
691 if(ffhd_energy) then
692 if(primitive) then
693 where(w(ixo^s,e_) < small_pressure) flag(ixo^s,e_) = .true.
694 else
695 tmp(ixo^s)=w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l)
696 where(tmp(ixo^s) < small_e) flag(ixo^s,e_) = .true.
697 end if
698 end if
699 end subroutine ffhd_check_w_origin
700
701 subroutine ffhd_to_conserved_origin(ixI^L,ixO^L,w,x)
703 integer, intent(in) :: ixi^l, ixo^l
704 double precision, intent(inout) :: w(ixi^s, nw)
705 double precision, intent(in) :: x(ixi^s, 1:ndim)
706
707 if(ffhd_energy) then
708 w(ixo^s,e_)=w(ixo^s,p_)*inv_gamma_1+half*w(ixo^s,mom(1))**2*w(ixo^s,rho_)
709 end if
710 w(ixo^s,mom(1))=w(ixo^s,rho_)*w(ixo^s,mom(1))
711 end subroutine ffhd_to_conserved_origin
712
713 subroutine ffhd_to_primitive_origin(ixI^L,ixO^L,w,x)
715 integer, intent(in) :: ixi^l, ixo^l
716 double precision, intent(inout) :: w(ixi^s, nw)
717 double precision, intent(in) :: x(ixi^s, 1:ndim)
718
719 if(fix_small_values) then
720 call ffhd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'ffhd_to_primitive_origin')
721 end if
722
723 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/w(ixo^s,rho_)
724 if(ffhd_energy) then
725 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-half*w(ixo^s,rho_)*w(ixo^s,mom(1))**2)
726 end if
727 end subroutine ffhd_to_primitive_origin
728
729 subroutine ffhd_ei_to_e(ixI^L,ixO^L,w,x)
731 integer, intent(in) :: ixi^l, ixo^l
732 double precision, intent(inout) :: w(ixi^s, nw)
733 double precision, intent(in) :: x(ixi^s, 1:ndim)
734
735 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
736 end subroutine ffhd_ei_to_e
737
738 subroutine ffhd_e_to_ei(ixI^L,ixO^L,w,x)
740 integer, intent(in) :: ixi^l, ixo^l
741 double precision, intent(inout) :: w(ixi^s, nw)
742 double precision, intent(in) :: x(ixi^s, 1:ndim)
743
744 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
745 if(fix_small_values) then
746 call ffhd_handle_small_ei(w,x,ixi^l,ixi^l,e_,'ffhd_e_to_ei')
747 end if
748 end subroutine ffhd_e_to_ei
749
750 subroutine ffhd_handle_small_values_origin(primitive, w, x, ixI^L, ixO^L, subname)
753 logical, intent(in) :: primitive
754 integer, intent(in) :: ixi^l,ixo^l
755 double precision, intent(inout) :: w(ixi^s,1:nw)
756 double precision, intent(in) :: x(ixi^s,1:ndim)
757 character(len=*), intent(in) :: subname
758
759 logical :: flag(ixi^s,1:nw)
760 double precision :: tmp2(ixi^s)
761
762 call phys_check_w(primitive, ixi^l, ixi^l, w, flag)
763
764 if(any(flag)) then
765 select case (small_values_method)
766 case ("replace")
767 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
768 if(small_values_fix_iw(mom(1))) then
769 where(flag(ixo^s,rho_)) w(ixo^s, mom(1)) = 0.0d0
770 end if
771 if(ffhd_energy) then
772 if(primitive) then
773 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
774 else
775 where(flag(ixo^s,e_))
776 w(ixo^s,e_) = small_e+ffhd_kin_en(w,ixi^l,ixo^l)
777 end where
778 end if
779 end if
780 case ("average")
781 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
782 if(ffhd_energy) then
783 if(primitive) then
784 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
785 else
786 w(ixi^s,e_)=w(ixi^s,e_)-ffhd_kin_en(w,ixi^l,ixi^l)
787 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
788 w(ixi^s,e_)=w(ixi^s,e_)+ffhd_kin_en(w,ixi^l,ixi^l)
789 end if
790 end if
791 case default
792 if(.not.primitive) then
793 if(ffhd_energy) then
794 w(ixo^s,p_)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
795 end if
796 w(ixo^s,mom(1))=w(ixo^s,mom(1))/w(ixo^s,rho_)
797 end if
798 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
799 end select
800 end if
801 end subroutine ffhd_handle_small_values_origin
802
803 subroutine ffhd_get_v_origin(w,x,ixI^L,ixO^L,v)
805 integer, intent(in) :: ixi^l, ixo^l
806 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
807 double precision, intent(out) :: v(ixi^s,ndir)
808 double precision :: rho(ixi^s)
809 integer :: idir
810
811 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
812 rho(ixo^s)=1.d0/rho(ixo^s)
813 do idir=1,ndir
814 v(ixo^s,ndir) = w(ixo^s,mom(1))*block%B0(ixo^s,idir,0)*rho(ixo^s)
815 end do
816 end subroutine ffhd_get_v_origin
817
818 subroutine ffhd_get_v_idim(w,x,ixI^L,ixO^L,idim,v)
820 integer, intent(in) :: ixi^l, ixo^l, idim
821 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:ndim)
822 double precision, intent(out) :: v(ixi^s)
823 double precision :: rho(ixi^s)
824
825 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
826 v(ixo^s) = (w(ixo^s, mom(1))*block%B0(ixo^s,idim,0)) / rho(ixo^s)
827 end subroutine ffhd_get_v_idim
828
829 subroutine ffhd_get_cmax_origin(wprim,x,ixI^L,ixO^L,idim,cmax)
831 integer, intent(in) :: ixi^l, ixo^l, idim
832 ! w in primitive form
833 double precision, intent(in) :: wprim(ixi^s, nw), x(ixi^s,1:ndim)
834 double precision, intent(inout) :: cmax(ixi^s)
835
836 if(ffhd_energy) then
837 cmax(ixo^s)=dsqrt(ffhd_gamma*wprim(ixo^s,p_)/wprim(ixo^s,rho_))
838 else
839 cmax(ixo^s)=dsqrt(ffhd_gamma*ffhd_adiab*wprim(ixo^s,rho_)**gamma_1)
840 end if
841 cmax(ixo^s)=dabs(wprim(ixo^s,mom(1))*block%B0(ixo^s,idim,0))+cmax(ixo^s)
842
843 end subroutine ffhd_get_cmax_origin
844
845 subroutine ffhd_get_a2max(w,x,ixI^L,ixO^L,a2max)
847 use mod_geometry
848 integer, intent(in) :: ixi^l, ixo^l
849 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
850 double precision, intent(inout) :: a2max(ndim)
851 double precision :: a2(ixi^s,ndim,nw)
852 integer :: gxo^l,hxo^l,jxo^l,kxo^l,i,j
853
854 if(.not.slab_uniform)then
855 call mpistop("subroutine get_a2max in mod_ffhd_phys adopts cartesian setting")
856 endif
857 a2=zero
858 do i = 1,ndim
859 !> 4th order
860 hxo^l=ixo^l-kr(i,^d);
861 gxo^l=hxo^l-kr(i,^d);
862 jxo^l=ixo^l+kr(i,^d);
863 kxo^l=jxo^l+kr(i,^d);
864 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
865 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
866 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/dxlevel(i)**2
867 end do
868 end subroutine ffhd_get_a2max
869
870 subroutine ffhd_get_tcutoff(ixI^L,ixO^L,w,x,Tco_local,Tmax_local)
872 use mod_geometry
873 integer, intent(in) :: ixi^l,ixo^l
874 double precision, intent(in) :: x(ixi^s,1:ndim)
875 double precision, intent(inout) :: w(ixi^s,1:nw)
876 double precision, intent(out) :: tco_local,tmax_local
877 double precision, parameter :: trac_delta=0.25d0
878 double precision :: te(ixi^s),lts(ixi^s)
879 double precision, dimension(ixI^S,1:ndim) :: gradt
880 double precision :: bdir(ndim)
881 double precision :: ltrc,ltrp,altr
882 integer :: idims,ix^d,jxo^l,hxo^l,ixa^d,ixb^d
883 integer :: jxp^l,hxp^l,ixp^l,ixq^l
884 logical :: lrlt(ixi^s)
885
886 call ffhd_get_temperature(w,x,ixi^l,ixi^l,te)
887 tco_local=zero
888 tmax_local=maxval(te(ixo^s))
889
890 {^ifoned
891 select case(ffhd_trac_type)
892 case(0)
893 !> test case, fixed cutoff temperature
894 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
895 case(1)
896 hxo^l=ixo^l-1;
897 jxo^l=ixo^l+1;
898 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
899 lrlt=.false.
900 where(lts(ixo^s) > trac_delta)
901 lrlt(ixo^s)=.true.
902 end where
903 if(any(lrlt(ixo^s))) then
904 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
905 end if
906 case(2)
907 !> iijima et al. 2021, LTRAC method
908 ltrc=1.5d0
909 ltrp=4.d0
910 ixp^l=ixo^l^ladd1;
911 hxo^l=ixo^l-1;
912 jxo^l=ixo^l+1;
913 hxp^l=ixp^l-1;
914 jxp^l=ixp^l+1;
915 lts(ixp^s)=0.5d0*dabs(te(jxp^s)-te(hxp^s))/te(ixp^s)
916 lts(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
917 lts(ixo^s)=0.25d0*(lts(jxo^s)+two*lts(ixo^s)+lts(hxo^s))
918 block%wextra(ixo^s,tcoff_)=te(ixo^s)*lts(ixo^s)**0.4d0
919 case default
920 call mpistop("ffhd_trac_type not allowed for 1D simulation")
921 end select
922 }
923 {^nooned
924 select case(ffhd_trac_type)
925 case(0)
926 !> test case, fixed cutoff temperature
927 if(slab_uniform) then
928 !> assume cgs units
929 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)
930 else
931 block%wextra(ixi^s,tcoff_)=2.5d5/unit_temperature
932 end if
933 case(1,4,6)
934 do idims=1,ndim
935 call gradient(te,ixi^l,ixo^l,idims,gradt(ixi^s,idims))
936 end do
937 if(ffhd_trac_type .gt. 1) then
938 ! B direction at cell center
939 bdir=zero
940 {do ixa^d=0,1\}
941 ixb^d=(ixomin^d+ixomax^d-1)/2+ixa^d;
942 bdir(1:ndim)=bdir(1:ndim)+block%B0(ixb^d,1:ndim,0)
943 {end do\}
944 if(sum(bdir(:)**2) .gt. zero) then
945 bdir(1:ndim)=bdir(1:ndim)/dsqrt(sum(bdir(:)**2))
946 end if
947 block%special_values(3:ndim+2)=bdir(1:ndim)
948 end if
949 {do ix^db=ixomin^db,ixomax^db\}
950 {^iftwod
951 ! temperature length scale inversed
952 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2))*&
953 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
954 }
955 {^ifthreed
956 ! temperature length scale inversed
957 lts(ix^d)=min(block%ds(ix^d,1),block%ds(ix^d,2),block%ds(ix^d,3))*&
958 abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
959 }
960 if(lts(ix^d)>trac_delta) then
961 block%special_values(1)=max(block%special_values(1),te(ix^d))
962 end if
963 {end do\}
964 block%special_values(2)=tmax_local
965 case(2)
966 !> iijima et al. 2021, LTRAC method
967 ltrc=1.5d0
968 ltrp=4.d0
969 ixp^l=ixo^l^ladd2;
970 do idims=1,ndim
971 ixq^l=ixp^l;
972 hxp^l=ixp^l;
973 jxp^l=ixp^l;
974 select case(idims)
975 {case(^d)
976 ixqmin^d=ixqmin^d+1
977 ixqmax^d=ixqmax^d-1
978 hxpmax^d=ixpmin^d
979 jxpmin^d=ixpmax^d
980 \}
981 end select
982 call gradient(te,ixi^l,ixq^l,idims,gradt(ixi^s,idims))
983 call gradientf(te,x,ixi^l,hxp^l,idims,gradt(ixi^s,idims),nghostcells,.true.)
984 call gradientf(te,x,ixi^l,jxp^l,idims,gradt(ixi^s,idims),nghostcells,.false.)
985 end do
986 {do ix^db=ixpmin^db,ixpmax^db\}
987 ! temperature length scale inversed
988 lts(ix^d)=abs(^d&gradt({ix^d},^d)*block%B0({ix^d},^d,0)+)/te(ix^d)
989 ! fraction of cells size to temperature length scale
990 lts(ix^d)=min(^d&block%ds({ix^d},^d))*lts(ix^d)
991 lts(ix^d)=max(one,(exp(lts(ix^d))/ltrc)**ltrp)
992 {end do\}
993
994 ! need one ghost layer for thermal conductivity
995 ixp^l=ixo^l^ladd1;
996 {do ix^db=ixpmin^db,ixpmax^db\}
997 {^iftwod
998 altr=0.25d0*((lts(ix1-1,ix2)+two*lts(ix^d)+lts(ix1+1,ix2))*block%B0(ix^d,1,0)**2+&
999 (lts(ix1,ix2-1)+two*lts(ix^d)+lts(ix1,ix2+1))*block%B0(ix^d,2,0)**2)
1000 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
1001 }
1002 {^ifthreed
1003 altr=0.25d0*((lts(ix1-1,ix2,ix3)+two*lts(ix^d)+lts(ix1+1,ix2,ix3))*block%B0(ix^d,1,0)**2+&
1004 (lts(ix1,ix2-1,ix3)+two*lts(ix^d)+lts(ix1,ix2+1,ix3))*block%B0(ix^d,2,0)**2+&
1005 (lts(ix1,ix2,ix3-1)+two*lts(ix^d)+lts(ix1,ix2,ix3+1))*block%B0(ix^d,3,0)**2)
1006 block%wextra(ix^d,tcoff_)=te(ix^d)*altr**0.4d0
1007 }
1008 {end do\}
1009 case(3,5)
1010 !> do nothing here
1011 case default
1012 call mpistop("unknown ffhd_trac_type")
1013 end select
1014 }
1015 end subroutine ffhd_get_tcutoff
1016
1017 subroutine ffhd_get_cbounds(wLC,wRC,wLp,wRp,x,ixI^L,ixO^L,idim,Hspeed,cmax,cmin)
1019 integer, intent(in) :: ixi^l, ixo^l, idim
1020 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1021 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1022 double precision, intent(in) :: x(ixi^s,1:ndim)
1023 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1024 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1025 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
1026 double precision :: wmean(ixi^s,nw)
1027 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1028
1029 select case (boundspeed)
1030 case (1)
1031 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1032 ! Methods for Fluid Dynamics" by Toro.
1033 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1034 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1035 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1036 umean(ixo^s)=(wlp(ixo^s,mom(1))*tmp1(ixo^s)&
1037 +wrp(ixo^s,mom(1))*tmp2(ixo^s))*tmp3(ixo^s)
1038 umean(ixo^s)=umean(ixo^s)*block%B0(ixo^s,idim,idim)
1039 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1040 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1041 dmean(ixo^s)=(tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1042 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1043 ((wrp(ixo^s,mom(1))-wlp(ixo^s,mom(1)))*block%B0(ixo^s,idim,idim))**2
1044 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1045 if(present(cmin)) then
1046 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1047 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1048 else
1049 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1050 end if
1051 case (2)
1052 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1053 tmp1(ixo^s)=wmean(ixo^s,mom(1))*block%B0(ixo^s,idim,idim)/wmean(ixo^s,rho_)
1054 call ffhd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1055 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1056 if(present(cmin)) then
1057 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1058 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1059 else
1060 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1061 end if
1062 case (3)
1063 ! Miyoshi 2005 JCP 208, 315 equation (67)
1064 call ffhd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1065 call ffhd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1066 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1067 if(present(cmin)) then
1068 cmin(ixo^s,1)=min(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1069 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))-csoundl(ixo^s)
1070 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1071 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1072 else
1073 cmax(ixo^s,1)=max(wlp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim),&
1074 wrp(ixo^s,mom(1))*block%B0(ixo^s,idim,idim))+csoundl(ixo^s)
1075 end if
1076 end select
1077 end subroutine ffhd_get_cbounds
1078
1079 subroutine ffhd_get_pthermal_iso(w,x,ixI^L,ixO^L,pth)
1081
1082 integer, intent(in) :: ixi^l, ixo^l
1083 double precision, intent(in) :: w(ixi^s,nw)
1084 double precision, intent(in) :: x(ixi^s,1:ndim)
1085 double precision, intent(out):: pth(ixi^s)
1086
1087 call ffhd_get_rho(w,x,ixi^l,ixo^l,pth)
1088 pth(ixo^s)=ffhd_adiab*pth(ixo^s)**ffhd_gamma
1089 end subroutine ffhd_get_pthermal_iso
1090
1091 subroutine ffhd_get_pthermal_origin(w,x,ixI^L,ixO^L,pth)
1094 integer, intent(in) :: ixi^l, ixo^l
1095 double precision, intent(in) :: w(ixi^s,nw)
1096 double precision, intent(in) :: x(ixi^s,1:ndim)
1097 double precision, intent(out):: pth(ixi^s)
1098 integer :: iw, ix^d
1099
1100 pth(ixo^s)=gamma_1*(w(ixo^s,e_)-ffhd_kin_en(w,ixi^l,ixo^l))
1101 if (fix_small_values) then
1102 {do ix^db= ixo^lim^db\}
1103 if(pth(ix^d)<small_pressure) then
1104 pth(ix^d)=small_pressure
1105 end if
1106 {end do^D&\}
1107 elseif(check_small_values) then
1108 {do ix^db= ixo^lim^db\}
1109 if(pth(ix^d)<small_pressure) then
1110 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1111 " encountered when call ffhd_get_pthermal"
1112 write(*,*) "Iteration: ", it, " Time: ", global_time
1113 write(*,*) "Location: ", x(ix^d,:)
1114 write(*,*) "Cell number: ", ix^d
1115 do iw=1,nw
1116 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1117 end do
1118 ! use erroneous arithmetic operation to crash the run
1119 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1120 write(*,*) "Saving status at the previous time step"
1121 crash=.true.
1122 end if
1123 {end do^D&\}
1124 end if
1125 end subroutine ffhd_get_pthermal_origin
1126
1127 subroutine ffhd_get_temperature_from_te(w, x, ixI^L, ixO^L, res)
1129 integer, intent(in) :: ixi^l, ixo^l
1130 double precision, intent(in) :: w(ixi^s, 1:nw)
1131 double precision, intent(in) :: x(ixi^s, 1:ndim)
1132 double precision, intent(out):: res(ixi^s)
1133
1134 res(ixo^s) = w(ixo^s, te_)
1135 end subroutine ffhd_get_temperature_from_te
1136
1137 subroutine ffhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1139 integer, intent(in) :: ixi^l, ixo^l
1140 double precision, intent(in) :: w(ixi^s, 1:nw)
1141 double precision, intent(in) :: x(ixi^s, 1:ndim)
1142 double precision, intent(out):: res(ixi^s)
1143 double precision :: r(ixi^s)
1144
1145 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1146 res(ixo^s) = gamma_1 * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1147 end subroutine ffhd_get_temperature_from_eint
1148
1149 subroutine ffhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1151 integer, intent(in) :: ixi^l, ixo^l
1152 double precision, intent(in) :: w(ixi^s, 1:nw)
1153 double precision, intent(in) :: x(ixi^s, 1:ndim)
1154 double precision, intent(out):: res(ixi^s)
1155
1156 double precision :: r(ixi^s)
1157
1158 call ffhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1159 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,res)
1160 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1161 end subroutine ffhd_get_temperature_from_etot
1162
1163 subroutine ffhd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1165 integer, intent(in) :: ixi^l, ixo^l
1166 double precision, intent(in) :: w(ixi^s,nw)
1167 double precision, intent(in) :: x(ixi^s,1:ndim)
1168 double precision, intent(out) :: csound2(ixi^s)
1169 double precision :: rho(ixi^s)
1170
1171 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1172 if(ffhd_energy) then
1173 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1174 csound2(ixo^s)=ffhd_gamma*csound2(ixo^s)/rho(ixo^s)
1175 else
1176 csound2(ixo^s)=ffhd_gamma*ffhd_adiab*rho(ixo^s)**gamma_1
1177 end if
1178 end subroutine ffhd_get_csound2
1179
1180 subroutine ffhd_get_flux(wC,w,x,ixI^L,ixO^L,idim,f)
1182 use mod_geometry
1183 integer, intent(in) :: ixi^l, ixo^l, idim
1184 ! conservative w
1185 double precision, intent(in) :: wc(ixi^s,nw)
1186 ! primitive w
1187 double precision, intent(in) :: w(ixi^s,nw)
1188 double precision, intent(in) :: x(ixi^s,1:ndim)
1189 double precision,intent(out) :: f(ixi^s,nwflux)
1190 double precision :: ptotal(ixo^s)
1191
1192 f(ixo^s,rho_)=w(ixo^s,mom(1))*w(ixo^s,rho_)*block%B0(ixo^s,idim,idim)
1193
1194 if(ffhd_energy) then
1195 ptotal(ixo^s)=w(ixo^s,p_)
1196 else
1197 ptotal(ixo^s)=ffhd_adiab*w(ixo^s,rho_)**ffhd_gamma
1198 end if
1199
1200 ! Get flux of momentum
1201 f(ixo^s,mom(1))=(wc(ixo^s,mom(1))*w(ixo^s,mom(1))+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1202
1203 ! Get flux of energy
1204 if(ffhd_energy) then
1205 f(ixo^s,e_)=w(ixo^s,mom(1))*(wc(ixo^s,e_)+ptotal(ixo^s))*block%B0(ixo^s,idim,idim)
1207 f(ixo^s,e_)=f(ixo^s,e_)+w(ixo^s,q_)*block%B0(ixo^s,idim,idim)
1208 f(ixo^s,q_)=zero
1209 end if
1210 end if
1211 end subroutine ffhd_get_flux
1212
1213 subroutine ffhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1217 integer, intent(in) :: ixi^l, ixo^l
1218 double precision, intent(in) :: qdt,dtfactor
1219 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw), x(ixi^s,1:ndim)
1220 double precision, intent(inout) :: w(ixi^s,1:nw)
1221 logical, intent(in) :: qsourcesplit
1222 logical, intent(inout) :: active
1223
1224 if (.not. qsourcesplit) then
1225 active = .true.
1226 call add_punitb(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1228 !!call add_hypertc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1229 call add_hypertc_source_orig(qdt,ixi^l,ixo^l,wct,w,x,wctprim)
1230 end if
1231 end if
1232
1233 if(ffhd_radiative_cooling) then
1234 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1235 w,x,qsourcesplit,active, rc_fl)
1236 end if
1237
1238 if(ffhd_gravity) then
1239 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,&
1240 w,x,ffhd_energy,qsourcesplit,active)
1241 end if
1242
1243 ! update temperature from new pressure, density, and old ionization degree
1245 if(.not.qsourcesplit) then
1246 active = .true.
1247 call ffhd_update_temperature(ixi^l,ixo^l,wct,w,x)
1248 end if
1249 end if
1250 end subroutine ffhd_add_source
1251
1252 subroutine add_punitb(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1254 use mod_geometry
1255 integer, intent(in) :: ixi^l,ixo^l
1256 double precision, intent(in) :: qdt
1257 double precision, intent(in) :: wct(ixi^s,1:nw),x(ixi^s,1:ndim)
1258 double precision, intent(in) :: wctprim(ixi^s,1:nw)
1259 double precision, intent(inout) :: w(ixi^s,1:nw)
1260
1261 integer :: idims,hxo^l
1262 double precision :: divb(ixi^s)
1263 double precision :: rhovpar(ixi^s),gradrhov(ixi^s)
1264
1265 divb=zero
1266 if(slab_uniform) then
1267 do idims=1,ndim
1268 hxo^l=ixo^l-kr(idims,^d);
1269 divb(ixo^s)=divb(ixo^s)+(block%B0(ixo^s,idims,idims)-block%B0(hxo^s,idims,idims))/dxlevel(idims)
1270 end do
1271 else
1272 call divvector(block%B0(ixi^s,1:ndir,0),ixi^l,ixo^l,divb)
1273 end if
1274 w(ixo^s,mom(1))=w(ixo^s,mom(1))+qdt*wctprim(ixo^s,p_)*divb(ixo^s)
1275 end subroutine add_punitb
1276
1277 subroutine ffhd_get_rho(w,x,ixI^L,ixO^L,rho)
1279 integer, intent(in) :: ixi^l, ixo^l
1280 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
1281 double precision, intent(out) :: rho(ixi^s)
1282
1283 rho(ixo^s) = w(ixo^s,rho_)
1284 end subroutine ffhd_get_rho
1285
1286 subroutine ffhd_handle_small_ei(w, x, ixI^L, ixO^L, ie, subname)
1289 integer, intent(in) :: ixi^l,ixo^l, ie
1290 double precision, intent(inout) :: w(ixi^s,1:nw)
1291 double precision, intent(in) :: x(ixi^s,1:ndim)
1292 character(len=*), intent(in) :: subname
1293 integer :: idir
1294 logical :: flag(ixi^s,1:nw)
1295 double precision :: rho(ixi^s)
1296
1297 flag=.false.
1298 where(w(ixo^s,ie)<small_e) flag(ixo^s,ie)=.true.
1299 if(any(flag(ixo^s,ie))) then
1300 select case (small_values_method)
1301 case ("replace")
1302 where(flag(ixo^s,ie)) w(ixo^s,ie)=small_e
1303 case ("average")
1304 call small_values_average(ixi^l, ixo^l, w, x, flag, ie)
1305 case default
1306 w(ixo^s,e_)=w(ixo^s,e_)*gamma_1
1307 call ffhd_get_rho(w,x,ixi^l,ixo^l,rho)
1308 w(ixo^s,mom(1)) = w(ixo^s,mom(1))/rho(ixo^s)
1309 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1310 end select
1311 end if
1312 end subroutine ffhd_handle_small_ei
1313
1314 subroutine ffhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1317 integer, intent(in) :: ixi^l, ixo^l
1318 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1319 double precision, intent(inout) :: w(ixi^s,1:nw)
1320 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1321
1322 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1323 call ffhd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1324 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1325 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1326 end subroutine ffhd_update_temperature
1327
1328 subroutine ffhd_get_dt(w,ixI^L,ixO^L,dtnew,dx^D,x)
1330 use mod_usr_methods
1331 use mod_gravity, only: gravity_get_dt
1332 integer, intent(in) :: ixi^l, ixo^l
1333 double precision, intent(inout) :: dtnew
1334 double precision, intent(in) :: dx^d
1335 double precision, intent(in) :: w(ixi^s,1:nw)
1336 double precision, intent(in) :: x(ixi^s,1:ndim)
1337
1338 dtnew = bigdouble
1339
1340 if(ffhd_gravity) then
1341 call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1342 end if
1343 end subroutine ffhd_get_dt
1344
1345 subroutine ffhd_add_source_geom(qdt,dtfactor,ixI^L,ixO^L,wCT,wprim,w,x)
1347 integer, intent(in) :: ixi^l, ixo^l
1348 double precision, intent(in) :: qdt, dtfactor,x(ixi^s,1:ndim)
1349 double precision, intent(inout) :: wct(ixi^s,1:nw), wprim(ixi^s,1:nw),w(ixi^s,1:nw)
1350
1351 ! no geometric source terms needed for ffhd, no divergences of tensors
1352 end subroutine ffhd_add_source_geom
1353
1354 function ffhd_kin_en_origin(w, ixI^L, ixO^L, inv_rho) result(ke)
1355 use mod_global_parameters, only: nw, ndim,block
1356 integer, intent(in) :: ixi^l, ixo^l
1357 double precision, intent(in) :: w(ixi^s, nw)
1358 double precision :: ke(ixo^s)
1359 double precision, intent(in), optional :: inv_rho(ixo^s)
1360
1361 if(present(inv_rho)) then
1362 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2*inv_rho(ixo^s)
1363 else
1364 ke(ixo^s)=0.5d0*w(ixo^s,mom(1))**2/w(ixo^s,rho_)
1365 end if
1366 end function ffhd_kin_en_origin
1367
1368 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1371 integer, intent(in) :: ixi^l, ixo^l
1372 double precision, intent(in) :: w(ixi^s,1:nw)
1373 double precision, intent(in) :: x(ixi^s,1:ndim)
1374 double precision, intent(out):: rfactor(ixi^s)
1375 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1376
1377 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1378 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)
1379 end subroutine rfactor_from_temperature_ionization
1380
1381 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1383 integer, intent(in) :: ixi^l, ixo^l
1384 double precision, intent(in) :: w(ixi^s,1:nw)
1385 double precision, intent(in) :: x(ixi^s,1:ndim)
1386 double precision, intent(out):: rfactor(ixi^s)
1387
1388 rfactor(ixo^s)=rr
1389 end subroutine rfactor_from_constant_ionization
1390
1391 subroutine add_hypertc_source(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1393 use mod_geometry
1394 integer, intent(in) :: ixi^l,ixo^l
1395 double precision, intent(in) :: qdt
1396 double precision, dimension(ixI^S,1:ndim), intent(in) :: x
1397 double precision, dimension(ixI^S,1:nw), intent(in) :: wct,wctprim
1398 double precision, dimension(ixI^S,1:nw), intent(inout) :: w
1399
1400 double precision, dimension(ixI^S) :: te,r,bgradt,gradt
1401 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1402 integer :: ix^d,idims
1403
1404 call ffhd_get_rfactor(wct,x,ixi^l,ixi^l,r)
1405 te(ixi^s)=wctprim(ixi^s,p_)/(r(ixi^s)*wct(ixi^s,rho_))
1406 {^ifoned
1407 call gradient(te,ixi^l,ixo^l,1,bgradt,2)
1408 }
1409 {^nooned
1410 bgradt(ixo^s)=zero
1411 do idims=1,ndim
1412 ! compute gradient conform the geometry, 4th order CD for uniform cartesian by setting 2
1413 call gradient(te,ixi^l,ixo^l,idims,gradt,2)
1414 bgradt(ixo^s)=bgradt(ixo^s)+(block%B0(ixo^s,idims,0))*gradt(ixo^s)
1415 enddo
1416 }
1417 {do ix^db=ixomin^db,ixomax^db\}
1418 if(ffhd_trac) then
1419 r(ix^d)=max(te(ix^d),block%wextra(ix^d,tcoff_))
1420 else
1421 r(ix^d)=te(ix^d)
1422 endif
1423 sigma_t5=hypertc_kappa*dsqrt(r(ix^d)**5)
1424 sigma_t7=sigma_t5*r(ix^d)
1425 sigmat5_bgradt=sigma_t5*bgradt(ix^d)
1426 if(ffhd_htc_sat) then
1427 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))
1428 tau=max(4.d0*dt, f_sat*sigma_t7/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1429 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1430 else
1431 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1432 max(4.d0*dt, sigma_t7/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1433 end if
1434 {end do\}
1435
1436 end subroutine add_hypertc_source
1437
1438 subroutine add_hypertc_source_orig(qdt,ixI^L,ixO^L,wCT,w,x,wCTprim)
1440
1441 integer, intent(in) :: ixi^l,ixo^l
1442 double precision, intent(in) :: qdt
1443 double precision, dimension(ixI^S,1:ndim), intent(in) :: x
1444 double precision, dimension(ixI^S,1:nw), intent(in) :: wct,wctprim
1445 double precision, dimension(ixI^S,1:nw), intent(inout) :: w
1446
1447 double precision, dimension(ixI^S) :: te
1448 double precision :: sigma_t5,sigma_t7,sigmat5_bgradt,f_sat,tau
1449 integer :: ix^d
1450
1451 te(ixi^s)=wctprim(ixi^s,p_)/wct(ixi^s,rho_)
1452 ! temperature on face T_(i+1/2)=(7(T_i+T_(i+1))-(T_(i-1)+T_(i+2)))/12
1453 ! T_(i+1/2)-T_(i-1/2)=(8(T_(i+1)-T_(i-1))-T_(i+2)+T_(i-2))/12
1454 {^iftwod
1455 do ix2=ixomin2,ixomax2
1456 do ix1=ixomin1,ixomax1
1457 if(ffhd_trac) then
1458 if(te(ix^d)<block%wextra(ix^d,tcoff_)) then
1459 sigma_t5=hypertc_kappa*sqrt(block%wextra(ix^d,tcoff_)**5)
1460 sigma_t7=sigma_t5*block%wextra(ix^d,tcoff_)
1461 else
1462 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1463 sigma_t7=sigma_t5*te(ix^d)
1464 end if
1465 else
1466 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1467 sigma_t7=sigma_t5*te(ix^d)
1468 end if
1469 sigmat5_bgradt=sigma_t5*(&
1470 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)&
1471 +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))
1472 if(ffhd_htc_sat) then
1473 f_sat=one/(one+abs(sigmat5_bgradt))/(1.5d0*wct(ix^d,rho_)*(ffhd_gamma*wctprim(ix^d,p_)/wct(ix^d,rho_))**1.5d0)
1474 tau=max(4.d0*dt, f_sat*sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1475 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1476 else
1477 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1478 max(4.d0*dt, sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1479 end if
1480 end do
1481 end do
1482 }
1483 {^ifthreed
1484 do ix3=ixomin3,ixomax3
1485 do ix2=ixomin2,ixomax2
1486 do ix1=ixomin1,ixomax1
1487 if(ffhd_trac) then
1488 if(te(ix^d)<block%wextra(ix^d,tcoff_)) then
1489 sigma_t5=hypertc_kappa*sqrt(block%wextra(ix^d,tcoff_)**5)
1490 sigma_t7=sigma_t5*block%wextra(ix^d,tcoff_)
1491 else
1492 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1493 sigma_t7=sigma_t5*te(ix^d)
1494 end if
1495 else
1496 sigma_t5=hypertc_kappa*sqrt(te(ix^d)**5)
1497 sigma_t7=sigma_t5*te(ix^d)
1498 end if
1499 sigmat5_bgradt=sigma_t5*(&
1500 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)&
1501 +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)&
1502 +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))
1503 if(ffhd_htc_sat) then
1504 f_sat=one/(one+abs(sigmat5_bgradt))/(1.5d0*wct(ix^d,rho_)*(ffhd_gamma*wctprim(ix^d,p_)/wct(ix^d,rho_))**1.5d0)
1505 tau=max(4.d0*dt, f_sat*sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1506 w(ix^d,q_)=w(ix^d,q_)-qdt*(f_sat*sigmat5_bgradt+wct(ix^d,q_))/tau
1507 else
1508 w(ix^d,q_)=w(ix^d,q_)-qdt*(sigmat5_bgradt+wct(ix^d,q_))/&
1509 max(4.d0*dt, sigma_t7*courantpar**2/(wctprim(ix^d,p_)*inv_gamma_1*cmax_global**2))
1510 end if
1511 end do
1512 end do
1513 end do
1514 }
1515 end subroutine add_hypertc_source_orig
1516
1517end 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