MPI-AMRVAC 3.2
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_hd_phys.t
Go to the documentation of this file.
1!> Hydrodynamics physics module
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 dust is added using and implicit update in IMEX
27 logical, public, protected :: hd_dust_implicit = .false.
28
29 !> Whether radiation-gas interaction is handled using flux limited diffusion
30 logical, public, protected :: hd_radiation_fld = .false.
31
32 !> Whether viscosity is added
33 logical, public, protected :: hd_viscosity = .false.
34
35 !> Whether gravity is added
36 logical, public, protected :: hd_gravity = .false.
37
38 !> Whether particles module is added
39 logical, public, protected :: hd_particles = .false.
40
41 !> Whether rotating frame is activated
42 logical, public, protected :: hd_rotating_frame = .false.
43
44 !> Whether CAK radiation line force is activated
45 logical, public, protected :: hd_cak_force = .false.
46
47 !> Number of tracer species
48 integer, public, protected :: hd_n_tracer = 0
49
50 !> Whether plasma is partially ionized
51 logical, public, protected :: hd_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 !> Indices of the momentum density for the form of better vectorization
60 integer, public, protected :: ^c&m^C_
61
62 !> Indices of the tracers
63 integer, allocatable, public, protected :: tracer(:)
64
65 !> Index of the energy density (-1 if not present)
66 integer, public, protected :: e_
67
68 !> Index of the gas pressure (-1 if not present) should equal e_
69 integer, public, protected :: p_
70
71 !> Index of the radiation energy (when fld active)
72 integer, public, protected :: r_e
73
74 !> Indices of temperature
75 integer, public, protected :: te_
76
77 !> Index of the cutoff temperature for the TRAC method
78 integer, public, protected :: tcoff_
79
80 !> The adiabatic index
81 double precision, public :: hd_gamma = 5.d0/3.0d0
82
83 !> gamma minus one and its inverse
84 double precision :: gamma_1, inv_gamma_1
85
86 !> The adiabatic constant
87 double precision, public :: hd_adiab = 1.0d0
88
89 !> Whether TRAC method is used
90 logical, public, protected :: hd_trac = .false.
91 integer, public, protected :: hd_trac_type = 1
92
93 !> Helium abundance over Hydrogen
94 double precision, public, protected :: he_abundance=0.1d0
95 !> Ionization fraction of H
96 !> H_ion_fr = H+/(H+ + H)
97 double precision, public, protected :: h_ion_fr=1d0
98 !> Ionization fraction of He
99 !> He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
100 double precision, public, protected :: he_ion_fr=1d0
101 !> Ratio of number He2+ / number He+ + He2+
102 !> He_ion_fr2 = He2+/(He2+ + He+)
103 double precision, public, protected :: he_ion_fr2=1d0
104 ! used for eq of state when it is not defined by units,
105 ! the units do not contain terms related to ionization fraction
106 ! and it is p = RR * rho * T
107 double precision, public, protected :: rr=1d0
108 ! remove the below flag and assume default value = .false.
109 ! when eq state properly implemented everywhere
110 ! and not anymore through units
111 logical, public, protected :: eq_state_units = .true.
112
113 procedure(sub_get_pthermal), pointer :: hd_get_rfactor => null()
114 ! Public methods
115 public :: hd_phys_init
116 public :: hd_kin_en
117 public :: hd_get_pthermal
118 public :: hd_get_csound2
119 public :: hd_to_conserved
120 public :: hd_to_primitive
121 public :: hd_check_params
122 public :: hd_check_w
123 public :: hd_e_to_ei
124 public :: hd_ei_to_e
125 ! in FLD used as phys_get_Rfactor
126 public :: hd_get_rfactor
127 ! Begin: following relevant for radiative hydro using FLD
128 ! first four are local and only of interest for mod_usr applications
129 ! where they can be used in diagnostics
130 ! NOTE those with _prim expect primitives on entry
132 public :: hd_get_csrad2
133 public :: hd_get_trad
135 ! the following used in FLD modules
136 ! as pointer phys_get_csrad2
137 public :: hd_get_csrad2_prim
138 ! the following used in FLD modules
139 ! as pointer phys_get_tgas
141 ! End: following relevant for radiative hydro using FLD
142
143contains
144
145 !> Read this module's parameters from a file
146 subroutine hd_read_params(files)
148 character(len=*), intent(in) :: files(:)
149 integer :: n
150
151 namelist /hd_list/ hd_energy, hd_n_tracer, hd_gamma, hd_adiab, &
157
158 do n = 1, size(files)
159 open(unitpar, file=trim(files(n)), status="old")
160 read(unitpar, hd_list, end=111)
161111 close(unitpar)
162 end do
163
164 end subroutine hd_read_params
165
166 !> Write this module's parameters to a snapsoht
167 subroutine hd_write_info(fh)
169 integer, intent(in) :: fh
170 integer, parameter :: n_par = 1
171 double precision :: values(n_par)
172 character(len=name_len) :: names(n_par)
173 integer, dimension(MPI_STATUS_SIZE) :: st
174 integer :: er
175
176 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
177
178 names(1) = "gamma"
179 values(1) = hd_gamma
180 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
181 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
182 end subroutine hd_write_info
183
184 !> Initialize the module
185 subroutine hd_phys_init()
189 use mod_dust, only: dust_init
191 use mod_gravity, only: gravity_init
193 use mod_cak_force, only: cak_init
198 use mod_fld
199
200 integer :: itr, idir
201
202 call hd_read_params(par_files)
203
204 physics_type = "hd"
205 phys_energy = hd_energy
206 phys_total_energy = hd_energy
207 phys_internal_e = .false.
208 phys_gamma = hd_gamma
209 phys_partial_ionization=hd_partial_ionization
210
212 if(phys_trac) then
213 if(ndim .eq. 1) then
214 if(hd_trac_type .gt. 2) then
216 if(mype==0) write(*,*) 'WARNING: set hd_trac_type=1'
217 end if
219 else
220 phys_trac=.false.
221 if(mype==0) write(*,*) 'WARNING: set hd_trac=F when ndim>=2'
222 end if
223 end if
224
225 ! set default gamma for polytropic/isothermal process
226 if(.not.hd_energy) then
227 if(hd_thermal_conduction) then
229 if(mype==0) write(*,*) 'WARNING: set hd_thermal_conduction=F when hd_energy=F'
230 end if
231 if(hd_radiative_cooling) then
233 if(mype==0) write(*,*) 'WARNING: set hd_radiative_cooling=F when hd_energy=F'
234 end if
235 end if
236 if(.not.eq_state_units) then
237 if(hd_partial_ionization) then
239 if(mype==0) write(*,*) 'WARNING: set hd_partial_ionization=F when eq_state_units=F'
240 end if
241 end if
243
244 allocate(start_indices(number_species),stop_indices(number_species))
245
246 ! set the index of the first flux variable for species 1
247 start_indices(1)=1
248 ! Determine flux variables
249 rho_ = var_set_rho()
250
251 allocate(mom(ndir))
252 mom(:) = var_set_momentum(ndir)
253 m^c_=mom(^c);
254
255 ! Set index of energy variable
256 if (hd_energy) then
257 e_ = var_set_energy()
258 p_ = e_
259 else
260 e_ = -1
261 p_ = -1
262 end if
263
264 if(hd_radiation_fld)then
265 if(hd_cak_force)then
266 if(mype==0) then
267 write(*,*)'Warning: CAK force addition together with FLD radiation'
268 endif
269 endif
271 if(mype==0) then
272 write(*,*)'Warning: Optically thin cooling together with FLD radiation'
273 endif
274 endif
275 if(hd_dust.and.hd_dust_implicit)then
276 call mpistop('implicit dust addition not compatible with FLD radiation')
277 endif
278 if(si_unit)then
279 call mpistop('using FLD implies the use of cgs units')
280 endif
281 if(.not.hd_energy)then
282 call mpistop('using FLD implies the use of an energy equation, set hd_energy=T')
283 else
284 !> set added variable and equation for radiation energy
285 r_e = var_set_radiation_energy()
286 phys_get_tgas => hd_get_temperature_from_prim
287 phys_get_csrad2 => hd_get_csrad2_prim
288 !> Initiate radiation-closure module
290 endif
291 else
292 r_e=-1
293 endif
294
295 phys_get_dt => hd_get_dt
296 phys_get_cmax => hd_get_cmax
297 phys_get_tcutoff => hd_get_tcutoff
298 phys_get_cbounds => hd_get_cbounds
299 phys_get_flux => hd_get_flux
300 phys_add_source_geom => hd_add_source_geom
301 phys_add_source => hd_add_source
302 phys_to_conserved => hd_to_conserved
303 phys_to_primitive => hd_to_primitive
304 phys_check_params => hd_check_params
305 phys_check_w => hd_check_w
306 phys_get_pthermal => hd_get_pthermal
307 phys_get_v => hd_get_v
308 phys_get_rho => hd_get_rho
309 phys_write_info => hd_write_info
310 phys_handle_small_values => hd_handle_small_values
311
312 ! derive units from basic units
313 call hd_physical_units()
314
315 if (hd_dust) then
316 call dust_init(rho_, mom(:), e_)
317 endif
318
319 allocate(tracer(hd_n_tracer))
320
321 ! Set starting index of tracers
322 do itr = 1, hd_n_tracer
323 tracer(itr) = var_set_fluxvar("trc", "trp", itr, need_bc=.false.)
324 end do
325
326 ! set temperature as an auxiliary variable to get ionization degree
327 if(hd_partial_ionization) then
328 te_ = var_set_auxvar('Te','Te')
329 else
330 te_ = -1
331 end if
332
333 ! set number of variables which need update ghostcells
334 nwgc=nwflux+nwaux
335
336 ! set the index of the last flux variable for species 1
337 stop_indices(1)=nwflux
338
339 if(hd_trac) then
340 tcoff_ = var_set_wextra()
341 iw_tcoff=tcoff_
342 else
343 tcoff_ = -1
344 end if
345
346 ! choose Rfactor in ideal gas law
347 if(hd_partial_ionization) then
348 hd_get_rfactor=>rfactor_from_temperature_ionization
349 phys_update_temperature => hd_update_temperature
350 else if(associated(usr_rfactor)) then
352 else
353 hd_get_rfactor=>rfactor_from_constant_ionization
354 end if
355
356 phys_get_rfactor => hd_get_rfactor
357
358 ! initialize thermal conduction module
359 if (hd_thermal_conduction) then
360 if (.not. hd_energy) &
361 call mpistop("thermal conduction needs hd_energy=T")
362
363 call sts_init()
365
366 allocate(tc_fl)
367 call tc_get_hd_params(tc_fl,tc_params_read_hd)
368 call add_sts_method(hd_get_tc_dt_hd,hd_sts_set_source_tc_hd,e_,1,e_,1,.false.)
370 call set_error_handling_to_head(hd_tc_handle_small_e)
371 tc_fl%get_temperature_from_conserved => hd_get_temperature_from_etot
372 tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
373 tc_fl%get_rho => hd_get_rho
374 tc_fl%e_ = e_
375 end if
376
377 ! Initialize radiative cooling module
378 if (hd_radiative_cooling) then
379 if (.not. hd_energy) &
380 call mpistop("radiative cooling needs hd_energy=T")
382 allocate(rc_fl)
383 call radiative_cooling_init(rc_fl,rc_params_read)
384 rc_fl%get_rho => hd_get_rho
385 rc_fl%get_pthermal => hd_get_pthermal
386 rc_fl%get_var_Rfactor => hd_get_rfactor
387 rc_fl%e_ = e_
388 rc_fl%Tcoff_ = tcoff_
389 end if
390 allocate(te_fl_hd)
391 te_fl_hd%get_rho=> hd_get_rho
392 te_fl_hd%get_pthermal=> hd_get_pthermal
393 te_fl_hd%get_var_Rfactor => hd_get_rfactor
394{^ifthreed
395 phys_te_images => hd_te_images
396}
397 ! Initialize viscosity module
398 if (hd_viscosity) call viscosity_init(phys_wider_stencil)
399
400 ! Initialize gravity module
401 if (hd_gravity) call gravity_init()
402
403 ! Initialize rotating_frame module
405
406 ! Initialize CAK radiation force module
408
409
410 ! Check whether custom flux types have been defined
411 if (.not. allocated(flux_type)) then
412 allocate(flux_type(ndir, nw))
413 flux_type = flux_default
414 else if (any(shape(flux_type) /= [ndir, nw])) then
415 call mpistop("phys_check error: flux_type has wrong shape")
416 end if
417
418 nvector = 1 ! No. vector vars
419 allocate(iw_vector(nvector))
420 iw_vector(1) = mom(1) - 1
421 ! initialize ionization degree table
423
424 end subroutine hd_phys_init
425
426{^ifthreed
427 subroutine hd_te_images
430 select case(convert_type)
431 case('EIvtiCCmpi','EIvtuCCmpi')
433 case('ESvtiCCmpi','ESvtuCCmpi')
435 case('SIvtiCCmpi','SIvtuCCmpi')
437 case('WIvtiCCmpi','WIvtuCCmpi')
439 case default
440 call mpistop("Error in synthesize emission: Unknown convert_type")
441 end select
442 end subroutine hd_te_images
443}
444!!start th cond
445 ! wrappers for STS functions in thermal_conductivity module
446 ! which take as argument the tc_fluid (defined in the physics module)
447 subroutine hd_sts_set_source_tc_hd(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_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
457 end subroutine hd_sts_set_source_tc_hd
458
459 function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
460 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para/rho)
461 !and tc_k_para can depend on T=p/rho
464
465 integer, intent(in) :: ixi^l, ixo^l
466 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
467 double precision, intent(in) :: w(ixi^s,1:nw)
468 double precision :: dtnew
469
470 dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
471 end function hd_get_tc_dt_hd
472
473 subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
474 ! move this in a different routine as in mhd if needed in more places
477
478 integer, intent(in) :: ixi^l,ixo^l
479 double precision, intent(inout) :: w(ixi^s,1:nw)
480 double precision, intent(in) :: x(ixi^s,1:ndim)
481 integer, intent(in) :: step
482
483 integer :: idir
484 logical :: flag(ixi^s,1:nw)
485 character(len=140) :: error_msg
486
487 flag=.false.
488 where(w(ixo^s,e_)<small_e) flag(ixo^s,e_)=.true.
489 if(any(flag(ixo^s,e_))) then
490 select case (small_values_method)
491 case ("replace")
492 where(flag(ixo^s,e_)) w(ixo^s,e_)=small_e
493 case ("average")
494 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
495 case default
496 ! small values error shows primitive variables
497 w(ixo^s,e_)=w(ixo^s,e_)*(hd_gamma - 1.0d0)
498 do idir = 1, ndir
499 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,rho_)
500 end do
501 write(error_msg,"(a,i3)") "Thermal conduction step ", step
502 call small_values_error(w, x, ixi^l, ixo^l, flag, error_msg)
503 end select
504 end if
505 end subroutine hd_tc_handle_small_e
506
507 ! fill in tc_fluid fields from namelist
508 subroutine tc_params_read_hd(fl)
510 type(tc_fluid), intent(inout) :: fl
511 integer :: n
512 logical :: tc_saturate=.false.
513 double precision :: tc_k_para=0d0
514
515 namelist /tc_list/ tc_saturate, tc_k_para
516
517 do n = 1, size(par_files)
518 open(unitpar, file=trim(par_files(n)), status="old")
519 read(unitpar, tc_list, end=111)
520111 close(unitpar)
521 end do
522 fl%tc_saturate = tc_saturate
523 fl%tc_k_para = tc_k_para
524
525 end subroutine tc_params_read_hd
526
527 subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
529 integer, intent(in) :: ixi^l, ixo^l
530 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
531 double precision, intent(out) :: rho(ixi^s)
532
533 rho(ixo^s) = w(ixo^s,rho_)
534
535 end subroutine hd_get_rho
536
537!!end th cond
538!!rad cool
539 subroutine rc_params_read(fl)
541 use mod_constants, only: bigdouble
542 use mod_basic_types, only: std_len
543 type(rc_fluid), intent(inout) :: fl
544 integer :: n
545 ! list parameters
546 integer :: ncool = 4000
547
548 !> Name of cooling curve
549 character(len=std_len) :: coolcurve='JCcorona'
550
551 !> Fixed temperature not lower than tlow
552 logical :: tfix=.false.
553
554 !> Lower limit of temperature
555 double precision :: tlow=bigdouble
556
557 !> Add cooling source in a split way (.true.) or un-split way (.false.)
558 logical :: rc_split=.false.
559
560
561 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split
562
563 do n = 1, size(par_files)
564 open(unitpar, file=trim(par_files(n)), status="old")
565 read(unitpar, rc_list, end=111)
566111 close(unitpar)
567 end do
568
569 fl%ncool=ncool
570 fl%coolcurve=coolcurve
571 fl%tlow=tlow
572 fl%Tfix=tfix
573 fl%rc_split=rc_split
574 end subroutine rc_params_read
575!! end rad cool
576
579 use mod_geometry, only: coordinate
582 use mod_particles, only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
583 use mod_fld
584
585 ! Initialize particles module, put here so additional gridvars and user payloads are known
586 if (hd_particles) then
587 call particles_init()
588 end if
589
590 if (.not. hd_energy) then
591 if (hd_gamma <= 0.0d0) call mpistop ("Error: hd_gamma <= 0")
592 if (hd_adiab < 0.0d0) call mpistop ("Error: hd_adiab < 0")
594 else
595 if (hd_gamma <= 0.0d0 .or. hd_gamma == 1.0d0) &
596 call mpistop ("Error: hd_gamma <= 0 or hd_gamma == 1.0")
598 inv_gamma_1=1.d0/(hd_gamma-1.d0)
600 end if
601
602 if (hd_dust) call dust_check_params()
603
604 if(hd_dust_implicit) then
605 if(.not.use_imex_scheme)then
606 call mpistop('select IMEX scheme for implicit dust update')
607 endif
608 ! implicit dust update
609 phys_implicit_update => dust_implicit_update
610 phys_evaluate_implicit => dust_evaluate_implicit
611 endif
612 if(hd_radiation_fld) then
613 if(.not.use_imex_scheme)then
614 call mpistop('select IMEX scheme for FLD radiation use')
615 endif
616 if(use_multigrid)then
617 call phys_set_mg_bounds()
618 else
619 if(.not.fld_no_mg)call mpistop('multigrid must have BCs for IMEX and FLD radiation use')
620 endif
621 if(mype==0)then
622 write(*,*)'========================'
623 write(*,*)'Using FLD with settings:'
624 write(*,*)'Using FLD with settings: hd_radiation_fld=',hd_radiation_fld
625 write(*,*)'Using FLD with settings: fld_fluxlimiter=',fld_fluxlimiter
626 write(*,*)'Using FLD with settings: fld_interaction_method=',fld_interaction_method
627 write(*,*)'Using FLD with settings: fld_opacity_law=',fld_opacity_law
628 write(*,*)'Using FLD with settings: fld_kappa0=',fld_kappa0
629 write(*,*)'Using FLD with settings: fld_opal_table=',fld_opal_table
630 write(*,*)'Using FLD with settings: fld_Radforce_split=',fld_radforce_split
631 write(*,*)'Using FLD with settings: fld_bisect_tol=',fld_bisect_tol
632 write(*,*)'Using FLD with settings: fld_diff_tol=',fld_diff_tol
633 write(*,*)'Using FLD with settings: nth_for_diff_mg=',nth_for_diff_mg
634 write(*,*)' FLD has use_imex_scheme and use_multigrid=',use_imex_scheme,use_multigrid
635 write(*,*)'========================'
636 endif
637 endif
638 if(mype==0)then
639 write(*,*)'====HD run with settings===================='
640 write(*,*)'Using mod_hd_phys with settings:'
641 write(*,*)'SI_unit=',si_unit
642 write(*,*)'Dimensionality :',ndim
643 write(*,*)'vector components:',ndir
644 write(*,*)'coordinate set to type,slab:',coordinate,slab
645 write(*,*)'number of variables nw=',nw
646 write(*,*)' start index iwstart=',iwstart
647 write(*,*)'number of vector variables=',nvector
648 write(*,*)'number of stagger variables nws=',nws
649 write(*,*)'number of variables with BCs=',nwgc
650 write(*,*)'number of vars with fluxes=',nwflux
651 write(*,*)'number of vars with flux + BC=',nwfluxbc
652 write(*,*)'number of auxiliary variables=',nwaux
653 write(*,*)'number of extra vars without flux=',nwextra
654 write(*,*)'number of extra vars for wextra=',nw_extra
655 write(*,*)'number of auxiliary I/O variables=',nwauxio
656 write(*,*)'number of hd_n_tracer=',hd_n_tracer
657 write(*,*)' hd_energy=',hd_energy
658 write(*,*)' hd_gravity=',hd_gravity
659 write(*,*)' hd_viscosity=',hd_viscosity
660 write(*,*)' hd_radiative_cooling=',hd_radiative_cooling
661 write(*,*)' hd_cak_force=',hd_cak_force
662 write(*,*)' hd_radiation_fld=',hd_radiation_fld
663 write(*,*)' hd_thermal_conduction=',hd_thermal_conduction
664 write(*,*)' hd_trac=',hd_trac
665 write(*,*)' hd_dust=',hd_dust
666 write(*,*)' hd_rotating_frame=',hd_rotating_frame
667 write(*,*)' hd_particles=',hd_particles
668 if(hd_particles) then
669 write(*,*) '*****Using particles: npayload,ngridvars :', npayload,ngridvars
670 write(*,*) '*****Using particles: nusrpayload :', nusrpayload
671 write(*,*) '*****Using particles: num_particles :', num_particles
672 write(*,*) '*****Using particles: physics_type_particles=',physics_type_particles
673 end if
674 write(*,*)'number of ghostcells=',nghostcells
675 write(*,*)'number due to phys_wider_stencil=',phys_wider_stencil
676 write(*,*)'==========================================='
677 endif
678
679 end subroutine hd_check_params
680
681 subroutine hd_physical_units
683 double precision :: mp,kb
684 double precision :: a,b
685 ! Derive scaling units
686 if(si_unit) then
687 mp=mp_si
688 kb=kb_si
689 else
690 mp=mp_cgs
691 kb=kb_cgs
692 end if
693 if(eq_state_units) then
694 a=1d0+4d0*he_abundance
695 if(hd_partial_ionization) then
696 b=1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0)
697 else
698 b=2d0+3d0*he_abundance
699 end if
700 rr=1d0
701 else
702 a=1d0
703 b=1d0
704 rr=(1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0))/(1d0+4d0*he_abundance)
705 end if
706 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
707 if(unit_density/=1.d0) then
709 else if(unit_numberdensity/=1.d0) then
711 end if
712 if(unit_temperature/=1.d0) then
715 if(unit_length/=1.d0) then
717 else if(unit_time/=1.d0) then
719 end if
720 else if(unit_pressure/=1.d0) then
723 if(unit_length/=1.d0) then
725 else if(unit_time/=1.d0) then
727 end if
728 else if(unit_velocity/=1.d0) then
731 if(unit_length/=1.d0) then
733 else if(unit_time/=1.d0) then
735 end if
736 else if(unit_time/=1.d0) then
740 end if
741 else if(unit_temperature/=1.d0) then
742 ! units of temperature and velocity are dependent
743 if(unit_pressure/=1.d0) then
747 if(unit_length/=1.d0) then
749 else if(unit_time/=1.d0) then
751 end if
752 end if
753 else if(unit_pressure/=1.d0) then
754 if(unit_velocity/=1.d0) then
758 if(unit_length/=1.d0) then
760 else if(unit_time/=1.d0) then
762 end if
763 else if(unit_time/=0.d0) then
768 end if
769 end if
771
772 !> Units for radiative flux and opacity, latter is used in FLD
775
776 end subroutine hd_physical_units
777
778 !> Returns logical argument flag where values are ok
779 subroutine hd_check_w(primitive, ixI^L, ixO^L, w, flag)
781 use mod_dust, only: dust_check_w
782
783 logical, intent(in) :: primitive
784 integer, intent(in) :: ixi^l, ixo^l
785 double precision, intent(in) :: w(ixi^s, nw)
786 logical, intent(inout) :: flag(ixi^s,1:nw)
787 double precision :: tmp(ixi^s)
788
789 flag=.false.
790
791 if (hd_energy) then
792 if (primitive) then
793 where(w(ixo^s, e_) < small_pressure) flag(ixo^s,e_) = .true.
794 else
795 tmp(ixo^s)=(hd_gamma-1.0d0)*(w(ixo^s,e_)-&
796 half*(^c&w(ixo^s,m^c_)**2+)/w(ixo^s,rho_))
797 where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_) = .true.
798 endif
799 if(hd_radiation_fld)then
800 where(w(ixo^s, r_e) < small_r_e) flag(ixo^s,r_e) = .true.
801 endif
802 end if
803
804 where(w(ixo^s, rho_) < small_density) flag(ixo^s,rho_) = .true.
805
806 if(hd_dust) call dust_check_w(ixi^l,ixo^l,w,flag)
807
808 end subroutine hd_check_w
809
810 !> Transform primitive variables into conservative ones
811 subroutine hd_to_conserved(ixI^L, ixO^L, w, x)
813 use mod_dust, only: dust_to_conserved
814 integer, intent(in) :: ixi^l, ixo^l
815 double precision, intent(inout) :: w(ixi^s, nw)
816 double precision, intent(in) :: x(ixi^s, 1:ndim)
817
818 integer :: ix^d
819
820 {do ix^db=ixomin^db,ixomax^db\}
821 if (hd_energy) then
822 ! Calculate total energy from pressure and kinetic energy
823 w(ix^d,e_)=w(ix^d, e_)*inv_gamma_1+&
824 half*(^c&w(ix^d,m^c_)**2+)*w(ix^d,rho_)
825 end if
826 ! Convert velocity to momentum
827 ^c&w(ix^d,m^c_)=w(ix^d,rho_)*w(ix^d,m^c_)\
828 {end do\}
829
830 if (hd_dust) then
831 call dust_to_conserved(ixi^l, ixo^l, w, x)
832 end if
833
834 end subroutine hd_to_conserved
835
836 !> Transform conservative variables into primitive ones
837 subroutine hd_to_primitive(ixI^L, ixO^L, w, x)
839 use mod_dust, only: dust_to_primitive
840 integer, intent(in) :: ixi^l, ixo^l
841 double precision, intent(inout) :: w(ixi^s, nw)
842 double precision, intent(in) :: x(ixi^s, 1:ndim)
843
844 double precision :: inv_rho
845 integer :: ix^d
846
847 if (fix_small_values) then
848 call hd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'hd_to_primitive')
849 end if
850
851 {do ix^db=ixomin^db,ixomax^db\}
852 inv_rho = 1.d0/w(ix^d,rho_)
853 ! Convert momentum to velocity
854 ^c&w(ix^d,m^c_)=w(ix^d,m^c_)*inv_rho\
855 ! Calculate pressure = (gamma-1) * (e-ek)
856 if(hd_energy) then
857 ! Compute pressure
858 w(ix^d,p_)=(hd_gamma-1.d0)*(w(ix^d,e_)&
859 -half*w(ix^d,rho_)*(^c&w(ix^d,m^c_)**2+))
860 end if
861 {end do\}
862
863 ! Convert dust momentum to dust velocity
864 if (hd_dust) then
865 call dust_to_primitive(ixi^l, ixo^l, w, x)
866 end if
867
868 end subroutine hd_to_primitive
869
870 !> Transform internal energy to total energy
871 subroutine hd_ei_to_e(ixI^L,ixO^L,w,x)
873 integer, intent(in) :: ixi^l, ixo^l
874 double precision, intent(inout) :: w(ixi^s, nw)
875 double precision, intent(in) :: x(ixi^s, 1:ndim)
876
877 ! Calculate total energy from internal and kinetic energy
878 w(ixo^s,e_)=w(ixo^s,e_)+half*(^c&w(ixo^s,m^c_)**2+)/w(ixo^s,rho_)
879
880 end subroutine hd_ei_to_e
881
882 !> Transform total energy to internal energy
883 subroutine hd_e_to_ei(ixI^L,ixO^L,w,x)
885 integer, intent(in) :: ixi^l, ixo^l
886 double precision, intent(inout) :: w(ixi^s, nw)
887 double precision, intent(in) :: x(ixi^s, 1:ndim)
888
889 ! Calculate ei = e - ek
890 w(ixo^s,e_)=w(ixo^s,e_)-half*(^c&w(ixo^s,m^c_)**2+)/w(ixo^s,rho_)
891
892 end subroutine hd_e_to_ei
893
894 !> Calculate v_i = m_i / rho within ixO^L
895 subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
897 integer, intent(in) :: ixi^l, ixo^l, idim
898 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
899 double precision, intent(out) :: v(ixi^s)
900
901 v(ixo^s) = w(ixo^s, mom(idim)) / w(ixo^s, rho_)
902 end subroutine hd_get_v_idim
903
904 !> Calculate velocity vector v_i = m_i / rho within ixO^L
905 subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
907
908 integer, intent(in) :: ixi^l, ixo^l
909 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:^nd)
910 double precision, intent(out) :: v(ixi^s,1:ndir)
911
912 integer :: idir
913
914 do idir=1,ndir
915 v(ixo^s,idir) = w(ixo^s, mom(idir)) / w(ixo^s, rho_)
916 end do
917
918 end subroutine hd_get_v
919
920 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
921 subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
925
926 integer, intent(in) :: ixi^l, ixo^l, idim
927 ! w in primitive form
928 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
929 double precision, intent(inout) :: cmax(ixi^s)
930
931 if(hd_energy) then
932 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(hd_gamma*w(ixo^s,p_)/w(ixo^s,rho_))
933 else
934 if (.not. associated(usr_set_pthermal)) then
935 cmax(ixo^s) = hd_adiab * w(ixo^s, rho_)**hd_gamma
936 else
937 call usr_set_pthermal(w,x,ixi^l,ixo^l,cmax)
938 end if
939 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(hd_gamma*cmax(ixo^s)/w(ixo^s,rho_))
940 end if
941
942 if (hd_dust) then
943 call dust_get_cmax_prim(w, x, ixi^l, ixo^l, idim, cmax)
944 end if
945 end subroutine hd_get_cmax
946
947 !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
948 subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
950 integer, intent(in) :: ixi^l,ixo^l
951 double precision, intent(in) :: x(ixi^s,1:ndim)
952 ! in primitive form
953 double precision, intent(inout) :: w(ixi^s,1:nw)
954 double precision, intent(out) :: tco_local, tmax_local
955
956 double precision, parameter :: trac_delta=0.25d0
957 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
958 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
959 integer :: jxo^l,hxo^l
960 integer :: jxp^l,hxp^l,ixp^l
961 logical :: lrlt(ixi^s)
962
963 {^ifoned
964 call hd_get_rfactor(w,x,ixi^l,ixi^l,te)
965 te(ixi^s)=w(ixi^s,p_)/(te(ixi^s)*w(ixi^s,rho_))
966
967 tco_local=zero
968 tmax_local=maxval(te(ixo^s))
969 select case(hd_trac_type)
970 case(0)
971 w(ixi^s,tcoff_)=3.d5/unit_temperature
972 case(1)
973 hxo^l=ixo^l-1;
974 jxo^l=ixo^l+1;
975 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
976 lrlt=.false.
977 where(lts(ixo^s) > trac_delta)
978 lrlt(ixo^s)=.true.
979 end where
980 if(any(lrlt(ixo^s))) then
981 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
982 end if
983 case(2)
984 !> iijima et al. 2021, LTRAC method
985 ltrc=1.5d0
986 ltrp=2.5d0
987 ixp^l=ixo^l^ladd1;
988 hxo^l=ixo^l-1;
989 jxo^l=ixo^l+1;
990 hxp^l=ixp^l-1;
991 jxp^l=ixp^l+1;
992 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
993 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
994 w(ixo^s,tcoff_)=te(ixo^s)*&
995 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
996 case default
997 call mpistop("hd_trac_type not allowed for 1D simulation")
998 end select
999 }
1000 end subroutine hd_get_tcutoff
1001
1002 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
1003 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
1005 use mod_dust, only: dust_get_cmax
1006 use mod_variables
1007
1008 integer, intent(in) :: ixi^l, ixo^l, idim
1009 ! conservative left and right status
1010 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1011 ! primitive left and right status
1012 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1013 double precision, intent(in) :: x(ixi^s, 1:ndim)
1014 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1015 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1016 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
1017
1018 double precision :: wmean(ixi^s,nw)
1019 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1020 integer :: ix^d
1021
1022 select case(boundspeed)
1023 case (1)
1024 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1025 ! Methods for Fluid Dynamics" by Toro.
1026
1027 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1028 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1029 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1030 umean(ixo^s)=(wlp(ixo^s,mom(idim))*tmp1(ixo^s)+wrp(ixo^s,mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1031
1032 if(hd_energy) then
1033 ! note usage of primitives here
1034 csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1035 csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1036 else
1037 ! note usage of conservatives here
1038 call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1039 call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1040 end if
1041
1042 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1043 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1044 (wrp(ixo^s,mom(idim))-wlp(ixo^s,mom(idim)))**2
1045
1046 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1047 if(present(cmin)) then
1048 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1049 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1050 if(h_correction) then
1051 {do ix^db=ixomin^db,ixomax^db\}
1052 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1053 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1054 {end do\}
1055 end if
1056 else
1057 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1058 end if
1059
1060 if (hd_dust) then
1061 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1062 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1063 end if
1064
1065 case (2)
1066 !if(hd_energy) then
1067 ! ! note usage of primitives here
1068 ! wmean(ixO^S,1:nwflux)=0.5d0*(wLp(ixO^S,1:nwflux)+wRp(ixO^S,1:nwflux))
1069 ! tmp1(ixO^S)=wmean(ixO^S,mom(idim))
1070 ! csoundR(ixO^S)=hd_gamma*wmean(ixO^S,p_)/wmean(ixO^S,rho_)
1071 !else
1072 ! note usage of conservatives here
1073 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1074 tmp1(ixo^s)=wmean(ixo^s,mom(idim))/wmean(ixo^s,rho_)
1075 call hd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1076 !endif
1077 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1078
1079 if(present(cmin)) then
1080 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1081 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1082 if(h_correction) then
1083 {do ix^db=ixomin^db,ixomax^db\}
1084 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1085 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1086 {end do\}
1087 end if
1088 else
1089 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1090 end if
1091
1092 if (hd_dust) then
1093 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1094 end if
1095 case (3)
1096 ! Miyoshi 2005 JCP 208, 315 equation (67)
1097 if(hd_energy) then
1098 ! note usage of primitives here
1099 csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1100 csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1101 else
1102 ! note usage of conservatives here
1103 call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1104 call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1105 end if
1106 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1107 if(present(cmin)) then
1108 cmin(ixo^s,1)=min(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))-csoundl(ixo^s)
1109 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1110 if(h_correction) then
1111 {do ix^db=ixomin^db,ixomax^db\}
1112 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1113 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1114 {end do\}
1115 end if
1116 else
1117 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1118 end if
1119 if (hd_dust) then
1120 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1121 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1122 end if
1123 end select
1124
1125 end subroutine hd_get_cbounds
1126
1127 !> Calculate the square of the thermal sound speed csound2 within ixO^L.
1128 !> csound2=gamma*p/rho
1129 subroutine hd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1131 integer, intent(in) :: ixi^l, ixo^l
1132 double precision, intent(in) :: w(ixi^s,nw)
1133 double precision, intent(in) :: x(ixi^s,1:ndim)
1134 double precision, intent(out) :: csound2(ixi^s)
1135
1136 call hd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1137 csound2(ixo^s)=hd_gamma*csound2(ixo^s)/w(ixo^s,rho_)
1138
1139 end subroutine hd_get_csound2
1140
1141 !> Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L
1142 subroutine hd_get_pthermal(w, x, ixI^L, ixO^L, pth)
1146
1147 integer, intent(in) :: ixi^l, ixo^l
1148 double precision, intent(in) :: w(ixi^s, 1:nw)
1149 double precision, intent(in) :: x(ixi^s, 1:ndim)
1150 double precision, intent(out):: pth(ixi^s)
1151 integer :: iw, ix^d
1152
1153 if (hd_energy) then
1154 pth(ixo^s) = (hd_gamma - 1.0d0) * (w(ixo^s, e_) - &
1155 hd_kin_en(w, ixi^l, ixo^l))
1156 else
1157 if (.not. associated(usr_set_pthermal)) then
1158 pth(ixo^s) = hd_adiab * w(ixo^s, rho_)**hd_gamma
1159 else
1160 call usr_set_pthermal(w,x,ixi^l,ixo^l,pth)
1161 end if
1162 end if
1163
1164 if (fix_small_values) then
1165 {do ix^db= ixo^lim^db\}
1166 if(pth(ix^d)<small_pressure) then
1167 pth(ix^d)=small_pressure
1168 endif
1169 {enddo^d&\}
1170 else if (check_small_values) then
1171 {do ix^db= ixo^lim^db\}
1172 if(pth(ix^d)<small_pressure) then
1173 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1174 " encountered when call hd_get_pthermal"
1175 write(*,*) "Iteration: ", it, " Time: ", global_time
1176 write(*,*) "Location: ", x(ix^d,:)
1177 write(*,*) "Cell number: ", ix^d
1178 do iw=1,nw
1179 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1180 end do
1181 ! use erroneous arithmetic operation to crash the run
1182 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1183 write(*,*) "Saving status at the previous time step"
1184 crash=.true.
1185 end if
1186 {enddo^d&\}
1187 end if
1188
1189 end subroutine hd_get_pthermal
1190
1191 !> Calculate modified squared sound speed for FLD
1192 !> NOTE: only for diagnostic purposes, unused subroutine
1193 subroutine hd_get_csrad2(w,x,ixI^L,ixO^L,csound)
1195
1196 integer, intent(in) :: ixi^l, ixo^l
1197 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
1198 double precision, intent(out):: csound(ixi^s)
1199
1200 double precision :: wprim(ixi^s, nw)
1201
1202 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1203 call hd_to_primitive(ixi^l,ixo^l,wprim,x)
1204 call hd_get_csrad2_prim(wprim,x,ixi^l,ixo^l,csound)
1205
1206 end subroutine hd_get_csrad2
1207
1208 !> Calculate modified squared sound speed for FLD
1209 !> NOTE: w is primitive on entry here!
1210 !> NOTE: used in FLD module as phys_get_csrad2
1211 subroutine hd_get_csrad2_prim(w,x,ixI^L,ixO^L,csound)
1213
1214 integer, intent(in) :: ixi^l, ixo^l
1215 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
1216 double precision, intent(out):: csound(ixi^s)
1217
1218 integer :: ix^d
1219 double precision :: inv_rho
1220 double precision :: prad_tensor(ixi^s, 1:ndim, 1:ndim)
1221 double precision :: prad_max(ixi^s)
1222
1223 call hd_get_pradiation_from_prim(w, x, ixi^l, ixo^l, prad_tensor)
1224
1225 {do ix^db=ixomin^db,ixomax^db \}
1226 inv_rho=1.d0/w(ix^d,rho_)
1227 prad_max(ix^d) = maxval(prad_tensor(ix^d,:,:))
1228 csound(ix^d)=(hd_gamma*w(ix^d,p_)+prad_max(ix^d))*inv_rho
1229 {end do\}
1230
1231 if(minval(csound(ixo^s))<smalldouble)then
1232 print *,'issue with squared speed and rad pressure'
1233 print *,minval(csound(ixo^s))
1234 print *,minval(prad_max(ixo^s))
1235 call mpistop("negative squared speed in get_csrad2 for dt")
1236 endif
1237
1238 end subroutine hd_get_csrad2_prim
1239
1240 !> Calculate radiation pressure within ixO^L
1241 !> NOTE: w is primitive on entry here!
1242 !> NOTE: used in FLD module as it is called from phys_get_csrad2
1243 subroutine hd_get_pradiation_from_prim(w, x, ixI^L, ixO^L, prad)
1245 use mod_fld
1246 integer, intent(in) :: ixi^l, ixo^l
1247 double precision, intent(in) :: w(ixi^s, 1:nw)
1248 double precision, intent(in) :: x(ixi^s, 1:ndim)
1249 double precision, intent(out):: prad(ixi^s, 1:ndim, 1:ndim)
1250
1251 call fld_get_radpress(w, x, ixi^l, ixo^l, prad)
1252
1253 end subroutine hd_get_pradiation_from_prim
1254
1255 !> calculates the sum of the gas pressure and max Prad tensor element
1256 !> NOTE: only for diagnostic purposes, unused subroutine
1257 subroutine hd_get_pthermal_plus_pradiation(w, x, ixI^L, ixO^L, pth_plus_prad)
1259 integer, intent(in) :: ixi^l, ixo^l
1260 double precision, intent(in) :: w(ixi^s, 1:nw)
1261 double precision, intent(in) :: x(ixi^s, 1:ndim)
1262 double precision, intent(out):: pth_plus_prad(ixi^s)
1263
1264 double precision :: wprim(ixi^s, 1:nw)
1265 double precision :: prad_tensor(ixi^s, 1:ndim, 1:ndim)
1266 double precision :: prad_max(ixi^s)
1267 integer :: ix^d
1268
1269 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1270 call hd_to_primitive(ixi^l,ixo^l,wprim,x)
1271 call hd_get_pradiation_from_prim(wprim, x, ixi^l, ixo^l, prad_tensor)
1272 {do ix^d = ixomin^d,ixomax^d\}
1273 prad_max(ix^d) = maxval(prad_tensor(ix^d,:,:))
1274 {enddo\}
1275 pth_plus_prad(ixo^s) = wprim(ixo^s,p_) + prad_max(ixo^s)
1276 end subroutine hd_get_pthermal_plus_pradiation
1277
1278 !> Calculates radiation temperature
1279 !> Note use of cgs units here in factor const_rad_a
1280 subroutine hd_get_trad(w, x, ixI^L, ixO^L, trad)
1282 use mod_constants
1283
1284 integer, intent(in) :: ixi^l, ixo^l
1285 double precision, intent(in) :: w(ixi^s, 1:nw)
1286 double precision, intent(in) :: x(ixi^s, 1:ndim)
1287 double precision, intent(out):: trad(ixi^s)
1288
1289 trad(ixi^s) = (w(ixi^s,r_e)*unit_pressure&
1290 /const_rad_a)**(1.d0/4.d0)/unit_temperature
1291
1292 end subroutine hd_get_trad
1293
1294 !> Calculate temperature=p/rho when in e_ the total energy is stored
1295 subroutine hd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1297 integer, intent(in) :: ixi^l, ixo^l
1298 double precision, intent(in) :: w(ixi^s, 1:nw)
1299 double precision, intent(in) :: x(ixi^s, 1:ndim)
1300 double precision, intent(out):: res(ixi^s)
1301
1302 double precision :: r(ixi^s)
1303
1304 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1305 call hd_get_pthermal(w, x, ixi^l, ixo^l, res)
1306 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1307 end subroutine hd_get_temperature_from_etot
1308
1309 !> Calculate temperature=p/rho when in e_ the internal energy is stored
1310 subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1312 integer, intent(in) :: ixi^l, ixo^l
1313 double precision, intent(in) :: w(ixi^s, 1:nw)
1314 double precision, intent(in) :: x(ixi^s, 1:ndim)
1315 double precision, intent(out):: res(ixi^s)
1316
1317 double precision :: r(ixi^s)
1318
1319 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1320 res(ixo^s) = (hd_gamma - 1.0d0) * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1321 end subroutine hd_get_temperature_from_eint
1322
1323 !> Calculate temperature=p/rho when in we work in primitives (hence e_ is p_)
1324 subroutine hd_get_temperature_from_prim(w, x, ixI^L, ixO^L, res)
1326 integer, intent(in) :: ixi^l, ixo^l
1327 double precision, intent(in) :: w(ixi^s, 1:nw)
1328 double precision, intent(in) :: x(ixi^s, 1:ndim)
1329 double precision, intent(out):: res(ixi^s)
1330
1331 double precision :: r(ixi^s)
1332
1333 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1334 res(ixo^s)=w(ixo^s,p_)/(r(ixo^s)*w(ixo^s,rho_))
1335
1336 end subroutine hd_get_temperature_from_prim
1337
1338 ! Calculate flux f_idim[iw]
1339 subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1341 use mod_dust, only: dust_get_flux_prim
1342
1343 integer, intent(in) :: ixi^l, ixo^l, idim
1344 ! conservative w
1345 double precision, intent(in) :: wc(ixi^s, 1:nw)
1346 ! primitive w
1347 double precision, intent(in) :: w(ixi^s, 1:nw)
1348 double precision, intent(in) :: x(ixi^s, 1:ndim)
1349 double precision, intent(out) :: f(ixi^s, nwflux)
1350
1351 double precision :: pth(ixi^s)
1352 integer :: ix^db
1353
1354 if (hd_energy) then
1355 {do ix^db=ixomin^db,ixomax^db\}
1356 f(ix^d,rho_)=w(ix^d,mom(idim))*w(ix^d,rho_)
1357 ! Momentum flux is v_i*m_i, +p in direction idim
1358 ^c&f(ix^d,m^c_)=w(ix^d,mom(idim))*wc(ix^d,m^c_)\
1359 f(ix^d,mom(idim))=f(ix^d,mom(idim))+w(ix^d,p_)
1360 ! Energy flux is v_i*(e + p)
1361 f(ix^d,e_)=w(ix^d,mom(idim))*(wc(ix^d,e_)+w(ix^d,p_))
1362 {end do\}
1363 else
1364 call hd_get_pthermal(wc, x, ixi^l, ixo^l, pth)
1365 {do ix^db=ixomin^db,ixomax^db\}
1366 f(ix^d,rho_)=w(ix^d,mom(idim))*w(ix^d,rho_)
1367 ! Momentum flux is v_i*m_i, +p in direction idim
1368 ^c&f(ix^d,m^c_)=w(ix^d,mom(idim))*wc(ix^d,m^c_)\
1369 f(ix^d,mom(idim))=f(ix^d,mom(idim))+pth(ix^d)
1370 {end do\}
1371 end if
1372
1373 if(hd_radiation_fld)then
1374 {do ix^db=ixomin^db,ixomax^db\}
1375 ! advection of radiation enery v_i*r_e
1376 f(ix^d,r_e)=w(ix^d,mom(idim))*wc(ix^d,r_e)
1377 {end do\}
1378 endif
1379
1380 do ix1 = 1, hd_n_tracer
1381 f(ixo^s, tracer(ix1)) = w(ixo^s,mom(idim)) * w(ixo^s, tracer(ix1))
1382 end do
1383
1384 ! Dust fluxes
1385 if (hd_dust) then
1386 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1387 end if
1388
1389 end subroutine hd_get_flux
1390
1391 !> Add geometrical source terms to w
1392 !>
1393 !> Notice that the expressions of the geometrical terms depend only on ndir,
1394 !> not ndim. Eg, they are the same in 2.5D and in 3D, for any geometry.
1395 !>
1396 subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1401 use mod_geometry
1402 integer, intent(in) :: ixi^l, ixo^l
1403 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:ndim)
1404 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
1405 double precision :: pth(ixi^s), source(ixi^s), minrho
1406 integer :: iw,idir, h1x^l{^nooned, h2x^l}
1407 integer :: mr_,mphi_ ! Polar var. names
1408 integer :: irho, ifluid, n_fluids
1409 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1410
1411 if (hd_dust) then
1412 n_fluids = 1 + dust_n_species
1413 else
1414 n_fluids = 1
1415 end if
1416
1417 select case (coordinate)
1418
1420 !the user provides the functions of exp_factor and del_exp_factor
1421 if(associated(usr_set_surface)) call usr_set_surface(ixi^l,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
1422 if(hd_energy) then
1423 source(ixo^s)=wprim(ixo^s, p_)
1424 else
1425 if(.not. associated(usr_set_pthermal)) then
1426 source(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1427 else
1428 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1429 end if
1430 end if
1431 source(ixo^s) = source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1432 w(ixo^s,mom(1)) = w(ixo^s,mom(1)) + qdt*source(ixo^s)
1433
1434 case (cylindrical)
1435 do ifluid = 0, n_fluids-1
1436 ! s[mr]=(pthermal+mphi**2/rho)/radius
1437 if (ifluid == 0) then
1438 ! gas
1439 irho = rho_
1440 mr_ = mom(r_)
1441 if(phi_>0) mphi_ = mom(phi_)
1442 if(hd_energy) then
1443 source(ixo^s)=wprim(ixo^s, p_)
1444 else
1445 if(.not. associated(usr_set_pthermal)) then
1446 source(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1447 else
1448 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1449 end if
1450 end if
1451 minrho = 0.0d0
1452 else
1453 ! dust : no pressure
1454 irho = dust_rho(ifluid)
1455 mr_ = dust_mom(r_, ifluid)
1456 if(phi_>0) mphi_ = dust_mom(phi_, ifluid)
1457 source(ixi^s) = zero
1458 minrho = 0.0d0
1459 end if
1460 if(phi_ > 0) then
1461 where (wct(ixo^s, irho) > minrho)
1462 source(ixo^s) = source(ixo^s) + wct(ixo^s,mphi_)*wprim(ixo^s,mphi_)
1463 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt*source(ixo^s)/x(ixo^s,r_)
1464 end where
1465 ! s[mphi]=(-mphi*vr)/radius
1466 where (wct(ixo^s, irho) > minrho)
1467 source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s, mr_)
1468 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * source(ixo^s) / x(ixo^s, r_)
1469 end where
1470 else
1471 ! s[mr]=2pthermal/radius
1472 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1473 end if
1474 end do
1475 case (spherical)
1476 if (hd_dust) then
1477 call mpistop("Dust geom source terms not implemented yet with spherical geometries")
1478 end if
1479 mr_ = mom(r_)
1480 if(phi_>0) mphi_ = mom(phi_)
1481 h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
1482 if(hd_energy) then
1483 pth(ixo^s)=wprim(ixo^s, p_)
1484 else
1485 if(.not. associated(usr_set_pthermal)) then
1486 pth(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1487 else
1488 call usr_set_pthermal(wct,x,ixi^l,ixo^l,pth)
1489 end if
1490 end if
1491 ! s[mr]=((vtheta**2+vphi**2)*rho+2*p)/r
1492 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1493 *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
1494 /block%dvolume(ixo^s)
1495 do idir = 2, ndir
1496 source(ixo^s) = source(ixo^s) + wprim(ixo^s, mom(idir))**2 * wprim(ixo^s, rho_)
1497 end do
1498 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, 1)
1499
1500 {^nooned
1501 ! s[mtheta]=-(vr*vtheta*rho)/r+cot(theta)*(vphi**2*rho+p)/r
1502 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1503 * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
1504 / block%dvolume(ixo^s)
1505 if (ndir == 3) then
1506 source(ixo^s) = source(ixo^s) + (wprim(ixo^s, mom(3))**2 * wprim(ixo^s, rho_)) / tan(x(ixo^s, 2))
1507 end if
1508 source(ixo^s) = source(ixo^s) - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)
1509 w(ixo^s, mom(2)) = w(ixo^s, mom(2)) + qdt * source(ixo^s) / x(ixo^s, 1)
1510
1511 if (ndir == 3) then
1512 ! s[mphi]=-(vphi*vr/rho)/r-cot(theta)*(vtheta*vphi/rho)/r
1513 source(ixo^s) = -(wprim(ixo^s, mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)&
1514 - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mom(3))) * wprim(ixo^s, rho_) / tan(x(ixo^s, 2))
1515 w(ixo^s, mom(3)) = w(ixo^s, mom(3)) + qdt * source(ixo^s) / x(ixo^s, 1)
1516 end if
1517 }
1518 end select
1519
1520 if (hd_rotating_frame) then
1521 if (hd_dust) then
1522 call mpistop("Rotating frame not implemented yet with dust")
1523 else
1524 call rotating_frame_add_source(qdt,dtfactor,ixi^l,ixo^l,wprim,w,x)
1525 end if
1526 end if
1527
1528 end subroutine hd_add_source_geom
1529
1530 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
1531 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1536 use mod_usr_methods, only: usr_gravity
1538 use mod_cak_force, only: cak_add_source
1539
1540 integer, intent(in) :: ixi^l, ixo^l
1541 double precision, intent(in) :: qdt, dtfactor
1542 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:ndim)
1543 double precision, intent(inout) :: w(ixi^s, 1:nw)
1544 logical, intent(in) :: qsourcesplit
1545 logical, intent(inout) :: active
1546
1547 double precision :: gravity_field(ixi^s, 1:ndim)
1548 integer :: idust, idim
1549
1550 if(hd_dust .and. .not. hd_dust_implicit) then
1551 call dust_add_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1552 end if
1553
1554 if(hd_radiative_cooling) then
1555 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1556 qsourcesplit,active, rc_fl)
1557 end if
1558
1559 if(hd_viscosity) then
1560 call viscosity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1561 hd_energy,qsourcesplit,active)
1562 end if
1563
1564 if (hd_gravity) then
1565 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1566 hd_energy,qsourcesplit,active)
1567
1568 if (hd_dust .and. qsourcesplit .eqv. grav_split) then
1569 active = .true.
1570
1571 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1572 do idust = 1, dust_n_species
1573 do idim = 1, ndim
1574 w(ixo^s, dust_mom(idim, idust)) = w(ixo^s, dust_mom(idim, idust)) &
1575 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s, dust_rho(idust))
1576 end do
1577 end do
1578 end if
1579 end if
1580
1581 if (hd_cak_force) then
1582 call cak_add_source(qdt,ixi^l,ixo^l,wct,w,x,hd_energy,qsourcesplit,active)
1583 end if
1584
1585 ! This is where the radiation force and heating/cooling are added
1586 if (hd_radiation_fld) then
1587 call hd_add_radiation_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,qsourcesplit,active)
1588 endif
1589
1590 if(hd_partial_ionization) then
1591 if(.not.qsourcesplit) then
1592 active = .true.
1593 call hd_update_temperature(ixi^l,ixo^l,wct,w,x)
1594 end if
1595 end if
1596
1597 end subroutine hd_add_source
1598
1599 subroutine hd_add_radiation_source(qdt,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1600 use mod_constants
1602 use mod_usr_methods
1603 use mod_fld
1604
1605 integer, intent(in) :: ixi^l, ixo^l
1606 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
1607 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
1608 double precision, intent(inout) :: w(ixi^s,1:nw)
1609 logical, intent(in) :: qsourcesplit
1610 logical, intent(inout) :: active
1611
1612 ! add radiation force and work done by it, changes momentum and gas energy
1613 ! handle photon tiring, heating and cooling exchange between gas and radiation field
1614 call add_fld_rad_force(qdt,ixi^l,ixo^l,wct,wctprim,w,x,qsourcesplit,active)
1615 if(fix_small_values) call hd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'fld_add_radiation')
1616
1617 end subroutine hd_add_radiation_source
1618
1619 subroutine hd_get_dt(wprim, ixI^L, ixO^L, dtnew, dx^D, x)
1621 use mod_dust, only: dust_get_dt
1623 use mod_gravity, only: gravity_get_dt
1624 use mod_cak_force, only: cak_get_dt
1625 use mod_fld, only: fld_radforce_get_dt
1626
1627 integer, intent(in) :: ixi^l, ixo^l
1628 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
1629 double precision, intent(in) :: wprim(ixi^s, 1:nw)
1630 double precision, intent(inout) :: dtnew
1631
1632 dtnew = bigdouble
1633
1634 if(hd_dust) then
1635 call dust_get_dt(wprim, ixi^l, ixo^l, dtnew, dx^d, x)
1636 end if
1637
1638 if(hd_viscosity) then
1639 call viscosity_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1640 end if
1641
1642 if(hd_gravity) then
1643 call gravity_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1644 end if
1645
1646 if (hd_cak_force) then
1647 call cak_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1648 end if
1649
1650 if(hd_radiation_fld) then
1651 call fld_radforce_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1652 endif
1653
1654 end subroutine hd_get_dt
1655
1656 function hd_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
1657 use mod_global_parameters, only: nw, ndim
1658 integer, intent(in) :: ixi^l, ixo^l
1659 double precision, intent(in) :: w(ixi^s, nw)
1660 double precision :: ke(ixo^s)
1661 double precision, intent(in), optional :: inv_rho(ixo^s)
1662
1663 if (present(inv_rho)) then
1664 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * inv_rho
1665 else
1666 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) / w(ixo^s, rho_)
1667 end if
1668 end function hd_kin_en
1669
1670 function hd_inv_rho(w, ixI^L, ixO^L) result(inv_rho)
1671 use mod_global_parameters, only: nw, ndim
1672 integer, intent(in) :: ixi^l, ixo^l
1673 double precision, intent(in) :: w(ixi^s, nw)
1674 double precision :: inv_rho(ixo^s)
1675
1676 ! Can make this more robust
1677 inv_rho = 1.0d0 / w(ixo^s, rho_)
1678 end function hd_inv_rho
1679
1680 subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1681 ! handles hydro (density,pressure,velocity) bootstrapping
1682 ! any negative dust density is flagged as well (and throws an error)
1683 ! small_values_method=replace also for dust
1687 logical, intent(in) :: primitive
1688 integer, intent(in) :: ixi^l,ixo^l
1689 double precision, intent(inout) :: w(ixi^s,1:nw)
1690 double precision, intent(in) :: x(ixi^s,1:ndim)
1691 character(len=*), intent(in) :: subname
1692
1693 integer :: n,idir
1694 logical :: flag(ixi^s,1:nw)
1695
1696 call hd_check_w(primitive, ixi^l, ixo^l, w, flag)
1697
1698 if (any(flag)) then
1699 select case (small_values_method)
1700 case ("replace")
1701 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
1702 do idir = 1, ndir
1703 if(small_values_fix_iw(mom(idir))) then
1704 where(flag(ixo^s,rho_)) w(ixo^s, mom(idir)) = 0.0d0
1705 end if
1706 end do
1707 if(hd_radiation_fld)then
1708 if (small_values_fix_iw(r_e)) then
1709 where(flag(ixo^s,r_e)) w(ixo^s,r_e) = small_r_e
1710 end if
1711 end if
1712 if(hd_energy)then
1713 if(small_values_fix_iw(e_)) then
1714 if(primitive) then
1715 where(flag(ixo^s,rho_)) w(ixo^s, p_) = small_pressure
1716 else
1717 where(flag(ixo^s,rho_)) w(ixo^s, e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1718 endif
1719 end if
1720 endif
1721
1722 if(hd_energy) then
1723 if(primitive) then
1724 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
1725 else
1726 where(flag(ixo^s,e_))
1727 ! Add kinetic energy
1728 w(ixo^s,e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1729 end where
1730 end if
1731 end if
1732
1733 if(hd_dust)then
1734 do n=1,dust_n_species
1735 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1736 do idir = 1, ndir
1737 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1738 enddo
1739 enddo
1740 endif
1741 case ("average")
1742 if(primitive)then
1743 ! averaging for all primitive fields, including dust
1744 call small_values_average(ixi^l, ixo^l, w, x, flag)
1745 else
1746 ! do averaging of density
1747 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
1748 if(hd_energy) then
1749 ! do averaging of pressure
1750 w(ixi^s,p_)=(hd_gamma-1.d0)*(w(ixi^s,e_) &
1751 -0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_))
1752 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
1753 w(ixi^s,e_)=w(ixi^s,p_)/(hd_gamma-1.d0) &
1754 +0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_)
1755 end if
1756 if(hd_radiation_fld) then
1757 ! do averaging of radiative energy density
1758 call small_values_average(ixi^l, ixo^l, w, x, flag, r_e)
1759 endif
1760 if(hd_dust)then
1761 do n=1,dust_n_species
1762 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1763 do idir = 1, ndir
1764 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1765 enddo
1766 enddo
1767 endif
1768 endif
1769 case default
1770 if(.not.primitive) then
1771 !convert w to primitive
1772 ! Calculate pressure = (gamma-1) * (e-ek)
1773 if(hd_energy) then
1774 w(ixo^s,p_)=(hd_gamma-1.d0)*(w(ixo^s,e_)-hd_kin_en(w,ixi^l,ixo^l))
1775 end if
1776 ! Convert gas momentum to velocity
1777 do idir = 1, ndir
1778 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir))/w(ixo^s,rho_)
1779 end do
1780 end if
1781 ! NOTE: dust entries may still have conserved values here
1782 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1783 end select
1784 end if
1785 end subroutine hd_handle_small_values
1786
1787 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1790 integer, intent(in) :: ixi^l, ixo^l
1791 double precision, intent(in) :: w(ixi^s,1:nw)
1792 double precision, intent(in) :: x(ixi^s,1:ndim)
1793 double precision, intent(out):: rfactor(ixi^s)
1794
1795 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1796
1797 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1798 ! assume the first and second ionization of Helium have the same degree
1799 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
1800
1801 end subroutine rfactor_from_temperature_ionization
1802
1803 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1805 integer, intent(in) :: ixi^l, ixo^l
1806 double precision, intent(in) :: w(ixi^s,1:nw)
1807 double precision, intent(in) :: x(ixi^s,1:ndim)
1808 double precision, intent(out):: rfactor(ixi^s)
1809
1810 rfactor(ixo^s)=rr
1811
1812 end subroutine rfactor_from_constant_ionization
1813
1814 subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1817
1818 integer, intent(in) :: ixi^l, ixo^l
1819 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1820 double precision, intent(inout) :: w(ixi^s,1:nw)
1821
1822 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1823
1824 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1825
1826 call hd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1827
1828 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1829 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1830
1831 end subroutine hd_update_temperature
1832
1833end 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...
subroutine cak_init(phys_gamma)
Initialize the module.
subroutine cak_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
Check time step for total radiation contribution.
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.
Module for physical and numeric constants.
double precision, parameter bigdouble
A very large real number.
double precision, parameter const_rad_a
Module for including dust species, which interact with the gas through a drag force.
Definition mod_dust.t:3
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:532
subroutine, public dust_evaluate_implicit(qtc, psa)
inplace update of psa==>F_im(psa)
Definition mod_dust.t:585
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:47
subroutine, public dust_to_conserved(ixil, ixol, w, x)
Definition mod_dust.t:208
integer, public, protected dust_n_species
The number of dust species.
Definition mod_dust.t:37
subroutine, public dust_get_flux_prim(w, x, ixil, ixol, idim, f)
Definition mod_dust.t:275
subroutine, public dust_check_w(ixil, ixol, w, flag)
Definition mod_dust.t:194
integer, dimension(:), allocatable, public, protected dust_rho
Indices of the dust densities.
Definition mod_dust.t:44
subroutine, public dust_get_cmax(w, x, ixil, ixol, idim, cmax, cmin)
Definition mod_dust.t:1018
subroutine, public dust_check_params()
Definition mod_dust.t:154
subroutine, public dust_get_cmax_prim(w, x, ixil, ixol, idim, cmax, cmin)
Definition mod_dust.t:1042
subroutine, public dust_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
Get dt related to dust and gas stopping time (Laibe 2011)
Definition mod_dust.t:900
subroutine, public dust_init(g_rho, g_mom, g_energy)
Definition mod_dust.t:95
subroutine, public dust_implicit_update(dtfactor, qdt, qtc, psb, psa)
Implicit solve of psb=psa+dtfactor*dt*F_im(psb)
Definition mod_dust.t:654
Module for flux conservation near refinement boundaries.
Module for flux limited diffusion (FLD)-approximation in Radiation-(Magneto)hydrodynamics simulations...
Definition mod_fld.t:13
logical fld_no_mg
Definition mod_fld.t:27
subroutine, public fld_get_radpress(w, x, ixil, ixol, rad_pressure)
Returns Radiation Pressure as tensor NOTE: w is primitive on entry.
Definition mod_fld.t:412
double precision, public fld_bisect_tol
Tolerance for bisection method for Energy sourceterms This is a percentage of the minimum of gas- and...
Definition mod_fld.t:23
double precision, public fld_diff_tol
Tolerance for radiative Energy diffusion.
Definition mod_fld.t:25
character(len=8) fld_opacity_law
switches for opacity
Definition mod_fld.t:29
character(len=16) fld_fluxlimiter
flux limiter choice
Definition mod_fld.t:32
character(len=50) fld_opal_table
Definition mod_fld.t:30
double precision, public fld_kappa0
Opacity per unit of unit_density.
Definition mod_fld.t:20
subroutine, public fld_init(he_abundance, r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
Definition mod_fld.t:78
subroutine, public add_fld_rad_force(qdt, ixil, ixol, wct, wctprim, w, x, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO This subroutine handles th...
Definition mod_fld.t:177
character(len=8) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
Definition mod_fld.t:38
logical fld_radforce_split
source split for energy interact and radforce:
Definition mod_fld.t:18
subroutine, public fld_radforce_get_dt(w, ixil, ixol, dtnew, dxd, x)
get dt limit for radiation force and FLD explicit source additions NOTE: w is primitive on entry
Definition mod_fld.t:304
integer nth_for_diff_mg
diffusion coefficient stencil control
Definition mod_fld.t:36
Module with geometry-related routines (e.g., divergence, curl)
Definition mod_geometry.t:2
integer coordinate
Definition mod_geometry.t:7
integer, parameter spherical
integer, parameter cylindrical
integer, parameter cartesian_expansion
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 unit_time
Physical scaling factor for time.
double precision unit_density
Physical scaling factor for density.
double precision unit_opacity
Physical scaling factor for Opacity.
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.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
logical slab
Cartesian geometry or not.
integer nwauxio
Number of auxiliary variables that are only included in output.
double precision unit_velocity
Physical scaling factor for velocity.
double precision unit_temperature
Physical scaling factor for temperature.
double precision unit_radflux
Physical scaling factor for radiation flux.
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
integer nghostcells
Number of ghost cells surrounding a grid.
logical phys_trac
Use TRAC for MHD or 1D HD.
logical fix_small_values
fix small values with average or replace methods
logical crash
Save a snapshot before crash a run met unphysical values.
logical use_multigrid
Use multigrid (only available in 2D and 3D)
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
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(wprim, 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
Hydrodynamics physics module.
Definition mod_hd_phys.t:2
integer, public, protected m
Definition mod_hd_phys.t:60
subroutine, public hd_check_params
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, public hd_get_pthermal_plus_pradiation(w, x, ixil, ixol, pth_plus_prad)
calculates the sum of the gas pressure and max Prad tensor element NOTE: only for diagnostic purposes...
subroutine, public hd_ei_to_e(ixil, ixol, w, x)
Transform internal energy to total energy.
subroutine, public hd_get_temperature_from_prim(w, x, ixil, ixol, res)
Calculate temperature=p/rho when in we work in primitives (hence e_ is p_)
integer, public, protected e_
Index of the energy density (-1 if not present)
Definition mod_hd_phys.t:66
logical, public, protected hd_radiative_cooling
Whether radiative cooling is added.
Definition mod_hd_phys.t:20
double precision, public, protected rr
double precision, public hd_gamma
The adiabatic index.
Definition mod_hd_phys.t:81
integer, public, protected hd_trac_type
Definition mod_hd_phys.t:91
logical, public, protected hd_particles
Whether particles module is added.
Definition mod_hd_phys.t:39
logical, public, protected hd_radiation_fld
Whether radiation-gas interaction is handled using flux limited diffusion.
Definition mod_hd_phys.t:30
type(tc_fluid), allocatable, public tc_fl
Definition mod_hd_phys.t:16
subroutine, public hd_check_w(primitive, ixil, ixol, w, flag)
Returns logical argument flag where values are ok.
logical, public, protected hd_viscosity
Whether viscosity is added.
Definition mod_hd_phys.t:33
integer, public, protected r_e
Index of the radiation energy (when fld active)
Definition mod_hd_phys.t:72
procedure(sub_get_pthermal), pointer, public hd_get_rfactor
integer, public, protected c
Indices of the momentum density for the form of better vectorization.
Definition mod_hd_phys.t:60
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.
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
Definition mod_hd_phys.t:78
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
integer, public, protected te_
Indices of temperature.
Definition mod_hd_phys.t:75
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
Definition mod_hd_phys.t:57
subroutine, public hd_get_pradiation_from_prim(w, x, ixil, ixol, prad)
Calculate radiation pressure within ixO^L NOTE: w is primitive on entry here! NOTE: used in FLD modul...
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
Definition mod_hd_phys.t:97
double precision function, dimension(ixo^s), public hd_kin_en(w, ixil, ixol, inv_rho)
subroutine, public hd_to_conserved(ixil, ixol, w, x)
Transform primitive variables into conservative ones.
logical, public, protected hd_cak_force
Whether CAK radiation line force is activated.
Definition mod_hd_phys.t:45
subroutine, public hd_phys_init()
Initialize the module.
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
Definition mod_hd_phys.t:63
subroutine, public hd_get_csrad2(w, x, ixil, ixol, csound)
Calculate modified squared sound speed for FLD NOTE: only for diagnostic purposes,...
logical, public, protected hd_thermal_conduction
Whether thermal conduction is added.
Definition mod_hd_phys.t:15
logical, public, protected eq_state_units
double precision, public hd_adiab
The adiabatic constant.
Definition mod_hd_phys.t:87
subroutine, public hd_get_trad(w, x, ixil, ixol, trad)
Calculates radiation temperature Note use of cgs units here in factor const_rad_a.
subroutine, public hd_to_primitive(ixil, ixol, w, x)
Transform conservative variables into primitive ones.
subroutine, public hd_get_csrad2_prim(w, x, ixil, ixol, csound)
Calculate modified squared sound speed for FLD NOTE: w is primitive on entry here!...
integer, public, protected rho_
Index of the density (in the w array)
Definition mod_hd_phys.t:54
logical, public, protected hd_partial_ionization
Whether plasma is partially ionized.
Definition mod_hd_phys.t:51
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
Definition mod_hd_phys.t:94
logical, public, protected hd_gravity
Whether gravity is added.
Definition mod_hd_phys.t:36
integer, public, protected c_
Definition mod_hd_phys.t:60
type(rc_fluid), allocatable, public rc_fl
Definition mod_hd_phys.t:21
logical, public, protected hd_dust_implicit
Whether dust is added using and implicit update in IMEX.
Definition mod_hd_phys.t:27
logical, public, protected hd_trac
Whether TRAC method is used.
Definition mod_hd_phys.t:90
integer, public, protected hd_n_tracer
Number of tracer species.
Definition mod_hd_phys.t:48
type(te_fluid), allocatable, public te_fl_hd
Definition mod_hd_phys.t:17
logical, public, protected hd_rotating_frame
Whether rotating frame is activated.
Definition mod_hd_phys.t:42
subroutine, public hd_get_pthermal(w, x, ixil, ixol, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
Definition mod_hd_phys.t:69
subroutine, public hd_e_to_ei(ixil, ixol, w, x)
Transform total energy to internal energy.
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.
subroutine particles_init()
Initialize particle data and parameters.
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 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.
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_hd(w, ixil, ixol, dxd, x, fl)
Get the explicit timestep for the TC (hd implementation) Note: also used in 1D MHD (or for neutrals i...
subroutine tc_init_params(phys_gamma)
subroutine, public sts_set_source_tc_hd(ixil, ixol, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux, fl)
subroutine 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
integer nw
Total number of variables.
integer number_species
number of species: each species has different characterictic speeds and should be used accordingly in...
The module add viscous source terms and check time step.
subroutine viscosity_init(phys_wider_stencil)
Initialize the module.
subroutine viscosity_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
procedure(sub_add_source), pointer, public viscosity_add_source