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