MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
mod_hd_phys.t
Go to the documentation of this file.
1 !> Hydrodynamics physics module
2 module mod_hd_phys
6  use mod_physics
7  use mod_comm_lib, only: mpistop
8  implicit none
9  private
10 
11  !> Whether an energy equation is used
12  logical, public, protected :: hd_energy = .true.
13 
14  !> Whether thermal conduction is added
15  logical, public, protected :: hd_thermal_conduction = .false.
16  type(tc_fluid), allocatable, public :: tc_fl
17  type(te_fluid), allocatable, public :: te_fl_hd
18 
19  !> Whether radiative cooling is added
20  logical, public, protected :: hd_radiative_cooling = .false.
21  type(rc_fluid), allocatable, public :: rc_fl
22 
23  !> Whether dust is added
24  logical, public, protected :: hd_dust = .false.
25 
26  !> Whether viscosity is added
27  logical, public, protected :: hd_viscosity = .false.
28 
29  !> Whether gravity is added
30  logical, public, protected :: hd_gravity = .false.
31 
32  !> Whether particles module is added
33  logical, public, protected :: hd_particles = .false.
34 
35  !> Whether rotating frame is activated
36  logical, public, protected :: hd_rotating_frame = .false.
37 
38  !> Whether CAK radiation line force is activated
39  logical, public, protected :: hd_cak_force = .false.
40 
41  !> Number of tracer species
42  integer, public, protected :: hd_n_tracer = 0
43 
44  !> Whether plasma is partially ionized
45  logical, public, protected :: hd_partial_ionization = .false.
46 
47  !> Index of the density (in the w array)
48  integer, public, protected :: rho_
49 
50  !> Indices of the momentum density
51  integer, allocatable, public, protected :: mom(:)
52 
53  !> Indices of the tracers
54  integer, allocatable, public, protected :: tracer(:)
55 
56  !> Index of the energy density (-1 if not present)
57  integer, public, protected :: e_
58 
59  !> Index of the gas pressure (-1 if not present) should equal e_
60  integer, public, protected :: p_
61 
62  !> Indices of temperature
63  integer, public, protected :: te_
64 
65  !> Index of the cutoff temperature for the TRAC method
66  integer, public, protected :: tcoff_
67 
68  !> The adiabatic index
69  double precision, public :: hd_gamma = 5.d0/3.0d0
70 
71  !> The adiabatic constant
72  double precision, public :: hd_adiab = 1.0d0
73 
74  !> The small_est allowed energy
75  double precision, protected :: small_e
76 
77  !> Whether TRAC method is used
78  logical, public, protected :: hd_trac = .false.
79  integer, public, protected :: hd_trac_type = 1
80 
81  !> Allows overruling default corner filling (for debug mode, since otherwise corner primitives fail)
82  logical, public, protected :: hd_force_diagonal = .false.
83 
84  !> Helium abundance over Hydrogen
85  double precision, public, protected :: he_abundance=0.1d0
86  !> Ionization fraction of H
87  !> H_ion_fr = H+/(H+ + H)
88  double precision, public, protected :: h_ion_fr=1d0
89  !> Ionization fraction of He
90  !> He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
91  double precision, public, protected :: he_ion_fr=1d0
92  !> Ratio of number He2+ / number He+ + He2+
93  !> He_ion_fr2 = He2+/(He2+ + He+)
94  double precision, public, protected :: he_ion_fr2=1d0
95  ! used for eq of state when it is not defined by units,
96  ! the units do not contain terms related to ionization fraction
97  ! and it is p = RR * rho * T
98  double precision, public, protected :: rr=1d0
99  ! remove the below flag and assume default value = .false.
100  ! when eq state properly implemented everywhere
101  ! and not anymore through units
102  logical, public, protected :: eq_state_units = .true.
103 
104  procedure(sub_get_pthermal), pointer :: hd_get_rfactor => null()
105  ! Public methods
106  public :: hd_phys_init
107  public :: hd_kin_en
108  public :: hd_get_pthermal
109  public :: hd_get_csound2
110  public :: hd_to_conserved
111  public :: hd_to_primitive
112  public :: hd_check_params
113  public :: hd_check_w
114 
115 contains
116 
117  !> Read this module's parameters from a file
118  subroutine hd_read_params(files)
120  character(len=*), intent(in) :: files(:)
121  integer :: n
122 
123  namelist /hd_list/ hd_energy, hd_n_tracer, hd_gamma, hd_adiab, &
128 
129  do n = 1, size(files)
130  open(unitpar, file=trim(files(n)), status="old")
131  read(unitpar, hd_list, end=111)
132 111 close(unitpar)
133  end do
134 
135  end subroutine hd_read_params
136 
137  !> Write this module's parameters to a snapsoht
138  subroutine hd_write_info(fh)
140  integer, intent(in) :: fh
141  integer, parameter :: n_par = 1
142  double precision :: values(n_par)
143  character(len=name_len) :: names(n_par)
144  integer, dimension(MPI_STATUS_SIZE) :: st
145  integer :: er
146 
147  call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
148 
149  names(1) = "gamma"
150  values(1) = hd_gamma
151  call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
152  call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
153  end subroutine hd_write_info
154 
155  !> Initialize the module
156  subroutine hd_phys_init()
160  use mod_dust, only: dust_init
161  use mod_viscosity, only: viscosity_init
162  use mod_gravity, only: gravity_init
163  use mod_particles, only: particles_init
165  use mod_cak_force, only: cak_init
169  use mod_usr_methods, only: usr_rfactor
170 
171  integer :: itr, idir
172 
173  call hd_read_params(par_files)
174 
175  physics_type = "hd"
176  phys_energy = hd_energy
177  phys_total_energy = hd_energy
178  phys_internal_e = .false.
179  phys_gamma = hd_gamma
180  phys_partial_ionization=hd_partial_ionization
181 
183  if(phys_trac) then
184  if(ndim .eq. 1) then
185  if(hd_trac_type .gt. 2) then
186  hd_trac_type=1
187  if(mype==0) write(*,*) 'WARNING: set hd_trac_type=1'
188  end if
190  else
191  phys_trac=.false.
192  if(mype==0) write(*,*) 'WARNING: set hd_trac=F when ndim>=2'
193  end if
194  end if
195 
196  ! set default gamma for polytropic/isothermal process
197  if(.not.hd_energy) then
198  if(hd_thermal_conduction) then
199  hd_thermal_conduction=.false.
200  if(mype==0) write(*,*) 'WARNING: set hd_thermal_conduction=F when hd_energy=F'
201  end if
202  if(hd_radiative_cooling) then
203  hd_radiative_cooling=.false.
204  if(mype==0) write(*,*) 'WARNING: set hd_radiative_cooling=F when hd_energy=F'
205  end if
206  end if
207  if(.not.eq_state_units) then
208  if(hd_partial_ionization) then
209  hd_partial_ionization=.false.
210  if(mype==0) write(*,*) 'WARNING: set hd_partial_ionization=F when eq_state_units=F'
211  end if
212  end if
214 
215  allocate(start_indices(number_species),stop_indices(number_species))
216 
217  ! set the index of the first flux variable for species 1
218  start_indices(1)=1
219  ! Determine flux variables
220  rho_ = var_set_rho()
221 
222  allocate(mom(ndir))
223  mom(:) = var_set_momentum(ndir)
224 
225  ! Set index of energy variable
226  if (hd_energy) then
227  e_ = var_set_energy()
228  p_ = e_
229  else
230  e_ = -1
231  p_ = -1
232  end if
233 
234  phys_get_dt => hd_get_dt
235  phys_get_cmax => hd_get_cmax
236  phys_get_a2max => hd_get_a2max
237  phys_get_tcutoff => hd_get_tcutoff
238  phys_get_cbounds => hd_get_cbounds
239  phys_get_flux => hd_get_flux
240  phys_add_source_geom => hd_add_source_geom
241  phys_add_source => hd_add_source
242  phys_to_conserved => hd_to_conserved
243  phys_to_primitive => hd_to_primitive
244  phys_check_params => hd_check_params
245  phys_check_w => hd_check_w
246  phys_get_pthermal => hd_get_pthermal
247  phys_get_v => hd_get_v
248  phys_write_info => hd_write_info
249  phys_handle_small_values => hd_handle_small_values
250 
251  ! Whether diagonal ghost cells are required for the physics
252  phys_req_diagonal = .false.
253 
254  ! derive units from basic units
255  call hd_physical_units()
256 
257  if (hd_dust) then
258  call dust_init(rho_, mom(:), e_)
259  endif
260 
261  if (hd_force_diagonal) then
262  ! ensure corners are filled, otherwise divide by zero when getting primitives
263  ! --> only for debug purposes
264  phys_req_diagonal = .true.
265  endif
266 
267  allocate(tracer(hd_n_tracer))
268 
269  ! Set starting index of tracers
270  do itr = 1, hd_n_tracer
271  tracer(itr) = var_set_fluxvar("trc", "trp", itr, need_bc=.false.)
272  end do
273 
274  ! set number of variables which need update ghostcells
275  nwgc=nwflux
276 
277  ! set the index of the last flux variable for species 1
278  stop_indices(1)=nwflux
279 
280  ! set temperature as an auxiliary variable to get ionization degree
281  if(hd_partial_ionization) then
282  te_ = var_set_auxvar('Te','Te')
283  else
284  te_ = -1
285  end if
286 
287  if(hd_trac) then
288  tcoff_ = var_set_wextra()
289  iw_tcoff=tcoff_
290  else
291  tcoff_ = -1
292  end if
293 
294  ! choose Rfactor in ideal gas law
295  if(hd_partial_ionization) then
297  phys_update_temperature => hd_update_temperature
298  else if(associated(usr_rfactor)) then
299  hd_get_rfactor=>usr_rfactor
300  else
301  hd_get_rfactor=>rfactor_from_constant_ionization
302  end if
303 
304  ! initialize thermal conduction module
305  if (hd_thermal_conduction) then
306  if (.not. hd_energy) &
307  call mpistop("thermal conduction needs hd_energy=T")
308  phys_req_diagonal = .true.
309 
310  call sts_init()
312 
313  allocate(tc_fl)
318  tc_fl%get_temperature_from_conserved => hd_get_temperature_from_etot
319  tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
320  tc_fl%get_rho => hd_get_rho
321  tc_fl%e_ = e_
322  end if
323 
324  ! Initialize radiative cooling module
325  if (hd_radiative_cooling) then
326  if (.not. hd_energy) &
327  call mpistop("radiative cooling needs hd_energy=T")
329  allocate(rc_fl)
331  rc_fl%get_rho => hd_get_rho
332  rc_fl%get_pthermal => hd_get_pthermal
333  rc_fl%get_var_Rfactor => hd_get_rfactor
334  rc_fl%e_ = e_
335  rc_fl%Tcoff_ = tcoff_
336  end if
337  allocate(te_fl_hd)
338  te_fl_hd%get_rho=> hd_get_rho
339  te_fl_hd%get_pthermal=> hd_get_pthermal
340  te_fl_hd%get_var_Rfactor => hd_get_rfactor
341 {^ifthreed
342  phys_te_images => hd_te_images
343 }
344  ! Initialize viscosity module
345  if (hd_viscosity) call viscosity_init(phys_wider_stencil,phys_req_diagonal)
346 
347  ! Initialize gravity module
348  if (hd_gravity) call gravity_init()
349 
350  ! Initialize rotating_frame module
352 
353  ! Initialize CAK radiation force module
354  if (hd_cak_force) call cak_init(hd_gamma)
355 
356  ! Initialize particles module
357  if (hd_particles) then
358  call particles_init()
359  phys_req_diagonal = .true.
360  end if
361 
362  ! Check whether custom flux types have been defined
363  if (.not. allocated(flux_type)) then
364  allocate(flux_type(ndir, nw))
365  flux_type = flux_default
366  else if (any(shape(flux_type) /= [ndir, nw])) then
367  call mpistop("phys_check error: flux_type has wrong shape")
368  end if
369 
370  nvector = 1 ! No. vector vars
371  allocate(iw_vector(nvector))
372  iw_vector(1) = mom(1) - 1
373  ! initialize ionization degree table
375 
376  end subroutine hd_phys_init
377 
378 {^ifthreed
379  subroutine hd_te_images
382  select case(convert_type)
383  case('EIvtiCCmpi','EIvtuCCmpi')
385  case('ESvtiCCmpi','ESvtuCCmpi')
387  case('SIvtiCCmpi','SIvtuCCmpi')
389  case('WIvtiCCmpi','WIvtuCCmpi')
391  case default
392  call mpistop("Error in synthesize emission: Unknown convert_type")
393  end select
394  end subroutine hd_te_images
395 }
396 !!start th cond
397  ! wrappers for STS functions in thermal_conductivity module
398  ! which take as argument the tc_fluid (defined in the physics module)
399  subroutine hd_sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
401  use mod_fix_conserve
403  integer, intent(in) :: ixI^L, ixO^L, igrid, nflux
404  double precision, intent(in) :: x(ixI^S,1:ndim)
405  double precision, intent(inout) :: wres(ixI^S,1:nw), w(ixI^S,1:nw)
406  double precision, intent(in) :: my_dt
407  logical, intent(in) :: fix_conserve_at_step
408  call sts_set_source_tc_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
409  end subroutine hd_sts_set_source_tc_hd
410 
411  function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
412  !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
413  !where tc_k_para_i=tc_k_para*B_i**2/B**2
414  !and T=p/rho
417 
418  integer, intent(in) :: ixi^l, ixo^l
419  double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
420  double precision, intent(in) :: w(ixi^s,1:nw)
421  double precision :: dtnew
422 
423  dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
424  end function hd_get_tc_dt_hd
425 
426  subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
427  ! move this in a different routine as in mhd if needed in more places
429  use mod_small_values
430 
431  integer, intent(in) :: ixI^L,ixO^L
432  double precision, intent(inout) :: w(ixI^S,1:nw)
433  double precision, intent(in) :: x(ixI^S,1:ndim)
434  integer, intent(in) :: step
435 
436  integer :: idir
437  logical :: flag(ixI^S,1:nw)
438  character(len=140) :: error_msg
439 
440  flag=.false.
441  where(w(ixo^s,e_)<small_e) flag(ixo^s,e_)=.true.
442  if(any(flag(ixo^s,e_))) then
443  select case (small_values_method)
444  case ("replace")
445  where(flag(ixo^s,e_)) w(ixo^s,e_)=small_e
446  case ("average")
447  call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
448  case default
449  ! small values error shows primitive variables
450  w(ixo^s,e_)=w(ixo^s,e_)*(hd_gamma - 1.0d0)
451  do idir = 1, ndir
452  w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,rho_)
453  end do
454  write(error_msg,"(a,i3)") "Thermal conduction step ", step
455  call small_values_error(w, x, ixi^l, ixo^l, flag, error_msg)
456  end select
457  end if
458  end subroutine hd_tc_handle_small_e
459 
460  ! fill in tc_fluid fields from namelist
461  subroutine tc_params_read_hd(fl)
463  type(tc_fluid), intent(inout) :: fl
464  integer :: n
465  logical :: tc_saturate=.false.
466  double precision :: tc_k_para=0d0
467 
468  namelist /tc_list/ tc_saturate, tc_k_para
469 
470  do n = 1, size(par_files)
471  open(unitpar, file=trim(par_files(n)), status="old")
472  read(unitpar, tc_list, end=111)
473 111 close(unitpar)
474  end do
475  fl%tc_saturate = tc_saturate
476  fl%tc_k_para = tc_k_para
477 
478  end subroutine tc_params_read_hd
479 
480  subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
482  integer, intent(in) :: ixI^L, ixO^L
483  double precision, intent(in) :: w(ixI^S,1:nw),x(ixI^S,1:ndim)
484  double precision, intent(out) :: rho(ixI^S)
485 
486  rho(ixo^s) = w(ixo^s,rho_)
487 
488  end subroutine hd_get_rho
489 
490 !!end th cond
491 !!rad cool
492  subroutine rc_params_read(fl)
494  use mod_constants, only: bigdouble
495  use mod_basic_types, only: std_len
496  type(rc_fluid), intent(inout) :: fl
497  integer :: n
498  ! list parameters
499  integer :: ncool = 4000
500  double precision :: cfrac=0.1d0
501 
502  !> Name of cooling curve
503  character(len=std_len) :: coolcurve='JCcorona'
504 
505  !> Name of cooling method
506  character(len=std_len) :: coolmethod='exact'
507 
508  !> Fixed temperature not lower than tlow
509  logical :: Tfix=.false.
510 
511  !> Lower limit of temperature
512  double precision :: tlow=bigdouble
513 
514  !> Add cooling source in a split way (.true.) or un-split way (.false.)
515  logical :: rc_split=.false.
516 
517 
518  namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
519 
520  do n = 1, size(par_files)
521  open(unitpar, file=trim(par_files(n)), status="old")
522  read(unitpar, rc_list, end=111)
523 111 close(unitpar)
524  end do
525 
526  fl%ncool=ncool
527  fl%coolcurve=coolcurve
528  fl%coolmethod=coolmethod
529  fl%tlow=tlow
530  fl%Tfix=tfix
531  fl%rc_split=rc_split
532  fl%cfrac=cfrac
533  end subroutine rc_params_read
534 !! end rad cool
535 
536  subroutine hd_check_params
539 
540  if (.not. hd_energy) then
541  if (hd_gamma <= 0.0d0) call mpistop ("Error: hd_gamma <= 0")
542  if (hd_adiab < 0.0d0) call mpistop ("Error: hd_adiab < 0")
544  else
545  if (hd_gamma <= 0.0d0 .or. hd_gamma == 1.0d0) &
546  call mpistop ("Error: hd_gamma <= 0 or hd_gamma == 1.0")
547  small_e = small_pressure/(hd_gamma - 1.0d0)
548  end if
549 
550  if (hd_dust) call dust_check_params()
551  if(use_imex_scheme) then
552  ! implicit dust update
553  phys_implicit_update => dust_implicit_update
554  phys_evaluate_implicit => dust_evaluate_implicit
555  endif
556 
557  end subroutine hd_check_params
558 
559  subroutine hd_physical_units
561  double precision :: mp,kB
562  double precision :: a,b
563  ! Derive scaling units
564  if(si_unit) then
565  mp=mp_si
566  kb=kb_si
567  else
568  mp=mp_cgs
569  kb=kb_cgs
570  end if
571  if(eq_state_units) then
572  a = 1d0 + 4d0 * he_abundance
573  if(hd_partial_ionization) then
574  b = 2.d0+3.d0*he_abundance
575  else
576  b = 1d0 + h_ion_fr + he_abundance*(he_ion_fr*(he_ion_fr2 + 1d0)+1d0)
577  end if
578  rr = 1d0
579  else
580  a = 1d0
581  b = 1d0
582  rr = (1d0 + h_ion_fr + he_abundance*(he_ion_fr*(he_ion_fr2 + 1d0)+1d0))/(1d0 + 4d0 * he_abundance)
583  end if
584  if(unit_density/=1.d0) then
586  else
587  ! unit of numberdensity is independent by default
589  end if
590  if(unit_velocity/=1.d0) then
593  else if(unit_pressure/=1.d0) then
596  else
597  ! unit of temperature is independent by default
600  end if
601  if(unit_time/=1.d0) then
603  else
604  ! unit of length is independent by default
606  end if
608 
609  end subroutine hd_physical_units
610 
611  !> Returns logical argument flag where values are ok
612  subroutine hd_check_w(primitive, ixI^L, ixO^L, w, flag)
614  use mod_dust, only: dust_check_w
615 
616  logical, intent(in) :: primitive
617  integer, intent(in) :: ixi^l, ixo^l
618  double precision, intent(in) :: w(ixi^s, nw)
619  logical, intent(inout) :: flag(ixi^s,1:nw)
620  double precision :: tmp(ixi^s)
621 
622  flag=.false.
623 
624  if (hd_energy) then
625  if (primitive) then
626  where(w(ixo^s, e_) < small_pressure) flag(ixo^s,e_) = .true.
627  else
628  tmp(ixo^s) = (hd_gamma - 1.0d0)*(w(ixo^s, e_) - &
629  hd_kin_en(w, ixi^l, ixo^l))
630  where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_) = .true.
631  endif
632  end if
633 
634  where(w(ixo^s, rho_) < small_density) flag(ixo^s,rho_) = .true.
635 
636  if(hd_dust) call dust_check_w(ixi^l,ixo^l,w,flag)
637 
638  end subroutine hd_check_w
639 
640  !> Transform primitive variables into conservative ones
641  subroutine hd_to_conserved(ixI^L, ixO^L, w, x)
643  use mod_dust, only: dust_to_conserved
644  integer, intent(in) :: ixi^l, ixo^l
645  double precision, intent(inout) :: w(ixi^s, nw)
646  double precision, intent(in) :: x(ixi^s, 1:ndim)
647  double precision :: invgam
648  integer :: idir, itr
649 
650  !!if (fix_small_values) then
651  !! call hd_handle_small_values(.true., w, x, ixI^L, ixO^L, 'hd_to_conserved')
652  !!end if
653 
654  if (hd_energy) then
655  invgam = 1.d0/(hd_gamma - 1.0d0)
656  ! Calculate total energy from pressure and kinetic energy
657  w(ixo^s, e_) = w(ixo^s, e_) * invgam + &
658  0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * w(ixo^s, rho_)
659  end if
660 
661  ! Convert velocity to momentum
662  do idir = 1, ndir
663  w(ixo^s, mom(idir)) = w(ixo^s, rho_) * w(ixo^s, mom(idir))
664  end do
665 
666  if (hd_dust) then
667  call dust_to_conserved(ixi^l, ixo^l, w, x)
668  end if
669 
670  end subroutine hd_to_conserved
671 
672  !> Transform conservative variables into primitive ones
673  subroutine hd_to_primitive(ixI^L, ixO^L, w, x)
675  use mod_dust, only: dust_to_primitive
676  integer, intent(in) :: ixi^l, ixo^l
677  double precision, intent(inout) :: w(ixi^s, nw)
678  double precision, intent(in) :: x(ixi^s, 1:ndim)
679  integer :: itr, idir
680  double precision :: inv_rho(ixo^s)
681 
682  if (fix_small_values) then
683  call hd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'hd_to_primitive')
684  end if
685 
686  inv_rho = 1.0d0 / w(ixo^s, rho_)
687 
688  if (hd_energy) then
689  ! Compute pressure
690  w(ixo^s, e_) = (hd_gamma - 1.0d0) * (w(ixo^s, e_) - &
691  hd_kin_en(w, ixi^l, ixo^l, inv_rho))
692  end if
693 
694  ! Convert momentum to velocity
695  do idir = 1, ndir
696  w(ixo^s, mom(idir)) = w(ixo^s, mom(idir)) * inv_rho
697  end do
698 
699  ! Convert dust momentum to dust velocity
700  if (hd_dust) then
701  call dust_to_primitive(ixi^l, ixo^l, w, x)
702  end if
703 
704  end subroutine hd_to_primitive
705 
706  !> Transform internal energy to total energy
707  subroutine hd_ei_to_e(ixI^L,ixO^L,w,x)
709  integer, intent(in) :: ixI^L, ixO^L
710  double precision, intent(inout) :: w(ixI^S, nw)
711  double precision, intent(in) :: x(ixI^S, 1:ndim)
712 
713  ! Calculate total energy from internal and kinetic energy
714  w(ixo^s,e_)=w(ixo^s,e_)&
715  +hd_kin_en(w,ixi^l,ixo^l)
716 
717  end subroutine hd_ei_to_e
718 
719  !> Transform total energy to internal energy
720  subroutine hd_e_to_ei(ixI^L,ixO^L,w,x)
722  integer, intent(in) :: ixI^L, ixO^L
723  double precision, intent(inout) :: w(ixI^S, nw)
724  double precision, intent(in) :: x(ixI^S, 1:ndim)
725 
726  ! Calculate ei = e - ek
727  w(ixo^s,e_)=w(ixo^s,e_)&
728  -hd_kin_en(w,ixi^l,ixo^l)
729 
730  end subroutine hd_e_to_ei
731 
732  subroutine e_to_rhos(ixI^L, ixO^L, w, x)
734 
735  integer, intent(in) :: ixI^L, ixO^L
736  double precision :: w(ixI^S, nw)
737  double precision, intent(in) :: x(ixI^S, 1:ndim)
738 
739  if (hd_energy) then
740  w(ixo^s, e_) = (hd_gamma - 1.0d0) * w(ixo^s, rho_)**(1.0d0 - hd_gamma) * &
741  (w(ixo^s, e_) - hd_kin_en(w, ixi^l, ixo^l))
742  else
743  call mpistop("energy from entropy can not be used with -eos = iso !")
744  end if
745  end subroutine e_to_rhos
746 
747  subroutine rhos_to_e(ixI^L, ixO^L, w, x)
749 
750  integer, intent(in) :: ixI^L, ixO^L
751  double precision :: w(ixI^S, nw)
752  double precision, intent(in) :: x(ixI^S, 1:ndim)
753 
754  if (hd_energy) then
755  w(ixo^s, e_) = w(ixo^s, rho_)**(hd_gamma - 1.0d0) * w(ixo^s, e_) &
756  / (hd_gamma - 1.0d0) + hd_kin_en(w, ixi^l, ixo^l)
757  else
758  call mpistop("entropy from energy can not be used with -eos = iso !")
759  end if
760  end subroutine rhos_to_e
761 
762  !> Calculate v_i = m_i / rho within ixO^L
763  subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
765  integer, intent(in) :: ixI^L, ixO^L, idim
766  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
767  double precision, intent(out) :: v(ixI^S)
768 
769  v(ixo^s) = w(ixo^s, mom(idim)) / w(ixo^s, rho_)
770  end subroutine hd_get_v_idim
771 
772  !> Calculate velocity vector v_i = m_i / rho within ixO^L
773  subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
775 
776  integer, intent(in) :: ixI^L, ixO^L
777  double precision, intent(in) :: w(ixI^S,nw), x(ixI^S,1:^ND)
778  double precision, intent(out) :: v(ixI^S,1:ndir)
779 
780  integer :: idir
781 
782  do idir=1,ndir
783  v(ixo^s,idir) = w(ixo^s, mom(idir)) / w(ixo^s, rho_)
784  end do
785 
786  end subroutine hd_get_v
787 
788  !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
789  subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
791  use mod_dust, only: dust_get_cmax
792 
793  integer, intent(in) :: ixI^L, ixO^L, idim
794  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S, 1:ndim)
795  double precision, intent(inout) :: cmax(ixI^S)
796  double precision :: csound(ixI^S)
797  double precision :: v(ixI^S)
798 
799  call hd_get_v_idim(w, x, ixi^l, ixo^l, idim, v)
800  call hd_get_csound2(w,x,ixi^l,ixo^l,csound)
801  csound(ixo^s) = dsqrt(csound(ixo^s))
802 
803  cmax(ixo^s) = dabs(v(ixo^s))+csound(ixo^s)
804 
805  if (hd_dust) then
806  call dust_get_cmax(w, x, ixi^l, ixo^l, idim, cmax)
807  end if
808  end subroutine hd_get_cmax
809 
810  subroutine hd_get_a2max(w,x,ixI^L,ixO^L,a2max)
812 
813  integer, intent(in) :: ixI^L, ixO^L
814  double precision, intent(in) :: w(ixI^S, nw), x(ixI^S,1:ndim)
815  double precision, intent(inout) :: a2max(ndim)
816  double precision :: a2(ixI^S,ndim,nw)
817  integer :: gxO^L,hxO^L,jxO^L,kxO^L,i,j
818 
819  a2=zero
820  do i = 1,ndim
821  !> 4th order
822  hxo^l=ixo^l-kr(i,^d);
823  gxo^l=hxo^l-kr(i,^d);
824  jxo^l=ixo^l+kr(i,^d);
825  kxo^l=jxo^l+kr(i,^d);
826  a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
827  -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
828  a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/dxlevel(i)**2
829  end do
830  end subroutine hd_get_a2max
831 
832  !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
833  subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
835  integer, intent(in) :: ixI^L,ixO^L
836  double precision, intent(in) :: x(ixI^S,1:ndim)
837  double precision, intent(inout) :: w(ixI^S,1:nw)
838  double precision, intent(out) :: tco_local, Tmax_local
839 
840  double precision, parameter :: trac_delta=0.25d0
841  double precision :: tmp1(ixI^S),Te(ixI^S),lts(ixI^S)
842  double precision :: ltr(ixI^S),ltrc,ltrp,tcoff(ixI^S)
843  integer :: jxO^L,hxO^L
844  integer :: jxP^L,hxP^L,ixP^L
845  logical :: lrlt(ixI^S)
846 
847  {^ifoned
848  call hd_get_temperature_from_etot(w,x,ixi^l,ixi^l,te)
849 
850  tco_local=zero
851  tmax_local=maxval(te(ixo^s))
852  select case(hd_trac_type)
853  case(0)
854  w(ixi^s,tcoff_)=3.d5/unit_temperature
855  case(1)
856  hxo^l=ixo^l-1;
857  jxo^l=ixo^l+1;
858  lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
859  lrlt=.false.
860  where(lts(ixo^s) > trac_delta)
861  lrlt(ixo^s)=.true.
862  end where
863  if(any(lrlt(ixo^s))) then
864  tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
865  end if
866  case(2)
867  !> iijima et al. 2021, LTRAC method
868  ltrc=1.5d0
869  ltrp=2.5d0
870  ixp^l=ixo^l^ladd1;
871  hxo^l=ixo^l-1;
872  jxo^l=ixo^l+1;
873  hxp^l=ixp^l-1;
874  jxp^l=ixp^l+1;
875  lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
876  ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
877  w(ixo^s,tcoff_)=te(ixo^s)*&
878  (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
879  case default
880  call mpistop("mhd_trac_type not allowed for 1D simulation")
881  end select
882  }
883  end subroutine hd_get_tcutoff
884 
885  !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
886  subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
888  use mod_dust, only: dust_get_cmax
889  use mod_variables
890 
891  integer, intent(in) :: ixI^L, ixO^L, idim
892  ! conservative left and right status
893  double precision, intent(in) :: wLC(ixI^S, nw), wRC(ixI^S, nw)
894  ! primitive left and right status
895  double precision, intent(in) :: wLp(ixI^S, nw), wRp(ixI^S, nw)
896  double precision, intent(in) :: x(ixI^S, 1:ndim)
897  double precision, intent(inout) :: cmax(ixI^S,1:number_species)
898  double precision, intent(inout), optional :: cmin(ixI^S,1:number_species)
899  double precision, intent(in) :: Hspeed(ixI^S,1:number_species)
900 
901  double precision :: wmean(ixI^S,nw)
902  double precision, dimension(ixI^S) :: umean, dmean, csoundL, csoundR, tmp1,tmp2,tmp3
903  integer :: ix^D
904 
905  select case(boundspeed)
906  case (1)
907  ! This implements formula (10.52) from "Riemann Solvers and Numerical
908  ! Methods for Fluid Dynamics" by Toro.
909 
910  tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
911  tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
912  tmp3(ixo^s)=1.d0/(dsqrt(wlp(ixo^s,rho_))+dsqrt(wrp(ixo^s,rho_)))
913  umean(ixo^s)=(wlp(ixo^s,mom(idim))*tmp1(ixo^s)+wrp(ixo^s,mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
914 
915  if(hd_energy) then
916  csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
917  csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
918  else
919  call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
920  call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
921  end if
922 
923  dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
924  tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
925  (wrp(ixo^s,mom(idim))-wlp(ixo^s,mom(idim)))**2
926 
927  dmean(ixo^s)=dsqrt(dmean(ixo^s))
928  if(present(cmin)) then
929  cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
930  cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
931  if(h_correction) then
932  {do ix^db=ixomin^db,ixomax^db\}
933  cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
934  cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
935  {end do\}
936  end if
937  else
938  cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
939  end if
940 
941  if (hd_dust) then
942  wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
943  call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
944  end if
945 
946  case (2)
947  wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
948  tmp1(ixo^s)=wmean(ixo^s,mom(idim))/wmean(ixo^s,rho_)
949  call hd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
950  csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
951 
952  if(present(cmin)) then
953  cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
954  cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
955  if(h_correction) then
956  {do ix^db=ixomin^db,ixomax^db\}
957  cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
958  cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
959  {end do\}
960  end if
961  else
962  cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
963  end if
964 
965  if (hd_dust) then
966  call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
967  end if
968  case (3)
969  ! Miyoshi 2005 JCP 208, 315 equation (67)
970  if(hd_energy) then
971  csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
972  csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
973  else
974  call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
975  call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
976  end if
977  csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
978  if(present(cmin)) then
979  cmin(ixo^s,1)=min(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))-csoundl(ixo^s)
980  cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
981  if(h_correction) then
982  {do ix^db=ixomin^db,ixomax^db\}
983  cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
984  cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
985  {end do\}
986  end if
987  else
988  cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
989  end if
990  if (hd_dust) then
991  wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
992  call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
993  end if
994  end select
995 
996  end subroutine hd_get_cbounds
997 
998  !> Calculate the square of the thermal sound speed csound2 within ixO^L.
999  !> csound2=gamma*p/rho
1000  subroutine hd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1002  integer, intent(in) :: ixi^l, ixo^l
1003  double precision, intent(in) :: w(ixi^s,nw)
1004  double precision, intent(in) :: x(ixi^s,1:ndim)
1005  double precision, intent(out) :: csound2(ixi^s)
1006 
1007  call hd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1008  csound2(ixo^s)=hd_gamma*csound2(ixo^s)/w(ixo^s,rho_)
1009 
1010  end subroutine hd_get_csound2
1011 
1012  !> Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L
1013  subroutine hd_get_pthermal(w, x, ixI^L, ixO^L, pth)
1017 
1018  integer, intent(in) :: ixi^l, ixo^l
1019  double precision, intent(in) :: w(ixi^s, 1:nw)
1020  double precision, intent(in) :: x(ixi^s, 1:ndim)
1021  double precision, intent(out):: pth(ixi^s)
1022  integer :: iw, ix^d
1023 
1024  if (hd_energy) then
1025  pth(ixo^s) = (hd_gamma - 1.0d0) * (w(ixo^s, e_) - &
1026  hd_kin_en(w, ixi^l, ixo^l))
1027  else
1028  if (.not. associated(usr_set_pthermal)) then
1029  pth(ixo^s) = hd_adiab * w(ixo^s, rho_)**hd_gamma
1030  else
1031  call usr_set_pthermal(w,x,ixi^l,ixo^l,pth)
1032  end if
1033  end if
1034 
1035  if (fix_small_values) then
1036  {do ix^db= ixo^lim^db\}
1037  if(pth(ix^d)<small_pressure) then
1038  pth(ix^d)=small_pressure
1039  endif
1040  {enddo^d&\}
1041  else if (check_small_values) then
1042  {do ix^db= ixo^lim^db\}
1043  if(pth(ix^d)<small_pressure) then
1044  write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1045  " encountered when call hd_get_pthermal"
1046  write(*,*) "Iteration: ", it, " Time: ", global_time
1047  write(*,*) "Location: ", x(ix^d,:)
1048  write(*,*) "Cell number: ", ix^d
1049  do iw=1,nw
1050  write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1051  end do
1052  ! use erroneous arithmetic operation to crash the run
1053  if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1054  write(*,*) "Saving status at the previous time step"
1055  crash=.true.
1056  end if
1057  {enddo^d&\}
1058  end if
1059 
1060  end subroutine hd_get_pthermal
1061 
1062  !> Calculate temperature=p/rho when in e_ the total energy is stored
1063  subroutine hd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1065  integer, intent(in) :: ixI^L, ixO^L
1066  double precision, intent(in) :: w(ixI^S, 1:nw)
1067  double precision, intent(in) :: x(ixI^S, 1:ndim)
1068  double precision, intent(out):: res(ixI^S)
1069 
1070  double precision :: R(ixI^S)
1071 
1072  call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1073  call hd_get_pthermal(w, x, ixi^l, ixo^l, res)
1074  res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1075  end subroutine hd_get_temperature_from_etot
1076 
1077  !> Calculate temperature=p/rho when in e_ the internal energy is stored
1078  subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1080  integer, intent(in) :: ixI^L, ixO^L
1081  double precision, intent(in) :: w(ixI^S, 1:nw)
1082  double precision, intent(in) :: x(ixI^S, 1:ndim)
1083  double precision, intent(out):: res(ixI^S)
1084 
1085  double precision :: R(ixI^S)
1086 
1087  call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1088  res(ixo^s) = (hd_gamma - 1.0d0) * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1089  end subroutine hd_get_temperature_from_eint
1090 
1091  ! Calculate flux f_idim[iw]
1092  subroutine hd_get_flux_cons(w, x, ixI^L, ixO^L, idim, f)
1094  use mod_dust, only: dust_get_flux
1095 
1096  integer, intent(in) :: ixI^L, ixO^L, idim
1097  double precision, intent(in) :: w(ixI^S, 1:nw), x(ixI^S, 1:ndim)
1098  double precision, intent(out) :: f(ixI^S, nwflux)
1099  double precision :: pth(ixI^S), v(ixI^S),frame_vel(ixI^S)
1100  integer :: idir, itr
1101 
1102  call hd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1103  call hd_get_v_idim(w, x, ixi^l, ixo^l, idim, v)
1104 
1105  f(ixo^s, rho_) = v(ixo^s) * w(ixo^s, rho_)
1106 
1107  ! Momentum flux is v_i*m_i, +p in direction idim
1108  do idir = 1, ndir
1109  f(ixo^s, mom(idir)) = v(ixo^s) * w(ixo^s, mom(idir))
1110  end do
1111 
1112  f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1113 
1114  if(hd_energy) then
1115  ! Energy flux is v_i*(e + p)
1116  f(ixo^s, e_) = v(ixo^s) * (w(ixo^s, e_) + pth(ixo^s))
1117  end if
1118 
1119  do itr = 1, hd_n_tracer
1120  f(ixo^s, tracer(itr)) = v(ixo^s) * w(ixo^s, tracer(itr))
1121  end do
1122 
1123  ! Dust fluxes
1124  if (hd_dust) then
1125  call dust_get_flux(w, x, ixi^l, ixo^l, idim, f)
1126  end if
1127 
1128  end subroutine hd_get_flux_cons
1129 
1130  ! Calculate flux f_idim[iw]
1131  subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1133  use mod_dust, only: dust_get_flux_prim
1134  use mod_viscosity, only: visc_get_flux_prim ! viscInDiv
1135 
1136  integer, intent(in) :: ixI^L, ixO^L, idim
1137  ! conservative w
1138  double precision, intent(in) :: wC(ixI^S, 1:nw)
1139  ! primitive w
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) :: f(ixI^S, nwflux)
1143  double precision :: pth(ixI^S),frame_vel(ixI^S)
1144  integer :: idir, itr
1145 
1146  if (hd_energy) then
1147  pth(ixo^s) = w(ixo^s,p_)
1148  else
1149  call hd_get_pthermal(wc, x, ixi^l, ixo^l, pth)
1150  end if
1151 
1152  f(ixo^s, rho_) = w(ixo^s,mom(idim)) * w(ixo^s, rho_)
1153 
1154  ! Momentum flux is v_i*m_i, +p in direction idim
1155  do idir = 1, ndir
1156  f(ixo^s, mom(idir)) = w(ixo^s,mom(idim)) * wc(ixo^s, mom(idir))
1157  end do
1158 
1159  f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1160 
1161  if(hd_energy) then
1162  ! Energy flux is v_i*(e + p)
1163  f(ixo^s, e_) = w(ixo^s,mom(idim)) * (wc(ixo^s, e_) + w(ixo^s,p_))
1164  end if
1165 
1166  do itr = 1, hd_n_tracer
1167  f(ixo^s, tracer(itr)) = w(ixo^s,mom(idim)) * w(ixo^s, tracer(itr))
1168  end do
1169 
1170  ! Dust fluxes
1171  if (hd_dust) then
1172  call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1173  end if
1174 
1175  ! Viscosity fluxes - viscInDiv
1176  if (hd_viscosity) then
1177  call visc_get_flux_prim(w, x, ixi^l, ixo^l, idim, f, hd_energy)
1178  endif
1179 
1180  end subroutine hd_get_flux
1181 
1182  !> Add geometrical source terms to w
1183  !>
1184  !> Notice that the expressions of the geometrical terms depend only on ndir,
1185  !> not ndim. Eg, they are the same in 2.5D and in 3D, for any geometry.
1186  !>
1187  !> Ileyk : to do :
1188  !> - address the source term for the dust in case (coordinate == spherical)
1189  subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, w, x)
1191  use mod_usr_methods, only: usr_set_surface
1192  use mod_viscosity, only: visc_add_source_geom ! viscInDiv
1195  use mod_geometry
1196  integer, intent(in) :: ixI^L, ixO^L
1197  double precision, intent(in) :: qdt, dtfactor, x(ixI^S, 1:ndim)
1198  double precision, intent(inout) :: wCT(ixI^S, 1:nw), w(ixI^S, 1:nw)
1199  ! to change and to set as a parameter in the parfile once the possibility to
1200  ! solve the equations in an angular momentum conserving form has been
1201  ! implemented (change tvdlf.t eg)
1202  double precision :: pth(ixI^S), source(ixI^S), minrho
1203  integer :: iw,idir, h1x^L{^NOONED, h2x^L}
1204  integer :: mr_,mphi_ ! Polar var. names
1205  integer :: irho, ifluid, n_fluids
1206  double precision :: exp_factor(ixI^S), del_exp_factor(ixI^S), exp_factor_primitive(ixI^S)
1207 
1208  if (hd_dust) then
1209  n_fluids = 1 + dust_n_species
1210  else
1211  n_fluids = 1
1212  end if
1213 
1214  select case (coordinate)
1215 
1216  case(cartesian_expansion)
1217  !the user provides the functions of exp_factor and del_exp_factor
1218  if(associated(usr_set_surface)) call usr_set_surface(ixi^l,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
1219  call hd_get_pthermal(wct, x, ixi^l, ixo^l, source)
1220  source(ixo^s) = source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1221  w(ixo^s,mom(1)) = w(ixo^s,mom(1)) + qdt*source(ixo^s)
1222 
1223  case (cylindrical)
1224  do ifluid = 0, n_fluids-1
1225  ! s[mr]=(pthermal+mphi**2/rho)/radius
1226  if (ifluid == 0) then
1227  ! gas
1228  irho = rho_
1229  mr_ = mom(r_)
1230  if(phi_>0) mphi_ = mom(phi_)
1231  call hd_get_pthermal(wct, x, ixi^l, ixo^l, source)
1232  minrho = 0.0d0
1233  else
1234  ! dust : no pressure
1235  irho = dust_rho(ifluid)
1236  mr_ = dust_mom(r_, ifluid)
1237  if(phi_>0) mphi_ = dust_mom(phi_, ifluid)
1238  source(ixi^s) = zero
1239  minrho = 0.0d0
1240  end if
1241  if (phi_ > 0) then
1242  where (wct(ixo^s, irho) > minrho)
1243  source(ixo^s) = source(ixo^s) + wct(ixo^s, mphi_)**2 / wct(ixo^s, irho)
1244  w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1245  end where
1246  ! s[mphi]=(-mphi*mr/rho)/radius
1247  where (wct(ixo^s, irho) > minrho)
1248  source(ixo^s) = -wct(ixo^s, mphi_) * wct(ixo^s, mr_) / wct(ixo^s, irho)
1249  w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * source(ixo^s) / x(ixo^s, r_)
1250  end where
1251  else
1252  ! s[mr]=2pthermal/radius
1253  w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1254  end if
1255  end do
1256  case (spherical)
1257  if (hd_dust) then
1258  call mpistop("Dust geom source terms not implemented yet with spherical geometries")
1259  end if
1260  mr_ = mom(r_)
1261  if(phi_>0) mphi_ = mom(phi_)
1262  h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
1263  ! s[mr]=((mtheta**2+mphi**2)/rho+2*p)/r
1264  call hd_get_pthermal(wct, x, ixi^l, ixo^l, pth)
1265  source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1266  *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
1267  /block%dvolume(ixo^s)
1268  if (ndir > 1) then
1269  do idir = 2, ndir
1270  source(ixo^s) = source(ixo^s) + wct(ixo^s, mom(idir))**2 / wct(ixo^s, rho_)
1271  end do
1272  end if
1273  w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, 1)
1274 
1275  {^nooned
1276  ! s[mtheta]=-(mr*mtheta/rho)/r+cot(theta)*(mphi**2/rho+p)/r
1277  source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1278  * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
1279  / block%dvolume(ixo^s)
1280  if (ndir == 3) then
1281  source(ixo^s) = source(ixo^s) + (wct(ixo^s, mom(3))**2 / wct(ixo^s, rho_)) / tan(x(ixo^s, 2))
1282  end if
1283  source(ixo^s) = source(ixo^s) - (wct(ixo^s, mom(2)) * wct(ixo^s, mr_)) / wct(ixo^s, rho_)
1284  w(ixo^s, mom(2)) = w(ixo^s, mom(2)) + qdt * source(ixo^s) / x(ixo^s, 1)
1285 
1286  if (ndir == 3) then
1287  ! s[mphi]=-(mphi*mr/rho)/r-cot(theta)*(mtheta*mphi/rho)/r
1288  source(ixo^s) = -(wct(ixo^s, mom(3)) * wct(ixo^s, mr_)) / wct(ixo^s, rho_)&
1289  - (wct(ixo^s, mom(2)) * wct(ixo^s, mom(3))) / wct(ixo^s, rho_) / tan(x(ixo^s, 2))
1290  w(ixo^s, mom(3)) = w(ixo^s, mom(3)) + qdt * source(ixo^s) / x(ixo^s, 1)
1291  end if
1292  }
1293  end select
1294 
1295  if (hd_viscosity) call visc_add_source_geom(qdt,ixi^l,ixo^l,wct,w,x)
1296 
1297  if (hd_rotating_frame) then
1298  if (hd_dust) then
1299  call mpistop("Rotating frame not implemented yet with dust")
1300  else
1301  call rotating_frame_add_source(qdt,dtfactor,ixi^l,ixo^l,wct,w,x)
1302  end if
1303  end if
1304 
1305  end subroutine hd_add_source_geom
1306 
1307  ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
1308  subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1313  use mod_usr_methods, only: usr_gravity
1315  use mod_cak_force, only: cak_add_source
1316 
1317  integer, intent(in) :: ixI^L, ixO^L
1318  double precision, intent(in) :: qdt, dtfactor
1319  double precision, intent(in) :: wCT(ixI^S, 1:nw),wCTprim(ixI^S,1:nw), x(ixI^S, 1:ndim)
1320  double precision, intent(inout) :: w(ixI^S, 1:nw)
1321  logical, intent(in) :: qsourcesplit
1322  logical, intent(inout) :: active
1323 
1324  double precision :: gravity_field(ixI^S, 1:ndim)
1325  integer :: idust, idim
1326 
1327  if(hd_dust .and. .not. use_imex_scheme) then
1328  call dust_add_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1329  end if
1330 
1331  if(hd_radiative_cooling) then
1332  call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1333  qsourcesplit,active, rc_fl)
1334  end if
1335 
1336  if(hd_viscosity) then
1337  call viscosity_add_source(qdt,ixi^l,ixo^l,wct,w,x,&
1338  hd_energy,qsourcesplit,active)
1339  end if
1340 
1341  if (hd_gravity) then
1342  call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1343  hd_energy,.false.,qsourcesplit,active)
1344 
1345  if (hd_dust .and. qsourcesplit .eqv. grav_split) then
1346  active = .true.
1347 
1348  call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1349  do idust = 1, dust_n_species
1350  do idim = 1, ndim
1351  w(ixo^s, dust_mom(idim, idust)) = w(ixo^s, dust_mom(idim, idust)) &
1352  + qdt * gravity_field(ixo^s, idim) * wct(ixo^s, dust_rho(idust))
1353  end do
1354  end do
1355  end if
1356  end if
1357 
1358  if (hd_cak_force) then
1359  call cak_add_source(qdt,ixi^l,ixo^l,wct,w,x,hd_energy,qsourcesplit,active)
1360  end if
1361 
1362  if(hd_partial_ionization) then
1363  if(.not.qsourcesplit) then
1364  active = .true.
1365  call hd_update_temperature(ixi^l,ixo^l,wct,w,x)
1366  end if
1367  end if
1368 
1369  end subroutine hd_add_source
1370 
1371  subroutine hd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
1373  use mod_dust, only: dust_get_dt
1375  use mod_viscosity, only: viscosity_get_dt
1376  use mod_gravity, only: gravity_get_dt
1377  use mod_cak_force, only: cak_get_dt
1378 
1379  integer, intent(in) :: ixI^L, ixO^L
1380  double precision, intent(in) :: dx^D, x(ixI^S, 1:^ND)
1381  double precision, intent(in) :: w(ixI^S, 1:nw)
1382  double precision, intent(inout) :: dtnew
1383 
1384  dtnew = bigdouble
1385 
1386  if(hd_dust) then
1387  call dust_get_dt(w, ixi^l, ixo^l, dtnew, dx^d, x)
1388  end if
1389 
1390  if(hd_radiative_cooling) then
1391  call cooling_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x,rc_fl)
1392  end if
1393 
1394  if(hd_viscosity) then
1395  call viscosity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1396  end if
1397 
1398  if(hd_gravity) then
1399  call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1400  end if
1401 
1402  if (hd_cak_force) then
1403  call cak_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1404  end if
1405 
1406  end subroutine hd_get_dt
1407 
1408  function hd_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
1409  use mod_global_parameters, only: nw, ndim
1410  integer, intent(in) :: ixi^l, ixo^l
1411  double precision, intent(in) :: w(ixi^s, nw)
1412  double precision :: ke(ixo^s)
1413  double precision, intent(in), optional :: inv_rho(ixo^s)
1414 
1415  if (present(inv_rho)) then
1416  ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * inv_rho
1417  else
1418  ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) / w(ixo^s, rho_)
1419  end if
1420  end function hd_kin_en
1421 
1422  function hd_inv_rho(w, ixI^L, ixO^L) result(inv_rho)
1423  use mod_global_parameters, only: nw, ndim
1424  integer, intent(in) :: ixi^l, ixo^l
1425  double precision, intent(in) :: w(ixi^s, nw)
1426  double precision :: inv_rho(ixo^s)
1427 
1428  ! Can make this more robust
1429  inv_rho = 1.0d0 / w(ixo^s, rho_)
1430  end function hd_inv_rho
1431 
1432  subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1433  ! handles hydro (density,pressure,velocity) bootstrapping
1434  ! any negative dust density is flagged as well (and throws an error)
1435  ! small_values_method=replace also for dust
1437  use mod_small_values
1439  logical, intent(in) :: primitive
1440  integer, intent(in) :: ixI^L,ixO^L
1441  double precision, intent(inout) :: w(ixI^S,1:nw)
1442  double precision, intent(in) :: x(ixI^S,1:ndim)
1443  character(len=*), intent(in) :: subname
1444 
1445  integer :: n,idir
1446  logical :: flag(ixI^S,1:nw)
1447 
1448  call hd_check_w(primitive, ixi^l, ixo^l, w, flag)
1449 
1450  if (any(flag)) then
1451  select case (small_values_method)
1452  case ("replace")
1453  where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
1454  do idir = 1, ndir
1455  if(small_values_fix_iw(mom(idir))) then
1456  where(flag(ixo^s,rho_)) w(ixo^s, mom(idir)) = 0.0d0
1457  end if
1458  end do
1459  if(hd_energy)then
1460  if(small_values_fix_iw(e_)) then
1461  if(primitive) then
1462  where(flag(ixo^s,rho_)) w(ixo^s, p_) = small_pressure
1463  else
1464  where(flag(ixo^s,rho_)) w(ixo^s, e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1465  endif
1466  end if
1467  endif
1468 
1469  if(hd_energy) then
1470  if(primitive) then
1471  where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
1472  else
1473  where(flag(ixo^s,e_))
1474  ! Add kinetic energy
1475  w(ixo^s,e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1476  end where
1477  end if
1478  end if
1479 
1480  if(hd_dust)then
1481  do n=1,dust_n_species
1482  where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1483  do idir = 1, ndir
1484  where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1485  enddo
1486  enddo
1487  endif
1488  case ("average")
1489  if(primitive)then
1490  ! averaging for all primitive fields, including dust
1491  call small_values_average(ixi^l, ixo^l, w, x, flag)
1492  else
1493  ! do averaging of density
1494  call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
1495  if(hd_energy) then
1496  ! do averaging of pressure
1497  w(ixi^s,p_)=(hd_gamma-1.d0)*(w(ixi^s,e_) &
1498  -0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_))
1499  call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
1500  w(ixi^s,e_)=w(ixi^s,p_)/(hd_gamma-1.d0) &
1501  +0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_)
1502  end if
1503  if(hd_dust)then
1504  do n=1,dust_n_species
1505  where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1506  do idir = 1, ndir
1507  where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1508  enddo
1509  enddo
1510  endif
1511  endif
1512  case default
1513  if(.not.primitive) then
1514  !convert w to primitive
1515  ! Calculate pressure = (gamma-1) * (e-ek)
1516  if(hd_energy) then
1517  w(ixo^s,p_)=(hd_gamma-1.d0)*(w(ixo^s,e_)-hd_kin_en(w,ixi^l,ixo^l))
1518  end if
1519  ! Convert gas momentum to velocity
1520  do idir = 1, ndir
1521  w(ixo^s, mom(idir)) = w(ixo^s, mom(idir))/w(ixo^s,rho_)
1522  end do
1523  end if
1524  ! NOTE: dust entries may still have conserved values here
1525  call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1526  end select
1527  end if
1528  end subroutine hd_handle_small_values
1529 
1530  subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1533  integer, intent(in) :: ixI^L, ixO^L
1534  double precision, intent(in) :: w(ixI^S,1:nw)
1535  double precision, intent(in) :: x(ixI^S,1:ndim)
1536  double precision, intent(out):: Rfactor(ixI^S)
1537 
1538  double precision :: iz_H(ixO^S),iz_He(ixO^S)
1539 
1540  call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1541  ! assume the first and second ionization of Helium have the same degree
1542  rfactor(ixo^s)=(1.d0+iz_h(ixo^s)+0.1d0*(1.d0+iz_he(ixo^s)*(1.d0+iz_he(ixo^s))))/2.3d0
1543 
1545 
1546  subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1548  integer, intent(in) :: ixI^L, ixO^L
1549  double precision, intent(in) :: w(ixI^S,1:nw)
1550  double precision, intent(in) :: x(ixI^S,1:ndim)
1551  double precision, intent(out):: Rfactor(ixI^S)
1552 
1553  rfactor(ixo^s)=rr
1554 
1555  end subroutine rfactor_from_constant_ionization
1556 
1557  subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1560 
1561  integer, intent(in) :: ixI^L, ixO^L
1562  double precision, intent(in) :: wCT(ixI^S,1:nw), x(ixI^S,1:ndim)
1563  double precision, intent(inout) :: w(ixI^S,1:nw)
1564 
1565  double precision :: iz_H(ixO^S),iz_He(ixO^S), pth(ixI^S)
1566 
1567  call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1568 
1569  call hd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1570 
1571  w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1572  he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1573 
1574  end subroutine hd_update_temperature
1575 
1576 end module mod_hd_phys
Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices iw=iwmin......
Module with basic data types used in amrvac.
integer, parameter std_len
Default length for strings.
Module to include CAK radiation line force in (magneto)hydrodynamic models Computes both the force fr...
Definition: mod_cak_force.t:27
subroutine cak_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
Check time step for total radiation contribution.
subroutine cak_init(phys_gamma)
Initialize the module.
Definition: mod_cak_force.t:98
subroutine cak_add_source(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Definition: mod_comm_lib.t:208
Module for physical and numeric constants.
Definition: mod_constants.t:2
double precision, parameter bigdouble
A very large real number.
Definition: mod_constants.t:11
Module for including dust species, which interact with the gas through a drag force.
Definition: mod_dust.t:3
subroutine, public dust_check_w(ixIL, ixOL, w, flag)
Definition: mod_dust.t:194
subroutine, public dust_get_flux(w, x, ixIL, ixOL, idim, f)
Definition: mod_dust.t:252
subroutine, public dust_implicit_update(dtfactor, qdt, qtC, psb, psa)
Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
Definition: mod_dust.t:651
subroutine, public dust_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
Get dt related to dust and gas stopping time (Laibe 2011)
Definition: mod_dust.t:897
subroutine, public dust_add_source(qdt, ixIL, ixOL, wCT, w, x, qsourcesplit, active)
w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
Definition: mod_dust.t:529
subroutine, public dust_to_primitive(ixIL, ixOL, w, x)
Definition: mod_dust.t:228
integer, dimension(:, :), allocatable, public, protected dust_mom
Indices of the dust momentum densities.
Definition: mod_dust.t:22
subroutine, public dust_get_cmax(w, x, ixIL, ixOL, idim, cmax, cmin)
Definition: mod_dust.t:1010
integer, public, protected dust_n_species
The number of dust species.
Definition: mod_dust.t:12
integer, dimension(:), allocatable, public, protected dust_rho
Indices of the dust densities.
Definition: mod_dust.t:19
subroutine, public dust_get_flux_prim(w, x, ixIL, ixOL, idim, f)
Definition: mod_dust.t:275
subroutine, public dust_check_params()
Definition: mod_dust.t:154
subroutine, public dust_evaluate_implicit(qtC, psa)
inplace update of psa==>F_im(psa)
Definition: mod_dust.t:582
subroutine, public dust_init(g_rho, g_mom, g_energy)
Definition: mod_dust.t:95
subroutine, public dust_to_conserved(ixIL, ixOL, w, x)
Definition: mod_dust.t:208
Module for flux conservation near refinement boundaries.
Module with geometry-related routines (e.g., divergence, curl)
Definition: mod_geometry.t:2
integer coordinate
Definition: mod_geometry.t:7
integer, parameter spherical
Definition: mod_geometry.t:11
integer, parameter cylindrical
Definition: mod_geometry.t:10
integer, parameter cartesian_expansion
Definition: mod_geometry.t:12
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.
logical h_correction
If true, do H-correction to fix the carbuncle problem at grid-aligned shocks.
double precision small_pressure
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 global_time
The global simulation time.
double precision unit_mass
Physical scaling factor for mass.
logical use_imex_scheme
whether IMEX in use or not
integer, dimension(3, 3) kr
Kronecker delta tensor.
integer it
Number of time steps taken.
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.
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
integer, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
double precision unit_velocity
Physical scaling factor for velocity.
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
logical phys_trac
Use TRAC (Johnston 2019 ApJL, 873, L22) for MHD or 1D HD.
logical crash
Save a snapshot before crash a run met unphysical values.
double precision small_density
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
double precision, dimension(ndim) dxlevel
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module for including gravity in (magneto)hydrodynamics simulations.
Definition: mod_gravity.t:2
logical grav_split
source split or not
Definition: mod_gravity.t:6
subroutine gravity_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
Definition: mod_gravity.t:81
subroutine gravity_add_source(qdt, ixIL, ixOL, wCT, wCTprim, w, x, energy, rhov, 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:45
subroutine gravity_init()
Initialize the module.
Definition: mod_gravity.t:27
Hydrodynamics physics module.
Definition: mod_hd_phys.t:2
subroutine, public hd_check_params
Definition: mod_hd_phys.t:537
logical, public, protected hd_energy
Whether an energy equation is used.
Definition: mod_hd_phys.t:12
logical, public, protected hd_dust
Whether dust is added.
Definition: mod_hd_phys.t:24
subroutine rfactor_from_temperature_ionization(w, x, ixIL, ixOL, Rfactor)
Definition: mod_hd_phys.t:1531
integer, public, protected e_
Index of the energy density (-1 if not present)
Definition: mod_hd_phys.t:57
logical, public, protected hd_radiative_cooling
Whether radiative cooling is added.
Definition: mod_hd_phys.t:20
double precision, public, protected rr
Definition: mod_hd_phys.t:98
subroutine hd_get_flux(wC, w, x, ixIL, ixOL, idim, f)
Definition: mod_hd_phys.t:1132
double precision, public hd_gamma
The adiabatic index.
Definition: mod_hd_phys.t:69
subroutine rhos_to_e(ixIL, ixOL, w, x)
Definition: mod_hd_phys.t:748
subroutine hd_get_a2max(w, x, ixIL, ixOL, a2max)
Definition: mod_hd_phys.t:811
subroutine hd_add_source_geom(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
Add geometrical source terms to w.
Definition: mod_hd_phys.t:1190
integer, public, protected hd_trac_type
Definition: mod_hd_phys.t:79
subroutine hd_get_flux_cons(w, x, ixIL, ixOL, idim, f)
Definition: mod_hd_phys.t:1093
logical, public, protected hd_particles
Whether particles module is added.
Definition: mod_hd_phys.t:33
subroutine hd_te_images
Definition: mod_hd_phys.t:380
type(tc_fluid), allocatable, public tc_fl
Definition: mod_hd_phys.t:16
logical, public, protected hd_viscosity
Whether viscosity is added.
Definition: mod_hd_phys.t:27
subroutine rc_params_read(fl)
Definition: mod_hd_phys.t:493
subroutine hd_get_rho(w, x, ixIL, ixOL, rho)
Definition: mod_hd_phys.t:481
subroutine, public hd_check_w(primitive, ixIL, ixOL, w, flag)
Returns logical argument flag where values are ok.
Definition: mod_hd_phys.t:613
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
Definition: mod_hd_phys.t:66
subroutine hd_e_to_ei(ixIL, ixOL, w, x)
Transform total energy to internal energy.
Definition: mod_hd_phys.t:721
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
Definition: mod_hd_phys.t:94
subroutine hd_sts_set_source_tc_hd(ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux)
Definition: mod_hd_phys.t:400
subroutine hd_ei_to_e(ixIL, ixOL, w, x)
Transform internal energy to total energy.
Definition: mod_hd_phys.t:708
subroutine hd_add_source(qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
Definition: mod_hd_phys.t:1309
subroutine hd_write_info(fh)
Write this module's parameters to a snapsoht.
Definition: mod_hd_phys.t:139
integer, public, protected te_
Indices of temperature.
Definition: mod_hd_phys.t:63
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
Definition: mod_hd_phys.t:51
subroutine hd_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
Definition: mod_hd_phys.t:1372
double precision function hd_get_tc_dt_hd(w, ixIL, ixOL, dxD, x)
Definition: mod_hd_phys.t:412
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
Definition: mod_hd_phys.t:88
logical, public, protected hd_cak_force
Whether CAK radiation line force is activated.
Definition: mod_hd_phys.t:39
subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
Calculate cmax_idim = csound + abs(v_idim) within ixO^L.
Definition: mod_hd_phys.t:887
subroutine, public hd_phys_init()
Initialize the module.
Definition: mod_hd_phys.t:157
subroutine hd_get_temperature_from_eint(w, x, ixIL, ixOL, res)
Calculate temperature=p/rho when in e_ the internal energy is stored.
Definition: mod_hd_phys.t:1079
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
Definition: mod_hd_phys.t:54
logical, public, protected hd_thermal_conduction
Whether thermal conduction is added.
Definition: mod_hd_phys.t:15
logical, public, protected eq_state_units
Definition: mod_hd_phys.t:102
logical, public, protected hd_force_diagonal
Allows overruling default corner filling (for debug mode, since otherwise corner primitives fail)
Definition: mod_hd_phys.t:82
subroutine, public hd_get_pthermal(w, x, ixIL, ixOL, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
Definition: mod_hd_phys.t:1014
double precision, public hd_adiab
The adiabatic constant.
Definition: mod_hd_phys.t:72
subroutine hd_get_tcutoff(ixIL, ixOL, w, x, tco_local, Tmax_local)
get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
Definition: mod_hd_phys.t:834
integer, public, protected rho_
Index of the density (in the w array)
Definition: mod_hd_phys.t:48
subroutine tc_params_read_hd(fl)
Definition: mod_hd_phys.t:462
subroutine hd_get_v_idim(w, x, ixIL, ixOL, idim, v)
Calculate v_i = m_i / rho within ixO^L.
Definition: mod_hd_phys.t:764
logical, public, protected hd_partial_ionization
Whether plasma is partially ionized.
Definition: mod_hd_phys.t:45
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
Definition: mod_hd_phys.t:91
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
Definition: mod_hd_phys.t:85
logical, public, protected hd_gravity
Whether gravity is added.
Definition: mod_hd_phys.t:30
subroutine hd_physical_units
Definition: mod_hd_phys.t:560
subroutine rfactor_from_constant_ionization(w, x, ixIL, ixOL, Rfactor)
Definition: mod_hd_phys.t:1547
subroutine, public hd_to_conserved(ixIL, ixOL, w, x)
Transform primitive variables into conservative ones.
Definition: mod_hd_phys.t:642
type(rc_fluid), allocatable, public rc_fl
Definition: mod_hd_phys.t:21
subroutine, public hd_to_primitive(ixIL, ixOL, w, x)
Transform conservative variables into primitive ones.
Definition: mod_hd_phys.t:674
logical, public, protected hd_trac
Whether TRAC method is used.
Definition: mod_hd_phys.t:78
integer, public, protected hd_n_tracer
Number of tracer species.
Definition: mod_hd_phys.t:42
type(te_fluid), allocatable, public te_fl_hd
Definition: mod_hd_phys.t:17
subroutine e_to_rhos(ixIL, ixOL, w, x)
Definition: mod_hd_phys.t:733
subroutine hd_get_cmax(w, x, ixIL, ixOL, idim, cmax)
Calculate cmax_idim = csound + abs(v_idim) within ixO^L.
Definition: mod_hd_phys.t:790
logical, public, protected hd_rotating_frame
Whether rotating frame is activated.
Definition: mod_hd_phys.t:36
subroutine hd_get_v(w, x, ixIL, ixOL, v)
Calculate velocity vector v_i = m_i / rho within ixO^L.
Definition: mod_hd_phys.t:774
double precision function, dimension(ixo^s), public hd_kin_en(w, ixIL, ixOL, inv_rho)
Definition: mod_hd_phys.t:1409
subroutine hd_update_temperature(ixIL, ixOL, wCT, w, x)
Definition: mod_hd_phys.t:1558
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
Definition: mod_hd_phys.t:60
subroutine hd_tc_handle_small_e(w, x, ixIL, ixOL, step)
Definition: mod_hd_phys.t:427
subroutine hd_get_temperature_from_etot(w, x, ixIL, ixOL, res)
Calculate temperature=p/rho when in e_ the total energy is stored.
Definition: mod_hd_phys.t:1064
subroutine hd_handle_small_values(primitive, w, x, ixIL, ixOL, subname)
Definition: mod_hd_phys.t:1433
subroutine, public hd_get_csound2(w, x, ixIL, ixOL, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. csound2=gamma*p/rho.
Definition: mod_hd_phys.t:1001
double precision function, dimension(ixo^s) hd_inv_rho(w, ixIL, ixOL)
Definition: mod_hd_phys.t:1423
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.
Definition: mod_particles.t:2
subroutine particles_init()
Initialize particle data and parameters.
Definition: mod_particles.t:15
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(fl, read_params)
subroutine cooling_get_dt(w, ixIL, ixOL, dtnew, dxD, x, fl)
subroutine radiative_cooling_init_params(phys_gamma, He_abund)
Radiative cooling initialization.
subroutine radiative_cooling_add_source(qdt, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active, fl)
Module for including rotating frame in (magneto)hydrodynamics simulations The rotation vector is assu...
subroutine rotating_frame_add_source(qdt, dtfactor, ixIL, ixOL, wCT, w, x)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
subroutine rotating_frame_init()
Initialize the module.
Module for handling problematic values in simulations, such as negative pressures.
logical, public trace_small_values
trace small values in the source file using traceback flag of compiler
subroutine, public small_values_average(ixIL, ixOL, w, x, w_flag, windex)
logical, dimension(:), allocatable, public small_values_fix_iw
Whether to apply small value fixes to certain variables.
subroutine, public small_values_error(wprim, x, ixIL, ixOL, w_flag, subname)
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 Adaptation of mod_thermal_conduction for the mod_supertimestepping ...
subroutine, public sts_set_source_tc_hd(ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
subroutine, public tc_get_hd_params(fl, read_hd_params)
Read tc module parameters from par file: MHD case.
subroutine tc_init_params(phys_gamma)
double precision function, public get_tc_dt_hd(w, ixIL, ixOL, dxD, x, fl)
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_surface), pointer usr_set_surface
procedure(phys_gravity), pointer usr_gravity
procedure(hd_pthermal), pointer usr_set_pthermal
The module add viscous source terms and check time step.
Definition: mod_viscosity.t:10
subroutine viscosity_get_dt(w, ixIL, ixOL, dtnew, dxD, x)
subroutine viscosity_add_source(qdt, ixIL, ixOL, wCT, w, x, energy, qsourcesplit, active)
Definition: mod_viscosity.t:89
subroutine viscosity_init(phys_wider_stencil, phys_req_diagonal)
Initialize the module.
Definition: mod_viscosity.t:58
subroutine visc_add_source_geom(qdt, ixIL, ixOL, wCT, w, x)
subroutine, public visc_get_flux_prim(w, x, ixIL, ixOL, idim, f, energy)