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 FIP passive scalar rho*fip in conserved form, fip in primitive form
78 integer, public, protected :: fip_ = -1
79
80 !> Whether FIP passive scalar is enabled
81 logical, public, protected :: hd_fip = .false.
82
83 !> Index of the cutoff temperature for the TRAC method
84 integer, public, protected :: tcoff_
85
86 !> The adiabatic index
87 double precision, public :: hd_gamma = 5.d0/3.0d0
88
89 !> gamma minus one and its inverse
90 double precision :: gamma_1, inv_gamma_1
91
92 !> The adiabatic constant
93 double precision, public :: hd_adiab = 1.0d0
94
95 !> Whether TRAC method is used
96 logical, public, protected :: hd_trac = .false.
97 integer, public, protected :: hd_trac_type = 1
98
99 !> Helium abundance over Hydrogen
100 double precision, public, protected :: he_abundance=0.1d0
101 !> Ionization fraction of H
102 !> H_ion_fr = H+/(H+ + H)
103 double precision, public, protected :: h_ion_fr=1d0
104 !> Ionization fraction of He
105 !> He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
106 double precision, public, protected :: he_ion_fr=1d0
107 !> Ratio of number He2+ / number He+ + He2+
108 !> He_ion_fr2 = He2+/(He2+ + He+)
109 double precision, public, protected :: he_ion_fr2=1d0
110 ! used for eq of state when it is not defined by units,
111 ! the units do not contain terms related to ionization fraction
112 ! and it is p = RR * rho * T
113 double precision, public, protected :: rr=1d0
114 ! remove the below flag and assume default value = .false.
115 ! when eq state properly implemented everywhere
116 ! and not anymore through units
117 logical, public, protected :: eq_state_units = .true.
118
119 procedure(sub_get_pthermal), pointer :: hd_get_rfactor => null()
120 ! Public methods
121 public :: hd_phys_init
122 public :: hd_kin_en
123 public :: hd_get_pthermal
124 public :: hd_get_csound2
125 public :: hd_to_conserved
126 public :: hd_to_primitive
127 public :: hd_check_params
128 public :: hd_check_w
129 public :: hd_e_to_ei
130 public :: hd_ei_to_e
131 ! in FLD used as phys_get_Rfactor
132 public :: hd_get_rfactor
133 ! Begin: following relevant for radiative hydro using FLD
134 ! first four are local and only of interest for mod_usr applications
135 ! where they can be used in diagnostics
136 ! NOTE those with _prim expect primitives on entry
138 public :: hd_get_csrad2
139 public :: hd_get_trad
141 ! the following used in FLD modules
142 ! as pointer phys_get_csrad2
143 public :: hd_get_csrad2_prim
144 ! the following used in FLD modules
145 ! as pointer phys_get_tgas
147 ! End: following relevant for radiative hydro using FLD
149
150contains
151
152 !> Read this module's parameters from a file
153 subroutine hd_read_params(files)
155 character(len=*), intent(in) :: files(:)
156 integer :: n
157
158 namelist /hd_list/ hd_energy, hd_n_tracer, hd_gamma, hd_adiab, &
164
165 do n = 1, size(files)
166 open(unitpar, file=trim(files(n)), status="old")
167 read(unitpar, hd_list, end=111)
168111 close(unitpar)
169 end do
170
171 end subroutine hd_read_params
172
173 !> Write this module's parameters to a snapsoht
174 subroutine hd_write_info(fh)
176 integer, intent(in) :: fh
177 integer, parameter :: n_par = 1
178 double precision :: values(n_par)
179 character(len=name_len) :: names(n_par)
180 integer, dimension(MPI_STATUS_SIZE) :: st
181 integer :: er
182
183 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
184
185 names(1) = "gamma"
186 values(1) = hd_gamma
187 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
188 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
189 end subroutine hd_write_info
190
191 !> Initialize the module
192 subroutine hd_phys_init()
196 use mod_dust, only: dust_init
198 use mod_gravity, only: gravity_init
200 use mod_cak_force, only: cak_init
205 use mod_fld
206
207 integer :: itr, idir
208
209 call hd_read_params(par_files)
210
211 physics_type = "hd"
212 phys_energy = hd_energy
213 phys_total_energy = hd_energy
214 phys_internal_e = .false.
215 phys_gamma = hd_gamma
216 phys_partial_ionization=hd_partial_ionization
217
219 if(phys_trac) then
220 if(ndim .eq. 1) then
221 if(hd_trac_type .gt. 2) then
223 if(mype==0) write(*,*) 'WARNING: set hd_trac_type=1'
224 end if
226 else
227 phys_trac=.false.
228 if(mype==0) write(*,*) 'WARNING: set hd_trac=F when ndim>=2'
229 end if
230 end if
231
232 ! set default gamma for polytropic/isothermal process
233 if(.not.hd_energy) then
234 if(hd_thermal_conduction) then
236 if(mype==0) write(*,*) 'WARNING: set hd_thermal_conduction=F when hd_energy=F'
237 end if
238 if(hd_radiative_cooling) then
240 if(mype==0) write(*,*) 'WARNING: set hd_radiative_cooling=F when hd_energy=F'
241 end if
242 end if
243 if(.not.eq_state_units) then
244 if(hd_partial_ionization) then
246 if(mype==0) write(*,*) 'WARNING: set hd_partial_ionization=F when eq_state_units=F'
247 end if
248 end if
250
251 allocate(start_indices(number_species),stop_indices(number_species))
252
253 ! set the index of the first flux variable for species 1
254 start_indices(1)=1
255 ! Determine flux variables
256 rho_ = var_set_rho()
257
258 allocate(mom(ndir))
259 mom(:) = var_set_momentum(ndir)
260 m^c_=mom(^c);
261
262 ! Set index of energy variable
263 if (hd_energy) then
264 e_ = var_set_energy()
265 p_ = e_
266 else
267 e_ = -1
268 p_ = -1
269 end if
270
271 if(hd_radiation_fld)then
272 if(hd_cak_force)then
273 if(mype==0) then
274 write(*,*)'Warning: CAK force addition together with FLD radiation'
275 endif
276 endif
278 if(mype==0) then
279 write(*,*)'Warning: Optically thin cooling together with FLD radiation'
280 endif
281 endif
282 if(hd_dust.and.hd_dust_implicit)then
283 call mpistop('implicit dust addition not compatible with FLD radiation')
284 endif
285 if(.not.hd_energy)then
286 call mpistop('using FLD implies the use of an energy equation, set hd_energy=T')
287 else
288 !> set added variable and equation for radiation energy
289 r_e = var_set_radiation_energy()
290 phys_get_tgas => hd_get_temperature_from_prim
291 phys_get_csrad2 => hd_get_csrad2_prim
292 !> Initiate radiation-closure module
293 call fld_init(hd_gamma)
294 endif
295 else
296 r_e=-1
297 endif
298
299 phys_get_dt => hd_get_dt
300 phys_get_cmax => hd_get_cmax
301 phys_get_tcutoff => hd_get_tcutoff
302 phys_get_cbounds => hd_get_cbounds
303 phys_get_flux => hd_get_flux
304 phys_add_source_geom => hd_add_source_geom
305 phys_add_source => hd_add_source
306 phys_to_conserved => hd_to_conserved
307 phys_to_primitive => hd_to_primitive
308 phys_check_params => hd_check_params
309 phys_check_w => hd_check_w
310 phys_get_pthermal => hd_get_pthermal
311 phys_get_v => hd_get_v
312 phys_get_rho => hd_get_rho
313 phys_write_info => hd_write_info
314 phys_handle_small_values => hd_handle_small_values
315
316 ! derive units from basic units
317 call hd_physical_units()
318
319 if (hd_dust) then
320 call dust_init(rho_, mom(:), e_)
321 endif
322
323 if (hd_fip) then
324 fip_ = var_set_fluxvar('rho_fip', 'fip', need_bc=.false.)
325 else
326 fip_ = -1
327 end if
328
329 allocate(tracer(hd_n_tracer))
330
331 ! Set starting index of tracers
332 do itr = 1, hd_n_tracer
333 tracer(itr) = var_set_fluxvar("trc", "trp", itr, need_bc=.false.)
334 end do
335
336 ! set temperature as an auxiliary variable to get ionization degree
337 if(hd_partial_ionization) then
338 te_ = var_set_auxvar('Te','Te')
339 else
340 te_ = -1
341 end if
342
343 ! set number of variables which need update ghostcells
344 nwgc=nwflux+nwaux
345
346 ! set the index of the last flux variable for species 1
347 stop_indices(1)=nwflux
348
349 if(hd_trac) then
350 tcoff_ = var_set_wextra()
351 iw_tcoff=tcoff_
352 else
353 tcoff_ = -1
354 end if
355
356 ! choose Rfactor in ideal gas law
357 if(hd_partial_ionization) then
358 hd_get_rfactor=>rfactor_from_temperature_ionization
359 phys_update_temperature => hd_update_temperature
360 else if(associated(usr_rfactor)) then
362 else
363 hd_get_rfactor=>rfactor_from_constant_ionization
364 end if
365
366 phys_get_rfactor => hd_get_rfactor
367
368 ! initialize thermal conduction module
369 if (hd_thermal_conduction) then
370 if (.not. hd_energy) &
371 call mpistop("thermal conduction needs hd_energy=T")
372
373 call sts_init()
375
376 allocate(tc_fl)
377 call tc_get_hd_params(tc_fl,tc_params_read_hd)
378 call add_sts_method(hd_get_tc_dt_hd,hd_sts_set_source_tc_hd,e_,1,e_,1,.false.)
380 call set_error_handling_to_head(hd_tc_handle_small_e)
381 tc_fl%get_temperature_from_conserved => hd_get_temperature_from_etot
382 tc_fl%get_temperature_from_eint => hd_get_temperature_from_eint
383 tc_fl%get_rho => hd_get_rho
384 tc_fl%e_ = e_
385 end if
386
387 ! Initialize radiative cooling module
388 if (hd_radiative_cooling) then
389 if (.not. hd_energy) &
390 call mpistop("radiative cooling needs hd_energy=T")
392 allocate(rc_fl)
393 rc_fl%fip_ = fip_
394 call radiative_cooling_init(rc_fl,rc_params_read)
395 rc_fl%get_rho => hd_get_rho
396 rc_fl%get_pthermal => hd_get_pthermal
397 rc_fl%get_var_Rfactor => hd_get_rfactor
398 rc_fl%e_ = e_
399 rc_fl%Tcoff_ = tcoff_
400 end if
401 allocate(te_fl_hd)
402 te_fl_hd%get_rho=> hd_get_rho
403 te_fl_hd%get_pthermal=> hd_get_pthermal
404 te_fl_hd%get_var_Rfactor => hd_get_rfactor
405{^ifthreed
406 phys_te_images => hd_te_images
407}
408 ! Initialize viscosity module
409 if (hd_viscosity) call viscosity_init(phys_wider_stencil)
410
411 ! Initialize gravity module
412 if (hd_gravity) call gravity_init()
413
414 ! Initialize rotating_frame module
416
417 ! Initialize CAK radiation force module
419
420
421 ! Check whether custom flux types have been defined
422 if (.not. allocated(flux_type)) then
423 allocate(flux_type(ndir, nw))
424 flux_type = flux_default
425 else if (any(shape(flux_type) /= [ndir, nw])) then
426 call mpistop("phys_check error: flux_type has wrong shape")
427 end if
428
429 nvector = 1 ! No. vector vars
430 allocate(iw_vector(nvector))
431 iw_vector(1) = mom(1) - 1
432 ! initialize ionization degree table
434
435 end subroutine hd_phys_init
436
437{^ifthreed
438 subroutine hd_te_images
441 select case(convert_type)
442 case('EIvtiCCmpi','EIvtuCCmpi')
444 case('ESvtiCCmpi','ESvtuCCmpi')
446 case('SIvtiCCmpi','SIvtuCCmpi')
448 case('WIvtiCCmpi','WIvtuCCmpi')
450 case default
451 call mpistop("Error in synthesize emission: Unknown convert_type")
452 end select
453 end subroutine hd_te_images
454}
455!!start th cond
456 ! wrappers for STS functions in thermal_conductivity module
457 ! which take as argument the tc_fluid (defined in the physics module)
458 subroutine hd_sts_set_source_tc_hd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
462 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
463 double precision, intent(in) :: x(ixi^s,1:ndim)
464 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
465 double precision, intent(in) :: my_dt
466 logical, intent(in) :: fix_conserve_at_step
467 call sts_set_source_tc_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
468 end subroutine hd_sts_set_source_tc_hd
469
470 function hd_get_tc_dt_hd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
471 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para/rho)
472 !and tc_k_para can depend on T=p/rho
475
476 integer, intent(in) :: ixi^l, ixo^l
477 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
478 double precision, intent(in) :: w(ixi^s,1:nw)
479 double precision :: dtnew
480
481 dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
482 end function hd_get_tc_dt_hd
483
484 subroutine hd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
485 ! move this in a different routine as in mhd if needed in more places
488
489 integer, intent(in) :: ixi^l,ixo^l
490 double precision, intent(inout) :: w(ixi^s,1:nw)
491 double precision, intent(in) :: x(ixi^s,1:ndim)
492 integer, intent(in) :: step
493
494 integer :: idir
495 logical :: flag(ixi^s,1:nw)
496 character(len=140) :: error_msg
497
498 flag=.false.
499 where(w(ixo^s,e_)<small_e) flag(ixo^s,e_)=.true.
500 if(any(flag(ixo^s,e_))) then
501 select case (small_values_method)
502 case ("replace")
503 where(flag(ixo^s,e_)) w(ixo^s,e_)=small_e
504 case ("average")
505 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
506 case default
507 ! small values error shows primitive variables
508 w(ixo^s,e_)=w(ixo^s,e_)*(hd_gamma - 1.0d0)
509 do idir = 1, ndir
510 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,rho_)
511 end do
512 write(error_msg,"(a,i3)") "Thermal conduction step ", step
513 call small_values_error(w, x, ixi^l, ixo^l, flag, error_msg)
514 end select
515 end if
516 end subroutine hd_tc_handle_small_e
517
518 ! fill in tc_fluid fields from namelist
519 subroutine tc_params_read_hd(fl)
521 type(tc_fluid), intent(inout) :: fl
522 integer :: n
523 logical :: tc_saturate=.false.
524 double precision :: tc_k_para=0d0
525
526 namelist /tc_list/ tc_saturate, tc_k_para
527
528 do n = 1, size(par_files)
529 open(unitpar, file=trim(par_files(n)), status="old")
530 read(unitpar, tc_list, end=111)
531111 close(unitpar)
532 end do
533 fl%tc_saturate = tc_saturate
534 fl%tc_k_para = tc_k_para
535
536 end subroutine tc_params_read_hd
537
538 subroutine hd_get_rho(w,x,ixI^L,ixO^L,rho)
540 integer, intent(in) :: ixi^l, ixo^l
541 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
542 double precision, intent(out) :: rho(ixi^s)
543
544 rho(ixo^s) = w(ixo^s,rho_)
545
546 end subroutine hd_get_rho
547
548!!end th cond
549!!rad cool
550 subroutine rc_params_read(fl)
552 use mod_constants, only: bigdouble
553 use mod_basic_types, only: std_len
554 type(rc_fluid), intent(inout) :: fl
555 integer :: n
556 ! list parameters
557 integer :: ncool = 4000
558
559 !> Name of cooling curve
560 character(len=std_len) :: coolcurve='JCcorona'
561
562 !> Fixed temperature not lower than tlow
563 logical :: tfix=.false.
564
565 !> Lower limit of temperature
566 double precision :: tlow=bigdouble
567
568 !> Add cooling source in a split way (.true.) or un-split way (.false.)
569 logical :: rc_split=.false.
570
571 logical :: rad_damp=.false.
572 double precision :: rad_damp_height=0.5d0
573 double precision :: rad_damp_scale=0.15d0
574 logical :: rad_newton = .false.
575 double precision :: rad_newton_trad = 0.006d0
576 double precision :: rad_newton_rhosurf = 1.d4
577 double precision :: rad_newton_pthick = 25.d0
578
579
580 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split, &
581 rad_newton, rad_newton_trad, rad_newton_rhosurf, &
582 rad_newton_pthick, rad_damp, rad_damp_height, rad_damp_scale
583
584 do n = 1, size(par_files)
585 open(unitpar, file=trim(par_files(n)), status="old")
586 read(unitpar, rc_list, end=111)
587111 close(unitpar)
588 end do
589
590 fl%ncool=ncool
591 fl%coolcurve=coolcurve
592 fl%tlow=tlow
593 fl%Tfix=tfix
594 fl%rc_split=rc_split
595 fl%rad_damp = rad_damp
596 fl%rad_damp_height = rad_damp_height
597 fl%rad_damp_scale = rad_damp_scale
598 fl%rad_newton = rad_newton
599 fl%rad_newton_trad = rad_newton_trad
600 fl%rad_newton_rhosurf = rad_newton_rhosurf
601 fl%rad_newton_pthick = rad_newton_pthick
602 end subroutine rc_params_read
603!! end rad cool
604
607 use mod_geometry, only: coordinate
610 use mod_particles, only: npayload,nusrpayload,ngridvars,num_particles,physics_type_particles
611 use mod_fld
612
613 double precision :: a,b,xfrac,yfrac
614
615 ! Initialize particles module, put here so additional gridvars and user payloads are known
616 if (hd_particles) then
617 call particles_init()
618 end if
619
620 if (.not. hd_energy) then
621 if (hd_gamma <= 0.0d0) call mpistop ("Error: hd_gamma <= 0")
622 if (hd_adiab < 0.0d0) call mpistop ("Error: hd_adiab < 0")
624 else
625 if (hd_gamma <= 0.0d0 .or. hd_gamma == 1.0d0) &
626 call mpistop ("Error: hd_gamma <= 0 or hd_gamma == 1.0")
628 inv_gamma_1=1.d0/(hd_gamma-1.d0)
630 end if
631
632 if (hd_dust) call dust_check_params()
633
634 if(hd_dust_implicit) then
635 if(.not.use_imex_scheme)then
636 call mpistop('select IMEX scheme for implicit dust update')
637 endif
638 ! implicit dust update
639 phys_implicit_update => dust_implicit_update
640 phys_evaluate_implicit => dust_evaluate_implicit
641 endif
642 if(hd_radiation_fld) then
643 if(.not.use_imex_scheme)then
644 call mpistop('select IMEX scheme for FLD radiation use')
645 endif
646 if(use_multigrid)then
647 call phys_set_mg_bounds()
648 else
649 if(.not.fld_no_mg)call mpistop('multigrid must have BCs for IMEX and FLD radiation use')
650 endif
651 if(mype==0)then
652 write(*,*)'==FLD SETUP======================'
653 write(*,*)'Using FLD with settings:'
654 write(*,*)'Using FLD with settings: hd_radiation_fld=',hd_radiation_fld
655 write(*,*)'Using FLD with settings: fld_fluxlimiter=',fld_fluxlimiter
656 write(*,*)'Using FLD with settings: fld_interaction_method=',fld_interaction_method
657 write(*,*)'Using FLD with settings: fld_opacity_law=',fld_opacity_law
658 write(*,*)'Using FLD with settings: fld_kappa0=',fld_kappa0
659 write(*,*)'Using FLD with settings: fld_opal_table=',fld_opal_table
660 write(*,*)'Using FLD with settings: fld_Radforce_split=',fld_radforce_split
661 write(*,*)'Using FLD with settings: fld_bisect_tol=',fld_bisect_tol
662 write(*,*)'Using FLD with settings: fld_diff_tol=',fld_diff_tol
663 write(*,*)'Using FLD with settings: nth_for_diff_mg=',nth_for_diff_mg
664 write(*,*)' FLD has use_imex_scheme and use_multigrid=',use_imex_scheme,use_multigrid
665 print *,'const_rad_a =',const_rad_a
666 print *,'NORMALIZED arad_norm=',arad_norm
667 print *,'NORMALIZED c_norm=',c_norm
668 print *,'const_kappae =',const_kappae
669 if(trim(fld_opacity_law).eq.'const_norm')then
670 print *,'NORMALIZED fld_kappa0 =',fld_kappa0
671 print *,'physical value (in cgs or SI) =',fld_kappa0*unit_opacity
672 endif
673 if(trim(fld_opacity_law).eq.'const')then
674 print *,'physical fld_kappa (in cgs or SI) =',fld_kappa0
675 print *,'NORMALIZED value =',fld_kappa0/unit_opacity
676 endif
677 if(fld_gamma/=hd_gamma)call mpistop("you must set fld_gamma and hd_gamma equal!")
678 write(*,*)'===FLD SETUP====================='
679 endif
680 endif
681 if(mype==0)then
682 write(*,*)'====HD run with settings===================='
683 write(*,*)'Using mod_hd_phys with settings:'
684 write(*,*)'SI_unit=',si_unit
685 write(*,*)'Dimensionality :',ndim
686 write(*,*)'vector components:',ndir
687 write(*,*)'coordinate set to type,slab:',coordinate,slab
688 write(*,*)'number of variables nw=',nw
689 write(*,*)' start index iwstart=',iwstart
690 write(*,*)'number of vector variables=',nvector
691 write(*,*)'number of stagger variables nws=',nws
692 write(*,*)'number of variables with BCs=',nwgc
693 write(*,*)'number of vars with fluxes=',nwflux
694 write(*,*)'number of vars with flux + BC=',nwfluxbc
695 write(*,*)'number of auxiliary variables=',nwaux
696 write(*,*)'number of extra vars without flux=',nwextra
697 write(*,*)'number of extra vars for wextra=',nw_extra
698 write(*,*)'number of auxiliary I/O variables=',nwauxio
699 write(*,*)'number of hd_n_tracer=',hd_n_tracer
700 write(*,*)' hd_energy=',hd_energy
701 write(*,*)' hd_gravity=',hd_gravity
702 write(*,*)' hd_viscosity=',hd_viscosity
703 write(*,*)' hd_radiative_cooling=',hd_radiative_cooling
704 write(*,*)' hd_cak_force=',hd_cak_force
705 write(*,*)' hd_radiation_fld=',hd_radiation_fld
706 write(*,*)' hd_thermal_conduction=',hd_thermal_conduction
707 write(*,*)' hd_trac=',hd_trac
708 write(*,*)' hd_dust=',hd_dust
709 write(*,*)' hd_rotating_frame=',hd_rotating_frame
710 write(*,*)' hd_particles=',hd_particles
711 if(hd_particles) then
712 write(*,*) '*****Using particles: npayload,ngridvars :', npayload,ngridvars
713 write(*,*) '*****Using particles: nusrpayload :', nusrpayload
714 write(*,*) '*****Using particles: num_particles :', num_particles
715 write(*,*) '*****Using particles: physics_type_particles=',physics_type_particles
716 end if
717 write(*,*)'number of ghostcells=',nghostcells
718 write(*,*)'number due to phys_wider_stencil=',phys_wider_stencil
719 write(*,*)'==========================================='
720 print *,'========EOS and UNITS==========='
721 print *,'SI_unit =',si_unit
722 print *,'gamma=',hd_gamma
723 print *,'eq_state_units=',eq_state_units
724 print *,'He_abundance =',he_abundance
725 print *,'RR =',rr
726 print *,'========EOS and UNITS==========='
727 print *,'unit_time =',unit_time
728 print *,'unit_length =',unit_length
729 print *,'unit_velocity =',unit_velocity
730 print *,'unit_pressure =',unit_pressure
731 print *,'unit_numberdensity =',unit_numberdensity
732 print *,'unit_density =',unit_density
733 print *,'unit_temperature =',unit_temperature
734 print *,'unit_mass =',unit_mass
735 print *,'unit_Erad =',unit_erad
736 print *,'unit_radflux =',unit_radflux
737 print *, 'CHECK that p_u ',unit_pressure,' equals ',unit_density*unit_velocity**2
738 print *, 'CHECK that L_u ',unit_length,' equals ',unit_velocity*unit_time
739 print *, 'CHECK that M_u',unit_mass,' equals ',unit_density*unit_length**3
740 print *, 'density to numberdensity has factor ',unit_density/unit_numberdensity
741 if(si_unit)then
742 print *, ' compare this to ',mp_si*(1.d0+4.d0*he_abundance)
743 else
744 print *, ' compare this to ',mp_cgs*(1.d0+4.d0*he_abundance)
745 endif
746 print *, 'pressure to n T has factor ',unit_pressure/(unit_numberdensity*unit_temperature)
747 if(si_unit)then
748 print *, ' compare this to ',kb_si*(2.d0+3.d0*he_abundance)
751 else
752 print *, ' compare this to ',kb_cgs*(2.d0+3.d0*he_abundance)
755 endif
756 if(eq_state_units)then
757 print *, 'mean molecular weight mu is =',a/b,' = ', (1.d0+4.d0*he_abundance)/(2.d0+3.d0*he_abundance)
758 xfrac=1.d0/a
759 yfrac=4.d0*he_abundance/(1.d0+4.d0*he_abundance)
760 print *, 'mass fraction hydrogen X is =',1/a,' and this equals ', 1.d0/(1.d0+4.d0*he_abundance)
761 print *, 'mass fraction helium Y is =',yfrac
762 print *, ' check that 1/mu', b/a,' is equal to 2X+3Y/4=',2.d0*xfrac+3.d0*yfrac/4.d0
763 print *, ' ratio n_e/n_p=',1.d0+2.0d0*he_abundance
764 endif
765 print *,'========UNITS==========='
766 endif
767
768 end subroutine hd_check_params
769
770 subroutine hd_physical_units
772 double precision :: mp,kb,c_lightspeed,xfrac,sigma_telectron
773 double precision :: a,b
774 ! Derive scaling units
775 if(si_unit) then
776 mp=mp_si
777 kb=kb_si
778 const_sigmasb=sigma_sb_si
779 c_lightspeed=c_si
780 sigma_telectron=sigma_te_si
781 else
782 mp=mp_cgs
783 kb=kb_cgs
784 const_sigmasb=sigma_sb_cgs
785 c_lightspeed=const_c
786 sigma_telectron=sigma_te_cgs
787 end if
788 if(eq_state_units) then
789 a=1d0+4d0*he_abundance
790 if(hd_partial_ionization) then
791 b=1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0)
792 else
793 b=2d0+3d0*he_abundance
794 end if
795 rr=1d0
796 xfrac=1.d0/a
797 else
798 a=1d0
799 b=1d0
800 rr=(1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0))/(1d0+4d0*he_abundance)
801 end if
802 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
803 if(unit_density/=1.d0) then
805 else if(unit_numberdensity/=1.d0) then
807 end if
808 if(unit_temperature/=1.d0) then
811 if(unit_length/=1.d0) then
813 else if(unit_time/=1.d0) then
815 end if
816 else if(unit_pressure/=1.d0) then
819 if(unit_length/=1.d0) then
821 else if(unit_time/=1.d0) then
823 end if
824 else if(unit_velocity/=1.d0) then
827 if(unit_length/=1.d0) then
829 else if(unit_time/=1.d0) then
831 end if
832 else if(unit_time/=1.d0) then
836 end if
837 else if(unit_temperature/=1.d0) then
838 ! units of temperature and velocity are dependent
839 if(unit_pressure/=1.d0) then
843 if(unit_length/=1.d0) then
845 else if(unit_time/=1.d0) then
847 end if
848 end if
849 else if(unit_pressure/=1.d0) then
850 if(unit_velocity/=1.d0) then
854 if(unit_length/=1.d0) then
856 else if(unit_time/=1.d0) then
858 end if
859 else if(unit_time/=0.d0) then
864 end if
865 end if
867
868 !> Units needed for radiative flux and opacity as used in FLD
869 ! normalized light speed
870 c_norm=c_lightspeed/unit_velocity
871 ! this is the radiation constant in either cgs or SI units
872 const_rad_a=4.d0*const_sigmasb/c_lightspeed
873 ! this is the dimensionless conversion factor for Erad to Trad
875 ! This is the Thomson scattering opacity in the correct units
876 ! note that the hydrogen mass fraction X=1/a in eq_state_units
877 if(eq_state_units) then
878 const_kappae=sigma_telectron*(1.d0+xfrac)/(2.0d0*mp)
879 else
880 const_kappae=0.34d0 ! specific value in cm^2/g for He=0.1 in cgs
881 endif
882 ! these are the units
886
887 end subroutine hd_physical_units
888
889 !> Returns logical argument flag where values are ok
890 subroutine hd_check_w(primitive, ixI^L, ixO^L, w, flag)
892 use mod_dust, only: dust_check_w
893
894 logical, intent(in) :: primitive
895 integer, intent(in) :: ixi^l, ixo^l
896 double precision, intent(in) :: w(ixi^s, nw)
897 logical, intent(inout) :: flag(ixi^s,1:nw)
898 double precision :: tmp(ixi^s)
899
900 flag=.false.
901
902 if (hd_energy) then
903 if (primitive) then
904 where(w(ixo^s, e_) < small_pressure) flag(ixo^s,e_) = .true.
905 else
906 tmp(ixo^s)=(hd_gamma-1.0d0)*(w(ixo^s,e_)-&
907 half*(^c&w(ixo^s,m^c_)**2+)/w(ixo^s,rho_))
908 where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_) = .true.
909 endif
910 if(hd_radiation_fld)then
911 where(w(ixo^s, r_e) < small_r_e) flag(ixo^s,r_e) = .true.
912 endif
913 end if
914
915 where(w(ixo^s, rho_) < small_density) flag(ixo^s,rho_) = .true.
916
917 if(hd_dust) call dust_check_w(ixi^l,ixo^l,w,flag)
918
919 end subroutine hd_check_w
920
921 subroutine hd_bound_fip(primitive, ixI^L, ixO^L, w)
923 logical, intent(in) :: primitive
924 integer, intent(in) :: ixi^l, ixo^l
925 double precision, intent(inout) :: w(ixi^s,1:nw)
926
927 double precision :: rho_safe(ixi^s), fip_prim(ixi^s)
928
929 if (.not. hd_fip) return
930
931 if (primitive) then
932 w(ixo^s,fip_) = min(maxfip, max(minfip, w(ixo^s,fip_)))
933 else
934 rho_safe(ixo^s) = max(w(ixo^s,rho_), small_density)
935 fip_prim(ixo^s) = w(ixo^s,fip_) / rho_safe(ixo^s)
936 fip_prim(ixo^s) = min(maxfip, max(minfip, fip_prim(ixo^s)))
937 w(ixo^s,fip_) = rho_safe(ixo^s) * fip_prim(ixo^s)
938 end if
939 end subroutine hd_bound_fip
940
941 !> Transform primitive variables into conservative ones
942 subroutine hd_to_conserved(ixI^L, ixO^L, w, x)
944 use mod_dust, only: dust_to_conserved
945 integer, intent(in) :: ixi^l, ixo^l
946 double precision, intent(inout) :: w(ixi^s, nw)
947 double precision, intent(in) :: x(ixi^s, 1:ndim)
948
949 integer :: ix^d
950
951 if (hd_fip) call hd_bound_fip(.true., ixi^l, ixo^l, w)
952
953 {do ix^db=ixomin^db,ixomax^db\}
954 if (hd_energy) then
955 ! Calculate total energy from pressure and kinetic energy
956 w(ix^d,e_)=w(ix^d, e_)*inv_gamma_1+&
957 half*(^c&w(ix^d,m^c_)**2+)*w(ix^d,rho_)
958 end if
959 ! Convert velocity to momentum
960 ^c&w(ix^d,m^c_)=w(ix^d,rho_)*w(ix^d,m^c_)\
961 if (hd_fip) w(ix^d,fip_) = w(ix^d,rho_) * w(ix^d,fip_)
962 {end do\}
963
964 if (hd_dust) then
965 call dust_to_conserved(ixi^l, ixo^l, w, x)
966 end if
967
968 end subroutine hd_to_conserved
969
970 !> Transform conservative variables into primitive ones
971 subroutine hd_to_primitive(ixI^L, ixO^L, w, x)
973 use mod_dust, only: dust_to_primitive
974 integer, intent(in) :: ixi^l, ixo^l
975 double precision, intent(inout) :: w(ixi^s, nw)
976 double precision, intent(in) :: x(ixi^s, 1:ndim)
977
978 double precision :: inv_rho
979 integer :: ix^d
980
981 if (fix_small_values) then
982 call hd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'hd_to_primitive')
983 end if
984
985 {do ix^db=ixomin^db,ixomax^db\}
986 inv_rho = 1.d0/w(ix^d,rho_)
987 if (hd_fip) w(ix^d,fip_) = w(ix^d,fip_) * inv_rho
988 ! Convert momentum to velocity
989 ^c&w(ix^d,m^c_)=w(ix^d,m^c_)*inv_rho\
990 ! Calculate pressure = (gamma-1) * (e-ek)
991 if(hd_energy) then
992 ! Compute pressure
993 w(ix^d,p_)=(hd_gamma-1.d0)*(w(ix^d,e_)&
994 -half*w(ix^d,rho_)*(^c&w(ix^d,m^c_)**2+))
995 end if
996 {end do\}
997
998 if (hd_fip) call hd_bound_fip(.true., ixi^l, ixo^l, w)
999
1000 ! Convert dust momentum to dust velocity
1001 if (hd_dust) then
1002 call dust_to_primitive(ixi^l, ixo^l, w, x)
1003 end if
1004
1005 end subroutine hd_to_primitive
1006
1007 !> Transform internal energy to total energy
1008 subroutine hd_ei_to_e(ixI^L,ixO^L,w,x)
1010 integer, intent(in) :: ixi^l, ixo^l
1011 double precision, intent(inout) :: w(ixi^s, nw)
1012 double precision, intent(in) :: x(ixi^s, 1:ndim)
1013
1014 ! Calculate total energy from internal and kinetic energy
1015 w(ixo^s,e_)=w(ixo^s,e_)+half*(^c&w(ixo^s,m^c_)**2+)/w(ixo^s,rho_)
1016
1017 end subroutine hd_ei_to_e
1018
1019 !> Transform total energy to internal energy
1020 subroutine hd_e_to_ei(ixI^L,ixO^L,w,x)
1022 integer, intent(in) :: ixi^l, ixo^l
1023 double precision, intent(inout) :: w(ixi^s, nw)
1024 double precision, intent(in) :: x(ixi^s, 1:ndim)
1025
1026 ! Calculate ei = e - ek
1027 w(ixo^s,e_)=w(ixo^s,e_)-half*(^c&w(ixo^s,m^c_)**2+)/w(ixo^s,rho_)
1028
1029 end subroutine hd_e_to_ei
1030
1031 !> Calculate v_i = m_i / rho within ixO^L
1032 subroutine hd_get_v_idim(w, x, ixI^L, ixO^L, idim, v)
1034 integer, intent(in) :: ixi^l, ixo^l, idim
1035 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
1036 double precision, intent(out) :: v(ixi^s)
1037
1038 v(ixo^s) = w(ixo^s, mom(idim)) / w(ixo^s, rho_)
1039 end subroutine hd_get_v_idim
1040
1041 !> Calculate velocity vector v_i = m_i / rho within ixO^L
1042 subroutine hd_get_v(w,x,ixI^L,ixO^L,v)
1044
1045 integer, intent(in) :: ixi^l, ixo^l
1046 double precision, intent(in) :: w(ixi^s,nw), x(ixi^s,1:^nd)
1047 double precision, intent(out) :: v(ixi^s,1:ndir)
1048
1049 integer :: idir
1050
1051 do idir=1,ndir
1052 v(ixo^s,idir) = w(ixo^s, mom(idir)) / w(ixo^s, rho_)
1053 end do
1054
1055 end subroutine hd_get_v
1056
1057 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
1058 subroutine hd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
1060 use mod_dust, only: dust_get_cmax_prim
1062
1063 integer, intent(in) :: ixi^l, ixo^l, idim
1064 ! w in primitive form
1065 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
1066 double precision, intent(inout) :: cmax(ixi^s)
1067
1068 if(hd_energy) then
1069 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(hd_gamma*w(ixo^s,p_)/w(ixo^s,rho_))
1070 else
1071 if (.not. associated(usr_set_pthermal)) then
1072 cmax(ixo^s) = hd_adiab * w(ixo^s, rho_)**hd_gamma
1073 else
1074 call usr_set_pthermal(w,x,ixi^l,ixo^l,cmax)
1075 end if
1076 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(hd_gamma*cmax(ixo^s)/w(ixo^s,rho_))
1077 end if
1078
1079 if (hd_dust) then
1080 call dust_get_cmax_prim(w, x, ixi^l, ixo^l, idim, cmax)
1081 end if
1082 end subroutine hd_get_cmax
1083
1084 !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
1085 subroutine hd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
1087 integer, intent(in) :: ixi^l,ixo^l
1088 double precision, intent(in) :: x(ixi^s,1:ndim)
1089 ! in primitive form
1090 double precision, intent(inout) :: w(ixi^s,1:nw)
1091 double precision, intent(out) :: tco_local, tmax_local
1092
1093 double precision, parameter :: trac_delta=0.25d0
1094 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
1095 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
1096 integer :: jxo^l,hxo^l
1097 integer :: jxp^l,hxp^l,ixp^l
1098 logical :: lrlt(ixi^s)
1099
1100 {^ifoned
1101 call hd_get_rfactor(w,x,ixi^l,ixi^l,te)
1102 te(ixi^s)=w(ixi^s,p_)/(te(ixi^s)*w(ixi^s,rho_))
1103
1104 tco_local=zero
1105 tmax_local=maxval(te(ixo^s))
1106 select case(hd_trac_type)
1107 case(0)
1108 w(ixi^s,tcoff_)=3.d5/unit_temperature
1109 case(1)
1110 hxo^l=ixo^l-1;
1111 jxo^l=ixo^l+1;
1112 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
1113 lrlt=.false.
1114 where(lts(ixo^s) > trac_delta)
1115 lrlt(ixo^s)=.true.
1116 end where
1117 if(any(lrlt(ixo^s))) then
1118 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
1119 end if
1120 case(2)
1121 !> iijima et al. 2021, LTRAC method
1122 ltrc=1.5d0
1123 ltrp=2.5d0
1124 ixp^l=ixo^l^ladd1;
1125 hxo^l=ixo^l-1;
1126 jxo^l=ixo^l+1;
1127 hxp^l=ixp^l-1;
1128 jxp^l=ixp^l+1;
1129 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
1130 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
1131 w(ixo^s,tcoff_)=te(ixo^s)*&
1132 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
1133 case default
1134 call mpistop("hd_trac_type not allowed for 1D simulation")
1135 end select
1136 }
1137 end subroutine hd_get_tcutoff
1138
1139 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
1140 subroutine hd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
1142 use mod_dust, only: dust_get_cmax
1143 use mod_variables
1144
1145 integer, intent(in) :: ixi^l, ixo^l, idim
1146 ! conservative left and right status
1147 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1148 ! primitive left and right status
1149 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1150 double precision, intent(in) :: x(ixi^s, 1:ndim)
1151 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1152 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1153 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
1154
1155 double precision :: wmean(ixi^s,nw)
1156 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1157 integer :: ix^d
1158
1159 select case(boundspeed)
1160 case (1)
1161 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1162 ! Methods for Fluid Dynamics" by Toro.
1163
1164 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1165 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1166 tmp3(ixo^s)=1.d0/(tmp1(ixo^s)+tmp2(ixo^s))
1167 umean(ixo^s)=(wlp(ixo^s,mom(idim))*tmp1(ixo^s)+wrp(ixo^s,mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1168
1169 if(hd_energy) then
1170 ! note usage of primitives here
1171 csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1172 csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1173 else
1174 ! note usage of conservatives here
1175 call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1176 call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1177 end if
1178
1179 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1180 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1181 (wrp(ixo^s,mom(idim))-wlp(ixo^s,mom(idim)))**2
1182
1183 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1184 if(present(cmin)) then
1185 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1186 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1187 if(h_correction) then
1188 {do ix^db=ixomin^db,ixomax^db\}
1189 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1190 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1191 {end do\}
1192 end if
1193 else
1194 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1195 end if
1196
1197 if (hd_dust) then
1198 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1199 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1200 end if
1201
1202 case (2)
1203 !if(hd_energy) then
1204 ! ! note usage of primitives here
1205 ! wmean(ixO^S,1:nwflux)=0.5d0*(wLp(ixO^S,1:nwflux)+wRp(ixO^S,1:nwflux))
1206 ! tmp1(ixO^S)=wmean(ixO^S,mom(idim))
1207 ! csoundR(ixO^S)=hd_gamma*wmean(ixO^S,p_)/wmean(ixO^S,rho_)
1208 !else
1209 ! note usage of conservatives here
1210 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1211 tmp1(ixo^s)=wmean(ixo^s,mom(idim))/wmean(ixo^s,rho_)
1212 call hd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1213 !endif
1214 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1215
1216 if(present(cmin)) then
1217 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1218 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1219 if(h_correction) then
1220 {do ix^db=ixomin^db,ixomax^db\}
1221 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1222 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1223 {end do\}
1224 end if
1225 else
1226 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1227 end if
1228
1229 if (hd_dust) then
1230 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1231 end if
1232 case (3)
1233 ! Miyoshi 2005 JCP 208, 315 equation (67)
1234 if(hd_energy) then
1235 ! note usage of primitives here
1236 csoundl(ixo^s)=hd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1237 csoundr(ixo^s)=hd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1238 else
1239 ! note usage of conservatives here
1240 call hd_get_csound2(wlc,x,ixi^l,ixo^l,csoundl)
1241 call hd_get_csound2(wrc,x,ixi^l,ixo^l,csoundr)
1242 end if
1243 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1244 if(present(cmin)) then
1245 cmin(ixo^s,1)=min(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))-csoundl(ixo^s)
1246 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1247 if(h_correction) then
1248 {do ix^db=ixomin^db,ixomax^db\}
1249 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1250 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1251 {end do\}
1252 end if
1253 else
1254 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1255 end if
1256 if (hd_dust) then
1257 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1258 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1259 end if
1260 end select
1261
1262 end subroutine hd_get_cbounds
1263
1264 !> Calculate the square of the thermal sound speed csound2 within ixO^L.
1265 !> csound2=gamma*p/rho
1266 subroutine hd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1268 integer, intent(in) :: ixi^l, ixo^l
1269 double precision, intent(in) :: w(ixi^s,nw)
1270 double precision, intent(in) :: x(ixi^s,1:ndim)
1271 double precision, intent(out) :: csound2(ixi^s)
1272
1273 call hd_get_pthermal(w,x,ixi^l,ixo^l,csound2)
1274 csound2(ixo^s)=hd_gamma*csound2(ixo^s)/w(ixo^s,rho_)
1275
1276 end subroutine hd_get_csound2
1277
1278 !> Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L
1279 subroutine hd_get_pthermal(w, x, ixI^L, ixO^L, pth)
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):: pth(ixi^s)
1288 integer :: iw, ix^d
1289
1290 if (hd_energy) then
1291 pth(ixo^s) = (hd_gamma - 1.0d0) * (w(ixo^s, e_) - &
1292 hd_kin_en(w, ixi^l, ixo^l))
1293 else
1294 if (.not. associated(usr_set_pthermal)) then
1295 pth(ixo^s) = hd_adiab * w(ixo^s, rho_)**hd_gamma
1296 else
1297 call usr_set_pthermal(w,x,ixi^l,ixo^l,pth)
1298 end if
1299 end if
1300
1301 if (fix_small_values) then
1302 {do ix^db= ixo^lim^db\}
1303 if(pth(ix^d)<small_pressure) then
1304 pth(ix^d)=small_pressure
1305 endif
1306 {enddo^d&\}
1307 else if (check_small_values) then
1308 {do ix^db= ixo^lim^db\}
1309 if(pth(ix^d)<small_pressure) then
1310 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1311 " encountered when call hd_get_pthermal"
1312 write(*,*) "Iteration: ", it, " Time: ", global_time
1313 write(*,*) "Location: ", x(ix^d,:)
1314 write(*,*) "Cell number: ", ix^d
1315 do iw=1,nw
1316 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1317 end do
1318 ! use erroneous arithmetic operation to crash the run
1319 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1320 write(*,*) "Saving status at the previous time step"
1321 crash=.true.
1322 end if
1323 {enddo^d&\}
1324 end if
1325
1326 end subroutine hd_get_pthermal
1327
1328 !> Calculate modified squared sound speed for FLD
1329 !> NOTE: only for diagnostic purposes, unused subroutine
1330 subroutine hd_get_csrad2(w,x,ixI^L,ixO^L,csound)
1332
1333 integer, intent(in) :: ixi^l, ixo^l
1334 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
1335 double precision, intent(out):: csound(ixi^s)
1336
1337 double precision :: wprim(ixi^s, nw)
1338
1339 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1340 call hd_to_primitive(ixi^l,ixo^l,wprim,x)
1341 call hd_get_csrad2_prim(wprim,x,ixi^l,ixo^l,csound)
1342
1343 end subroutine hd_get_csrad2
1344
1345 !> Calculate modified squared sound speed for FLD
1346 !> NOTE: w is primitive on entry here!
1347 !> NOTE: used in FLD module as phys_get_csrad2
1348 subroutine hd_get_csrad2_prim(w,x,ixI^L,ixO^L,csound)
1350
1351 integer, intent(in) :: ixi^l, ixo^l
1352 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
1353 double precision, intent(out):: csound(ixi^s)
1354
1355 integer :: ix^d
1356 double precision :: inv_rho
1357 double precision :: prad_tensor(ixi^s, 1:ndim, 1:ndim)
1358 double precision :: prad_max(ixi^s)
1359
1360 call hd_get_pradiation_from_prim(w, x, ixi^l, ixo^l, prad_tensor)
1361
1362 {do ix^db=ixomin^db,ixomax^db \}
1363 inv_rho=1.d0/w(ix^d,rho_)
1364 prad_max(ix^d) = maxval(prad_tensor(ix^d,:,:))
1365 csound(ix^d)=(hd_gamma*w(ix^d,p_)+prad_max(ix^d))*inv_rho
1366 {end do\}
1367
1368 if(minval(csound(ixo^s))<smalldouble)then
1369 print *,'issue with squared speed and rad pressure'
1370 print *,minval(csound(ixo^s))
1371 print *,minval(prad_max(ixo^s))
1372 call mpistop("negative squared speed in get_csrad2 for dt")
1373 endif
1374
1375 end subroutine hd_get_csrad2_prim
1376
1377 !> Calculate radiation pressure within ixO^L
1378 !> NOTE: w is primitive on entry here!
1379 !> NOTE: used in FLD module as it is called from phys_get_csrad2
1380 subroutine hd_get_pradiation_from_prim(w, x, ixI^L, ixO^L, prad)
1382 use mod_fld
1383 integer, intent(in) :: ixi^l, ixo^l
1384 double precision, intent(in) :: w(ixi^s, 1:nw)
1385 double precision, intent(in) :: x(ixi^s, 1:ndim)
1386 double precision, intent(out):: prad(ixi^s, 1:ndim, 1:ndim)
1387
1388 call fld_get_radpress(w, x, ixi^l, ixo^l, prad)
1389
1390 end subroutine hd_get_pradiation_from_prim
1391
1392 !> calculates the sum of the gas pressure and max Prad tensor element
1393 !> NOTE: only for diagnostic purposes, unused subroutine
1394 subroutine hd_get_pthermal_plus_pradiation(w, x, ixI^L, ixO^L, pth_plus_prad)
1396 integer, intent(in) :: ixi^l, ixo^l
1397 double precision, intent(in) :: w(ixi^s, 1:nw)
1398 double precision, intent(in) :: x(ixi^s, 1:ndim)
1399 double precision, intent(out):: pth_plus_prad(ixi^s)
1400
1401 double precision :: wprim(ixi^s, 1:nw)
1402 double precision :: prad_tensor(ixi^s, 1:ndim, 1:ndim)
1403 double precision :: prad_max(ixi^s)
1404 integer :: ix^d
1405
1406 wprim(ixi^s,1:nw)=w(ixi^s,1:nw)
1407 call hd_to_primitive(ixi^l,ixo^l,wprim,x)
1408 call hd_get_pradiation_from_prim(wprim, x, ixi^l, ixo^l, prad_tensor)
1409 {do ix^d = ixomin^d,ixomax^d\}
1410 prad_max(ix^d) = maxval(prad_tensor(ix^d,:,:))
1411 {enddo\}
1412 pth_plus_prad(ixo^s) = wprim(ixo^s,p_) + prad_max(ixo^s)
1413 end subroutine hd_get_pthermal_plus_pradiation
1414
1415 !> Calculates radiation temperature
1416 subroutine hd_get_trad(w, x, ixI^L, ixO^L, trad)
1418 use mod_constants
1419
1420 integer, intent(in) :: ixi^l, ixo^l
1421 double precision, intent(in) :: w(ixi^s, 1:nw)
1422 double precision, intent(in) :: x(ixi^s, 1:ndim)
1423 double precision, intent(out):: trad(ixi^s)
1424
1425 trad(ixi^s) = (w(ixi^s,r_e)/arad_norm)**(1.d0/4.d0)
1426
1427 end subroutine hd_get_trad
1428
1429 !> Calculate temperature=p/rho when in e_ the total energy is stored
1430 subroutine hd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1432 integer, intent(in) :: ixi^l, ixo^l
1433 double precision, intent(in) :: w(ixi^s, 1:nw)
1434 double precision, intent(in) :: x(ixi^s, 1:ndim)
1435 double precision, intent(out):: res(ixi^s)
1436
1437 double precision :: r(ixi^s)
1438
1439 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1440 call hd_get_pthermal(w, x, ixi^l, ixo^l, res)
1441 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1442 end subroutine hd_get_temperature_from_etot
1443
1444 !> Calculate temperature=p/rho when in e_ the internal energy is stored
1445 subroutine hd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1447 integer, intent(in) :: ixi^l, ixo^l
1448 double precision, intent(in) :: w(ixi^s, 1:nw)
1449 double precision, intent(in) :: x(ixi^s, 1:ndim)
1450 double precision, intent(out):: res(ixi^s)
1451
1452 double precision :: r(ixi^s)
1453
1454 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1455 res(ixo^s) = (hd_gamma - 1.0d0) * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1456 end subroutine hd_get_temperature_from_eint
1457
1458 !> Calculate temperature=p/rho when in we work in primitives (hence e_ is p_)
1459 subroutine hd_get_temperature_from_prim(w, x, ixI^L, ixO^L, res)
1461 integer, intent(in) :: ixi^l, ixo^l
1462 double precision, intent(in) :: w(ixi^s, 1:nw)
1463 double precision, intent(in) :: x(ixi^s, 1:ndim)
1464 double precision, intent(out):: res(ixi^s)
1465
1466 double precision :: r(ixi^s)
1467
1468 call hd_get_rfactor(w,x,ixi^l,ixo^l,r)
1469 res(ixo^s)=w(ixo^s,p_)/(r(ixo^s)*w(ixo^s,rho_))
1470
1471 end subroutine hd_get_temperature_from_prim
1472
1473 ! Calculate flux f_idim[iw]
1474 subroutine hd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1476 use mod_dust, only: dust_get_flux_prim
1477
1478 integer, intent(in) :: ixi^l, ixo^l, idim
1479 ! conservative w
1480 double precision, intent(in) :: wc(ixi^s, 1:nw)
1481 ! primitive w
1482 double precision, intent(in) :: w(ixi^s, 1:nw)
1483 double precision, intent(in) :: x(ixi^s, 1:ndim)
1484 double precision, intent(out) :: f(ixi^s, nwflux)
1485
1486 double precision :: pth(ixi^s)
1487 integer :: ix^db
1488
1489 if (hd_energy) then
1490 {do ix^db=ixomin^db,ixomax^db\}
1491 f(ix^d,rho_)=w(ix^d,mom(idim))*w(ix^d,rho_)
1492 ! Momentum flux is v_i*m_i, +p in direction idim
1493 ^c&f(ix^d,m^c_)=w(ix^d,mom(idim))*wc(ix^d,m^c_)\
1494 f(ix^d,mom(idim))=f(ix^d,mom(idim))+w(ix^d,p_)
1495 ! Energy flux is v_i*(e + p)
1496 f(ix^d,e_)=w(ix^d,mom(idim))*(wc(ix^d,e_)+w(ix^d,p_))
1497 {end do\}
1498 else
1499 call hd_get_pthermal(wc, x, ixi^l, ixo^l, pth)
1500 {do ix^db=ixomin^db,ixomax^db\}
1501 f(ix^d,rho_)=w(ix^d,mom(idim))*w(ix^d,rho_)
1502 ! Momentum flux is v_i*m_i, +p in direction idim
1503 ^c&f(ix^d,m^c_)=w(ix^d,mom(idim))*wc(ix^d,m^c_)\
1504 f(ix^d,mom(idim))=f(ix^d,mom(idim))+pth(ix^d)
1505 {end do\}
1506 end if
1507
1508 if(hd_radiation_fld)then
1509 {do ix^db=ixomin^db,ixomax^db\}
1510 ! advection of radiation enery v_i*r_e
1511 f(ix^d,r_e)=w(ix^d,mom(idim))*wc(ix^d,r_e)
1512 {end do\}
1513 endif
1514
1515 do ix1 = 1, hd_n_tracer
1516 f(ixo^s, tracer(ix1)) = w(ixo^s,mom(idim)) * w(ixo^s, tracer(ix1))
1517 end do
1518
1519 if (hd_fip) then
1520 f(ixo^s,fip_) = w(ixo^s,mom(idim)) * wc(ixo^s,fip_)
1521 end if
1522
1523 ! Dust fluxes
1524 if (hd_dust) then
1525 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1526 end if
1527
1528 end subroutine hd_get_flux
1529
1530 !> Add geometrical source terms to w
1531 !>
1532 !> Notice that the expressions of the geometrical terms depend only on ndir,
1533 !> not ndim. Eg, they are the same in 2.5D and in 3D, for any geometry.
1534 !>
1535 subroutine hd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1540 use mod_geometry
1541 integer, intent(in) :: ixi^l, ixo^l
1542 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:ndim)
1543 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw),w(ixi^s, 1:nw)
1544 double precision :: pth(ixi^s), source(ixi^s), minrho
1545 integer :: iw,idir, h1x^l{^nooned, h2x^l}
1546 integer :: mr_,mphi_ ! Polar var. names
1547 integer :: irho, ifluid, n_fluids
1548 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1549
1550 if (hd_dust) then
1551 n_fluids = 1 + dust_n_species
1552 else
1553 n_fluids = 1
1554 end if
1555
1556 select case (coordinate)
1557
1559 !the user provides the functions of exp_factor and del_exp_factor
1560 if(associated(usr_set_surface)) call usr_set_surface(ixi^l,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
1561 if(hd_energy) then
1562 source(ixo^s)=wprim(ixo^s, p_)
1563 else
1564 if(.not. associated(usr_set_pthermal)) then
1565 source(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1566 else
1567 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1568 end if
1569 end if
1570 source(ixo^s) = source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1571 w(ixo^s,mom(1)) = w(ixo^s,mom(1)) + qdt*source(ixo^s)
1572
1573 case (cylindrical)
1574 do ifluid = 0, n_fluids-1
1575 ! s[mr]=(pthermal+mphi**2/rho)/radius
1576 if (ifluid == 0) then
1577 ! gas
1578 irho = rho_
1579 mr_ = mom(r_)
1580 if(phi_>0) mphi_ = mom(phi_)
1581 if(hd_energy) then
1582 source(ixo^s)=wprim(ixo^s, p_)
1583 else
1584 if(.not. associated(usr_set_pthermal)) then
1585 source(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1586 else
1587 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1588 end if
1589 end if
1590 minrho = 0.0d0
1591 else
1592 ! dust : no pressure
1593 irho = dust_rho(ifluid)
1594 mr_ = dust_mom(r_, ifluid)
1595 if(phi_>0) mphi_ = dust_mom(phi_, ifluid)
1596 source(ixi^s) = zero
1597 minrho = 0.0d0
1598 end if
1599 if(phi_ > 0) then
1600 where (wct(ixo^s, irho) > minrho)
1601 source(ixo^s) = source(ixo^s) + wct(ixo^s,mphi_)*wprim(ixo^s,mphi_)
1602 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt*source(ixo^s)/x(ixo^s,r_)
1603 end where
1604 ! s[mphi]=(-mphi*vr)/radius
1605 where (wct(ixo^s, irho) > minrho)
1606 source(ixo^s) = -wct(ixo^s, mphi_) * wprim(ixo^s, mr_)
1607 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * source(ixo^s) / x(ixo^s, r_)
1608 end where
1609 else
1610 ! s[mr]=2pthermal/radius
1611 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1612 end if
1613 end do
1614 case (spherical)
1615 if (hd_dust) then
1616 call mpistop("Dust geom source terms not implemented yet with spherical geometries")
1617 end if
1618 mr_ = mom(r_)
1619 if(phi_>0) mphi_ = mom(phi_)
1620 h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
1621 if(hd_energy) then
1622 pth(ixo^s)=wprim(ixo^s, p_)
1623 else
1624 if(.not. associated(usr_set_pthermal)) then
1625 pth(ixo^s)=hd_adiab * wprim(ixo^s, rho_)**hd_gamma
1626 else
1627 call usr_set_pthermal(wct,x,ixi^l,ixo^l,pth)
1628 end if
1629 end if
1630 ! s[mr]=((vtheta**2+vphi**2)*rho+2*p)/r
1631 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1632 *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
1633 /block%dvolume(ixo^s)
1634 do idir = 2, ndir
1635 source(ixo^s) = source(ixo^s) + wprim(ixo^s, mom(idir))**2 * wprim(ixo^s, rho_)
1636 end do
1637 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, 1)
1638
1639 {^nooned
1640 ! s[mtheta]=-(vr*vtheta*rho)/r+cot(theta)*(vphi**2*rho+p)/r
1641 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1642 * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
1643 / block%dvolume(ixo^s)
1644 if (ndir == 3) then
1645 source(ixo^s) = source(ixo^s) + (wprim(ixo^s, mom(3))**2 * wprim(ixo^s, rho_)) / tan(x(ixo^s, 2))
1646 end if
1647 source(ixo^s) = source(ixo^s) - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)
1648 w(ixo^s, mom(2)) = w(ixo^s, mom(2)) + qdt * source(ixo^s) / x(ixo^s, 1)
1649
1650 if (ndir == 3) then
1651 ! s[mphi]=-(vphi*vr/rho)/r-cot(theta)*(vtheta*vphi/rho)/r
1652 source(ixo^s) = -(wprim(ixo^s, mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)&
1653 - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mom(3))) * wprim(ixo^s, rho_) / tan(x(ixo^s, 2))
1654 w(ixo^s, mom(3)) = w(ixo^s, mom(3)) + qdt * source(ixo^s) / x(ixo^s, 1)
1655 end if
1656 }
1657 end select
1658
1659 if (hd_rotating_frame) then
1660 if (hd_dust) then
1661 call mpistop("Rotating frame not implemented yet with dust")
1662 else
1663 call rotating_frame_add_source(qdt,dtfactor,ixi^l,ixo^l,wprim,w,x)
1664 end if
1665 end if
1666
1667 end subroutine hd_add_source_geom
1668
1669 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
1670 subroutine hd_add_source(qdt,dtfactor, ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1675 use mod_usr_methods, only: usr_gravity
1677 use mod_cak_force, only: cak_add_source
1678
1679 integer, intent(in) :: ixi^l, ixo^l
1680 double precision, intent(in) :: qdt, dtfactor
1681 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw), x(ixi^s, 1:ndim)
1682 double precision, intent(inout) :: w(ixi^s, 1:nw)
1683 logical, intent(in) :: qsourcesplit
1684 logical, intent(inout) :: active
1685
1686 double precision :: gravity_field(ixi^s, 1:ndim)
1687 integer :: idust, idim
1688
1689 if(hd_dust .and. .not. hd_dust_implicit) then
1690 call dust_add_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1691 end if
1692
1693 if(hd_radiative_cooling) then
1694 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1695 qsourcesplit,active, rc_fl)
1696 end if
1697
1698 if(hd_viscosity) then
1699 call viscosity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1700 hd_energy,qsourcesplit,active)
1701 end if
1702
1703 if (hd_gravity) then
1704 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1705 hd_energy,qsourcesplit,active)
1706
1707 if (hd_dust .and. qsourcesplit .eqv. grav_split) then
1708 active = .true.
1709
1710 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1711 do idust = 1, dust_n_species
1712 do idim = 1, ndim
1713 w(ixo^s, dust_mom(idim, idust)) = w(ixo^s, dust_mom(idim, idust)) &
1714 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s, dust_rho(idust))
1715 end do
1716 end do
1717 end if
1718 end if
1719
1720 if (hd_cak_force) then
1721 call cak_add_source(qdt,ixi^l,ixo^l,wct,w,x,hd_energy,qsourcesplit,active)
1722 end if
1723
1724 ! This is where the radiation force and heating/cooling are added
1725 if (hd_radiation_fld) then
1726 call hd_add_radiation_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,qsourcesplit,active)
1727 endif
1728
1729 if(hd_partial_ionization) then
1730 if(.not.qsourcesplit) then
1731 active = .true.
1732 call hd_update_temperature(ixi^l,ixo^l,wct,w,x)
1733 end if
1734 end if
1735
1736 end subroutine hd_add_source
1737
1738 subroutine hd_add_radiation_source(qdt,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1739 use mod_constants
1741 use mod_usr_methods
1742 use mod_fld
1743
1744 integer, intent(in) :: ixi^l, ixo^l
1745 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
1746 double precision, intent(in) :: wct(ixi^s,1:nw),wctprim(ixi^s,1:nw)
1747 double precision, intent(inout) :: w(ixi^s,1:nw)
1748 logical, intent(in) :: qsourcesplit
1749 logical, intent(inout) :: active
1750
1751 ! add radiation force and work done by it, changes momentum and gas energy
1752 ! handle photon tiring, heating and cooling exchange between gas and radiation field
1753 call add_fld_rad_force(qdt,ixi^l,ixo^l,wct,wctprim,w,x,qsourcesplit,active)
1754
1755 end subroutine hd_add_radiation_source
1756
1757 subroutine hd_get_dt(wprim, ixI^L, ixO^L, dtnew, dx^D, x)
1759 use mod_dust, only: dust_get_dt
1761 use mod_gravity, only: gravity_get_dt
1762 use mod_cak_force, only: cak_get_dt
1763 use mod_fld, only: fld_radforce_get_dt
1764
1765 integer, intent(in) :: ixi^l, ixo^l
1766 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
1767 double precision, intent(in) :: wprim(ixi^s, 1:nw)
1768 double precision, intent(inout) :: dtnew
1769
1770 dtnew = bigdouble
1771
1772 if(hd_dust) then
1773 call dust_get_dt(wprim, ixi^l, ixo^l, dtnew, dx^d, x)
1774 end if
1775
1776 if(hd_viscosity) then
1777 call viscosity_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1778 end if
1779
1780 if(hd_gravity) then
1781 call gravity_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1782 end if
1783
1784 if (hd_cak_force) then
1785 call cak_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1786 end if
1787
1788 if(hd_radiation_fld) then
1789 call fld_radforce_get_dt(wprim,ixi^l,ixo^l,dtnew,dx^d,x)
1790 endif
1791
1792 end subroutine hd_get_dt
1793
1794 function hd_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
1795 use mod_global_parameters, only: nw, ndim
1796 integer, intent(in) :: ixi^l, ixo^l
1797 double precision, intent(in) :: w(ixi^s, nw)
1798 double precision :: ke(ixo^s)
1799 double precision, intent(in), optional :: inv_rho(ixo^s)
1800
1801 if (present(inv_rho)) then
1802 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * inv_rho
1803 else
1804 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) / w(ixo^s, rho_)
1805 end if
1806 end function hd_kin_en
1807
1808 function hd_inv_rho(w, ixI^L, ixO^L) result(inv_rho)
1809 use mod_global_parameters, only: nw, ndim
1810 integer, intent(in) :: ixi^l, ixo^l
1811 double precision, intent(in) :: w(ixi^s, nw)
1812 double precision :: inv_rho(ixo^s)
1813
1814 ! Can make this more robust
1815 inv_rho = 1.0d0 / w(ixo^s, rho_)
1816 end function hd_inv_rho
1817
1818 subroutine hd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1819 ! handles hydro (density,pressure,velocity) bootstrapping
1820 ! any negative dust density is flagged as well (and throws an error)
1821 ! small_values_method=replace also for dust
1825 logical, intent(in) :: primitive
1826 integer, intent(in) :: ixi^l,ixo^l
1827 double precision, intent(inout) :: w(ixi^s,1:nw)
1828 double precision, intent(in) :: x(ixi^s,1:ndim)
1829 character(len=*), intent(in) :: subname
1830
1831 integer :: n,idir
1832 logical :: flag(ixi^s,1:nw)
1833
1834 call hd_check_w(primitive, ixi^l, ixo^l, w, flag)
1835
1836 if (any(flag)) then
1837 select case (small_values_method)
1838 case ("replace")
1839 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
1840 do idir = 1, ndir
1841 if(small_values_fix_iw(mom(idir))) then
1842 where(flag(ixo^s,rho_)) w(ixo^s, mom(idir)) = 0.0d0
1843 end if
1844 end do
1845 if(hd_radiation_fld)then
1846 if (small_values_fix_iw(r_e)) then
1847 where(flag(ixo^s,r_e)) w(ixo^s,r_e) = small_r_e
1848 end if
1849 end if
1850 if(hd_energy)then
1851 if(small_values_fix_iw(e_)) then
1852 if(primitive) then
1853 where(flag(ixo^s,rho_)) w(ixo^s, p_) = small_pressure
1854 else
1855 where(flag(ixo^s,rho_)) w(ixo^s, e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1856 endif
1857 end if
1858 endif
1859
1860 if(hd_energy) then
1861 if(primitive) then
1862 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
1863 else
1864 where(flag(ixo^s,e_))
1865 ! Add kinetic energy
1866 w(ixo^s,e_) = small_e + hd_kin_en(w,ixi^l,ixo^l)
1867 end where
1868 end if
1869 end if
1870
1871 if(hd_dust)then
1872 do n=1,dust_n_species
1873 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1874 do idir = 1, ndir
1875 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1876 enddo
1877 enddo
1878 endif
1879 case ("average")
1880 if(primitive)then
1881 ! averaging for all primitive fields, including dust
1882 call small_values_average(ixi^l, ixo^l, w, x, flag)
1883 else
1884 ! do averaging of density
1885 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
1886 if(hd_energy) then
1887 ! do averaging of pressure
1888 w(ixi^s,p_)=(hd_gamma-1.d0)*(w(ixi^s,e_) &
1889 -0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_))
1890 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
1891 w(ixi^s,e_)=w(ixi^s,p_)/(hd_gamma-1.d0) &
1892 +0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_)
1893 end if
1894 if(hd_radiation_fld) then
1895 ! do averaging of radiative energy density
1896 call small_values_average(ixi^l, ixo^l, w, x, flag, r_e)
1897 endif
1898 if(hd_dust)then
1899 do n=1,dust_n_species
1900 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1901 do idir = 1, ndir
1902 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1903 enddo
1904 enddo
1905 endif
1906 endif
1907 case default
1908 if(.not.primitive) then
1909 !convert w to primitive
1910 ! Calculate pressure = (gamma-1) * (e-ek)
1911 if(hd_energy) then
1912 w(ixo^s,p_)=(hd_gamma-1.d0)*(w(ixo^s,e_)-hd_kin_en(w,ixi^l,ixo^l))
1913 end if
1914 ! Convert gas momentum to velocity
1915 do idir = 1, ndir
1916 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir))/w(ixo^s,rho_)
1917 end do
1918 end if
1919 ! NOTE: dust entries may still have conserved values here
1920 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1921 end select
1922 end if
1923 if (hd_fip) call hd_bound_fip(primitive, ixi^l, ixo^l, w)
1924 end subroutine hd_handle_small_values
1925
1926 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1929 integer, intent(in) :: ixi^l, ixo^l
1930 double precision, intent(in) :: w(ixi^s,1:nw)
1931 double precision, intent(in) :: x(ixi^s,1:ndim)
1932 double precision, intent(out):: rfactor(ixi^s)
1933
1934 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1935
1936 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1937 ! assume the first and second ionization of Helium have the same degree
1938 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
1939
1940 end subroutine rfactor_from_temperature_ionization
1941
1942 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1944 integer, intent(in) :: ixi^l, ixo^l
1945 double precision, intent(in) :: w(ixi^s,1:nw)
1946 double precision, intent(in) :: x(ixi^s,1:ndim)
1947 double precision, intent(out):: rfactor(ixi^s)
1948
1949 rfactor(ixo^s)=rr
1950
1951 end subroutine rfactor_from_constant_ionization
1952
1953 subroutine hd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1956
1957 integer, intent(in) :: ixi^l, ixo^l
1958 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1959 double precision, intent(inout) :: w(ixi^s,1:nw)
1960
1961 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1962
1963 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1964
1965 call hd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1966
1967 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1968 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1969
1970 end subroutine hd_update_temperature
1971
1972end 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.
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:33
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:476
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:25
double precision, public fld_diff_tol
Tolerance for radiative Energy diffusion.
Definition mod_fld.t:27
double precision, public fld_gamma
A copy of (m)hd_gamma.
Definition mod_fld.t:46
character(len=40) fld_fluxlimiter
flux limiter choice
Definition mod_fld.t:38
character(len=40) fld_opal_table
Definition mod_fld.t:36
double precision, public fld_kappa0
Opacity value when using constant opacity.
Definition mod_fld.t:22
character(len=40) fld_opacity_law
switches for opacity
Definition mod_fld.t:35
character(len=40) fld_interaction_method
Which method to find the root for the energy interaction polynomial.
Definition mod_fld.t:44
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:197
logical fld_radforce_split
source split for energy interact and radforce:
Definition mod_fld.t:18
subroutine, public fld_init(r_gamma)
Initialising FLD-module Read opacities Initialise Multigrid and adimensionalise kappa.
Definition mod_fld.t:86
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:334
integer nth_for_diff_mg
diffusion coefficient stencil control
Definition mod_fld.t:42
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 arad_norm
Normalised radiation constant.
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.
double precision const_rad_a
Physical factors useful for radiation fld.
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 c_norm
Normalised speed of light.
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
double precision unit_erad
Physical scaling factor for radiation energy density.
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:87
subroutine, public hd_get_temperature_from_etot(w, x, ixil, ixol, res)
Calculate temperature=p/rho when in e_ the total energy is stored.
integer, public, protected hd_trac_type
Definition mod_hd_phys.t:97
logical, public, protected hd_particles
Whether particles module is added.
Definition mod_hd_phys.t:39
logical, public, protected hd_fip
Whether FIP passive scalar is enabled.
Definition mod_hd_phys.t:81
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:84
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)
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:93
subroutine, public hd_get_trad(w, x, ixil, ixol, trad)
Calculates radiation temperature.
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.
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:96
integer, public, protected fip_
Index of the FIP passive scalar rho*fip in conserved form, fip in primitive form.
Definition mod_hd_phys.t:78
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, public viscosity_get_dt(wprim, ixil, ixol, dtnew, dxd, x)
procedure(sub_add_source), pointer, public viscosity_add_source
subroutine, public viscosity_init(phys_wider_stencil)
Initialize the module.