MPI-AMRVAC 3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code (development version)
Loading...
Searching...
No Matches
mod_rhd_phys.t
Go to the documentation of this file.
1!> Radiation-Hydrodynamics physics module
2!> Module aims at solving the Hydrodynamic equations toghether with
3!> the zeroth moment of the radiative transfer equation. A closure is
4!> provided by the flux limited diffusion (FLD)-approximation in the mod_fld.t module. See
5!> [1]Moens, N., Sundqvist, J. O., El Mellah, I., Poniatowski, L., Teunissen, J., and Keppens, R.,
6!> “Radiation-hydrodynamics with MPI-AMRVAC . Flux-limited diffusion”,
7!> <i>Astronomy and Astrophysics</i>, vol. 657, 2022. doi:10.1051/0004-6361/202141023.
8!> For more information.
9!> Another possible closure in the works is the anisotropic flux limited diffusion approximation (AFLD) described in mod_afld.t.
10
15 use mod_physics
16 use mod_comm_lib, only: mpistop
17 implicit none
18 private
19
20 !> Whether an energy equation is used
21 logical, public, protected :: rhd_energy = .true.
22
23 !> Whether thermal conduction is added
24 logical, public, protected :: rhd_thermal_conduction = .false.
25 type(tc_fluid), allocatable, public :: tc_fl
26 type(te_fluid), allocatable, public :: te_fl_rhd
27
28 !> Whether radiative cooling is added
29 logical, public, protected :: rhd_radiative_cooling = .false.
30 type(rc_fluid), allocatable, public :: rc_fl
31
32 !> Whether dust is added
33 logical, public, protected :: rhd_dust = .false.
34
35 !> Whether viscosity is added
36 logical, public, protected :: rhd_viscosity = .false.
37
38 !> Whether gravity is added
39 logical, public, protected :: rhd_gravity = .false.
40
41 !> Whether particles module is added
42 logical, public, protected :: rhd_particles = .false.
43
44 !> Whether rotating frame is activated
45 logical, public, protected :: rhd_rotating_frame = .false.
46
47 !> Number of tracer species
48 integer, public, protected :: rhd_n_tracer = 0
49
50 !> Index of the density (in the w array)
51 integer, public, protected :: rho_
52
53 !> Indices of the momentum density
54 integer, allocatable, public, protected :: mom(:)
55
56 !> Indices of the tracers
57 integer, allocatable, public, protected :: tracer(:)
58
59 !> Index of the energy density (-1 if not present)
60 integer, public, protected :: e_
61
62 !> Index of the gas pressure (-1 if not present) should equal e_
63 integer, public, protected :: p_
64
65 !> Index of the radiation energy
66 integer, public, protected :: r_e
67
68 !> Indices of temperature
69 integer, public, protected :: te_
70
71 !> Index of the cutoff temperature for the TRAC method
72 integer, public, protected :: tcoff_
73
74 !> The adiabatic index
75 double precision, public :: rhd_gamma = 5.d0/3.0d0
76
77 !> The adiabatic constant
78 double precision, public :: rhd_adiab = 1.0d0
79
80 !> The small_est allowed energy
81 double precision, protected :: small_e
82
83 !> The smallest allowed radiation energy
84 double precision, public, protected :: small_r_e = 0.d0
85
86 !> Whether TRAC method is used
87 logical, public, protected :: rhd_trac = .false.
88 integer, public, protected :: rhd_trac_type = 1
89
90 !> Helium abundance over Hydrogen
91 double precision, public, protected :: he_abundance=0.1d0
92
93 !> Formalism to treat radiation
94 character(len=8), public :: rhd_radiation_formalism = 'fld'
95
96 !> In the case of no rhd_energy, how to compute pressure
97 character(len=8), public :: rhd_pressure = 'Trad'
98
99 !> Treat radiation fld_Rad_force
100 logical, public, protected :: rhd_radiation_force = .true.
101
102 !> Treat radiation-gas energy interaction
103 logical, public, protected :: rhd_energy_interact = .true.
104
105 !> Treat radiation energy diffusion
106 logical, public, protected :: rhd_radiation_diffusion = .true.
107
108 !> Treat radiation advection
109 logical, public, protected :: rhd_radiation_advection = .true.
110
111 !> Whether plasma is partially ionized
112 logical, public, protected :: rhd_partial_ionization = .false.
113
114 !> Do a running mean over the radiation pressure when determining dt
115 logical, protected :: radio_acoustic_filter = .false.
116 integer, protected :: size_ra_filter = 1
117
118 !> kb/(m_p mu)* 1/a_rad**4,
119 double precision, public :: kbmpmua4
120
121 !> Use the speed of light to calculate the timestep, usefull for debugging
122 logical :: dt_c = .false.
123
124 !> Ionization fraction of H
125 !> H_ion_fr = H+/(H+ + H)
126 double precision, public, protected :: h_ion_fr=1d0
127 !> Ionization fraction of He
128 !> He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
129 double precision, public, protected :: he_ion_fr=1d0
130 !> Ratio of number He2+ / number He+ + He2+
131 !> He_ion_fr2 = He2+/(He2+ + He+)
132 double precision, public, protected :: he_ion_fr2=1d0
133 ! used for eq of state when it is not defined by units,
134 ! the units do not contain terms related to ionization fraction
135 ! and it is p = RR * rho * T
136 double precision, public, protected :: rr=1d0
137 ! remove the below flag and assume default value = .false.
138 ! when eq state properly implemented everywhere
139 ! and not anymore through units
140 logical, public, protected :: eq_state_units = .true.
141
142 procedure(sub_get_pthermal), pointer :: rhd_get_rfactor => null()
143 ! Public methods
144 public :: rhd_phys_init
145 public :: rhd_kin_en
146 public :: rhd_get_pthermal
147 public :: rhd_get_pradiation
148 public :: rhd_get_ptot
149 public :: rhd_get_csound2
150 public :: rhd_to_conserved
151 public :: rhd_to_primitive
152 public :: rhd_check_params
153 public :: rhd_check_w
154 public :: rhd_get_tgas
155 public :: rhd_get_trad
156 public :: rhd_set_mg_bounds
157
158contains
159
160 !> Read this module's parameters from a file
161 subroutine rhd_read_params(files)
163 character(len=*), intent(in) :: files(:)
164 integer :: n
165
166 namelist /rhd_list/ rhd_energy, rhd_pressure, rhd_n_tracer, rhd_gamma, rhd_adiab, &
172 rhd_radiation_advection, radio_acoustic_filter, size_ra_filter, dt_c, rhd_partial_ionization
173
174 do n = 1, size(files)
175 open(unitpar, file=trim(files(n)), status="old")
176 read(unitpar, rhd_list, end=111)
177111 close(unitpar)
178 end do
179
180 end subroutine rhd_read_params
181
182 !> Write this module's parameters to a snapsoht
183 subroutine rhd_write_info(fh)
185 integer, intent(in) :: fh
186 integer, parameter :: n_par = 1
187 double precision :: values(n_par)
188 character(len=name_len) :: names(n_par)
189 integer, dimension(MPI_STATUS_SIZE) :: st
190 integer :: er
191
192 call mpi_file_write(fh, n_par, 1, mpi_integer, st, er)
193
194 names(1) = "gamma"
195 values(1) = rhd_gamma
196 call mpi_file_write(fh, values, n_par, mpi_double_precision, st, er)
197 call mpi_file_write(fh, names, n_par * name_len, mpi_character, st, er)
198 end subroutine rhd_write_info
199
200 !> Initialize the module
201 subroutine rhd_phys_init()
205 use mod_dust, only: dust_init
207 use mod_gravity, only: gravity_init
210 use mod_fld
211 use mod_afld
216
217 integer :: itr, idir
218
219 call rhd_read_params(par_files)
220
221 physics_type = "rhd"
222 phys_energy = rhd_energy
223 phys_total_energy = rhd_energy
224 phys_gamma = rhd_gamma
225
227 if(phys_trac) then
228 if(ndim .eq. 1) then
229 if(rhd_trac_type .gt. 2) then
231 if(mype==0) write(*,*) 'WARNING: set rhd_trac_type=1'
232 end if
234 else
235 phys_trac=.false.
236 if(mype==0) write(*,*) 'WARNING: set rhd_trac=F when ndim>=2'
237 end if
238 end if
239
240 ! set default gamma for polytropic/isothermal process
241 if(.not.rhd_energy) then
244 if(mype==0) write(*,*) 'WARNING: set rhd_thermal_conduction=F when rhd_energy=F'
245 end if
246 if(rhd_radiative_cooling) then
248 if(mype==0) write(*,*) 'WARNING: set rhd_radiative_cooling=F when rhd_energy=F'
249 end if
250 end if
251 if(.not.eq_state_units) then
254 if(mype==0) write(*,*) 'WARNING: set rhd_partial_ionization=F when eq_state_units=F'
255 end if
256 end if
258
259 allocate(start_indices(number_species),stop_indices(number_species))
260
261 ! set the index of the first flux variable for species 1
262 start_indices(1)=1
263 ! Determine flux variables
264 rho_ = var_set_rho()
265
266 allocate(mom(ndir))
267 mom(:) = var_set_momentum(ndir)
268
269 ! Set index of energy variable
270 if (rhd_energy) then
271 e_ = var_set_energy()
272 p_ = e_
273 else
274 e_ = -1
275 p_ = -1
276 end if
277
278 !> set radiation energy
279 r_e = var_set_radiation_energy()
280
281 ! set temperature as an auxiliary variable to get ionization degree
283 te_ = var_set_auxvar('Te','Te')
284 else
285 te_ = -1
286 end if
287
288 phys_get_dt => rhd_get_dt
289 phys_get_cmax => rhd_get_cmax
290 phys_get_a2max => rhd_get_a2max
291 phys_get_tcutoff => rhd_get_tcutoff
292 phys_get_cbounds => rhd_get_cbounds
293 phys_get_flux => rhd_get_flux
294 phys_add_source_geom => rhd_add_source_geom
295 phys_add_source => rhd_add_source
296 phys_to_conserved => rhd_to_conserved
297 phys_to_primitive => rhd_to_primitive
298 ! phys_ei_to_e => rhd_ei_to_e
299 ! phys_e_to_ei => rhd_e_to_ei
300 phys_check_params => rhd_check_params
301 phys_check_w => rhd_check_w
302 phys_get_pthermal => rhd_get_pthermal
303 phys_write_info => rhd_write_info
304 phys_handle_small_values => rhd_handle_small_values
305 phys_set_mg_bounds => rhd_set_mg_bounds
306 phys_get_trad => rhd_get_trad
307 phys_get_tgas => rhd_get_tgas
308
309 ! derive units from basic units
310 call rhd_physical_units()
311
312 if (rhd_dust) then
313 call dust_init(rho_, mom(:), e_)
314 endif
315
316 !> Initiate radiation-closure module
317 select case (rhd_radiation_formalism)
318 case('fld')
320 case('afld')
321 call afld_init(he_abundance, rhd_radiation_diffusion, rhd_gamma)
322 case default
323 call mpistop('Radiation formalism unknown')
324 end select
325
326 allocate(tracer(rhd_n_tracer))
327
328 ! Set starting index of tracers
329 do itr = 1, rhd_n_tracer
330 tracer(itr) = var_set_fluxvar("trc", "trp", itr, need_bc=.false.)
331 end do
332
333 ! set number of variables which need update ghostcells
334 nwgc=nwflux+nwaux
335
336 ! set the index of the last flux variable for species 1
337 stop_indices(1)=nwflux
338
339 ! Number of variables need reconstruction in w
340 nw_recon=nwflux
341
342 if(rhd_trac) then
343 tcoff_ = var_set_wextra()
344 iw_tcoff=tcoff_
345 else
346 tcoff_ = -1
347 end if
348
349 ! choose Rfactor in ideal gas law
351 rhd_get_rfactor=>rfactor_from_temperature_ionization
352 phys_update_temperature => rhd_update_temperature
353 else if(associated(usr_rfactor)) then
354 rhd_get_rfactor=>usr_rfactor
355 else
356 rhd_get_rfactor=>rfactor_from_constant_ionization
357 end if
358
359 ! initialize thermal conduction module
360 if (rhd_thermal_conduction) then
361 if (.not. rhd_energy) &
362 call mpistop("thermal conduction needs rhd_energy=T")
363
364 call sts_init()
366
367 allocate(tc_fl)
368 call tc_get_hd_params(tc_fl,tc_params_read_rhd)
369 tc_fl%get_temperature_from_conserved => rhd_get_temperature_from_etot
370 call add_sts_method(rhd_get_tc_dt_rhd,rhd_sts_set_source_tc_rhd,e_,1,e_,1,.false.)
371 call set_conversion_methods_to_head(rhd_e_to_ei, rhd_ei_to_e)
372 call set_error_handling_to_head(rhd_tc_handle_small_e)
373 tc_fl%get_temperature_from_eint => rhd_get_temperature_from_eint
374 tc_fl%get_rho => rhd_get_rho
375 tc_fl%e_ = e_
376 end if
377
378 ! Initialize radiative cooling module
379 if (rhd_radiative_cooling) then
380 if (.not. rhd_energy) &
381 call mpistop("radiative cooling needs rhd_energy=T")
383 allocate(rc_fl)
384 call radiative_cooling_init(rc_fl,rc_params_read)
385 rc_fl%get_rho => rhd_get_rho
386 rc_fl%get_pthermal => rhd_get_pthermal
387 rc_fl%get_var_Rfactor => rhd_get_rfactor
388 rc_fl%e_ = e_
389 rc_fl%Tcoff_ = tcoff_
390 end if
391 allocate(te_fl_rhd)
392 te_fl_rhd%get_rho=> rhd_get_rho
393 te_fl_rhd%get_pthermal=> rhd_get_pthermal
394 te_fl_rhd%get_var_Rfactor => rhd_get_rfactor
395{^ifthreed
396 phys_te_images => rhd_te_images
397}
398 ! Initialize viscosity module
399 if (rhd_viscosity) call viscosity_init(phys_wider_stencil)
400
401 ! Initialize gravity module
402 if (rhd_gravity) call gravity_init()
403
404 ! Initialize rotating_frame module
406
407 ! Initialize particles module
408 if (rhd_particles) then
409 call particles_init()
410 end if
411
412 ! Check whether custom flux types have been defined
413 if (.not. allocated(flux_type)) then
414 allocate(flux_type(ndir, nw))
415 flux_type = flux_default
416 else if (any(shape(flux_type) /= [ndir, nw])) then
417 call mpistop("phys_check error: flux_type has wrong shape")
418 end if
419
420 nvector = 1 ! No. vector vars
421 allocate(iw_vector(nvector))
422 iw_vector(1) = mom(1) - 1
423
424 !> Usefull constante
425 kbmpmua4 = unit_pressure**(-3.d0/4.d0)*unit_density*const_kb/(const_mp*fld_mu)*const_rad_a**(-1.d0/4.d0)
426 ! initialize ionization degree table
428 end subroutine rhd_phys_init
429
430{^ifthreed
431 subroutine rhd_te_images()
434 select case(convert_type)
435 case('EIvtiCCmpi','EIvtuCCmpi')
437 case('ESvtiCCmpi','ESvtuCCmpi')
439 case('SIvtiCCmpi','SIvtuCCmpi')
441 case default
442 call mpistop("Error in synthesize emission: Unknown convert_type")
443 end select
444 end subroutine rhd_te_images
445}
446!!start th cond
447 ! wrappers for STS functions in thermal_conductivity module
448 ! which take as argument the tc_fluid (defined in the physics module)
449 subroutine rhd_sts_set_source_tc_rhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
453 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
454 double precision, intent(in) :: x(ixi^s,1:ndim)
455 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
456 double precision, intent(in) :: my_dt
457 logical, intent(in) :: fix_conserve_at_step
458 call sts_set_source_tc_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
459 end subroutine rhd_sts_set_source_tc_rhd
460
461
462 function rhd_get_tc_dt_rhd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
463 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
464 !where tc_k_para_i=tc_k_para*B_i**2/B**2
465 !and T=p/rho
468
469 integer, intent(in) :: ixi^l, ixo^l
470 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
471 double precision, intent(in) :: w(ixi^s,1:nw)
472 double precision :: dtnew
473
474 dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
475 end function rhd_get_tc_dt_rhd
476
477
478 subroutine rhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
479 ! move this in a different routine as in mhd if needed in more places
482
483 integer, intent(in) :: ixi^l,ixo^l
484 double precision, intent(inout) :: w(ixi^s,1:nw)
485 double precision, intent(in) :: x(ixi^s,1:ndim)
486 integer, intent(in) :: step
487
488 integer :: idir
489 logical :: flag(ixi^s,1:nw)
490 character(len=140) :: error_msg
491
492 flag=.false.
493 where(w(ixo^s,e_)<small_e) flag(ixo^s,e_)=.true.
494 if(any(flag(ixo^s,e_))) then
495 select case (small_values_method)
496 case ("replace")
497 where(flag(ixo^s,e_)) w(ixo^s,e_)=small_e
498 case ("average")
499 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
500 case default
501 ! small values error shows primitive variables
502 w(ixo^s,e_)=w(ixo^s,e_)*(rhd_gamma - 1.0d0)
503 do idir = 1, ndir
504 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,rho_)
505 end do
506 write(error_msg,"(a,i3)") "Thermal conduction step ", step
507 call small_values_error(w, x, ixi^l, ixo^l, flag, error_msg)
508 end select
509 end if
510 end subroutine rhd_tc_handle_small_e
511
512 ! fill in tc_fluid fields from namelist
513 subroutine tc_params_read_rhd(fl)
516 type(tc_fluid), intent(inout) :: fl
517 integer :: n
518 logical :: tc_saturate=.false.
519 double precision :: tc_k_para=0d0
520
521 namelist /tc_list/ tc_saturate, tc_k_para
522
523 do n = 1, size(par_files)
524 open(unitpar, file=trim(par_files(n)), status="old")
525 read(unitpar, tc_list, end=111)
526111 close(unitpar)
527 end do
528 fl%tc_saturate = tc_saturate
529 fl%tc_k_para = tc_k_para
530
531 end subroutine tc_params_read_rhd
532
533 subroutine rhd_get_rho(w,x,ixI^L,ixO^L,rho)
535 integer, intent(in) :: ixi^l, ixo^l
536 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
537 double precision, intent(out) :: rho(ixi^s)
538
539 rho(ixo^s) = w(ixo^s,rho_)
540
541 end subroutine rhd_get_rho
542
543!!end th cond
544!!rad cool
545 subroutine rc_params_read(fl)
547 use mod_constants, only: bigdouble
548 use mod_basic_types, only: std_len
549 type(rc_fluid), intent(inout) :: fl
550 integer :: n
551 ! list parameters
552 integer :: ncool = 4000
553 double precision :: cfrac=0.1d0
554
555 !> Name of cooling curve
556 character(len=std_len) :: coolcurve='JCcorona'
557
558 !> Name of cooling method
559 character(len=std_len) :: coolmethod='exact'
560
561 !> Fixed temperature not lower than tlow
562 logical :: tfix=.false.
563
564 !> Lower limit of temperature
565 double precision :: tlow=bigdouble
566
567 !> Add cooling source in a split way (.true.) or un-split way (.false.)
568 logical :: rc_split=.false.
569
570
571 namelist /rc_list/ coolcurve, coolmethod, ncool, cfrac, tlow, tfix, rc_split
572
573 do n = 1, size(par_files)
574 open(unitpar, file=trim(par_files(n)), status="old")
575 read(unitpar, rc_list, end=111)
576111 close(unitpar)
577 end do
578
579 fl%ncool=ncool
580 fl%coolcurve=coolcurve
581 fl%coolmethod=coolmethod
582 fl%tlow=tlow
583 fl%Tfix=tfix
584 fl%rc_split=rc_split
585 fl%cfrac=cfrac
586 end subroutine rc_params_read
587!! end rad cool
588
591 use mod_dust, only: dust_check_params
592
593 if (.not. rhd_energy .and. rhd_pressure == 'adiabatic') then
594 if (rhd_gamma <= 0.0d0) call mpistop ("Error: rhd_gamma <= 0")
595 if (rhd_adiab < 0.0d0) call mpistop ("Error: rhd_adiab < 0")
597 elseif (rhd_pressure == 'Tcond') then
598 small_pressure = smalldouble
599 else
600 if (rhd_gamma <= 0.0d0 .or. rhd_gamma == 1.0d0) &
601 call mpistop ("Error: rhd_gamma <= 0 or rhd_gamma == 1.0")
602 small_e = small_pressure/(rhd_gamma - 1.0d0)
603 end if
604
606
607 if (rhd_dust) call dust_check_params()
608
609 ! if (rhd_radiation_diffusion .and. (.not. use_imex_scheme) ) &
610 if (rhd_radiation_diffusion .and. (.not. use_imex_scheme) .and. ((rhd_radiation_formalism .eq. 'fld') .or. (rhd_radiation_formalism .eq. 'afld')) ) &
611 call mpistop("Use an IMEX scheme when doing FLD")
612
614
615 end subroutine rhd_check_params
616
617 !> Set the boundaries for the diffusion of E
622
623 integer :: ib
624
625 ! Set boundary conditions for the multigrid solver
626 do ib = 1, 2*ndim
627 select case (typeboundary(r_e, ib))
628 case (bc_symm)
629 ! d/dx u = 0
630 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
631 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
632 case (bc_asymm)
633 ! u = 0
634 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
635 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
636 case (bc_cont)
637 ! d/dx u = 0
638 ! mg%bc(iB, mg_iphi)%bc_type = mg_bc_continuous
639 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
640 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
641 case (bc_periodic)
642 ! Nothing to do here
643 case (bc_noinflow)
644 call usr_special_mg_bc(ib)
645 case (bc_special)
646 call usr_special_mg_bc(ib)
647 case default
648 call mpistop("divE_multigrid warning: unknown b.c. ")
649 end select
650 end do
651 end subroutine rhd_set_mg_bounds
652
653 subroutine rhd_physical_units
655 double precision :: mp,kb
656 double precision :: a,b
657 ! Derive scaling units
658 if(si_unit) then
659 mp=mp_si
660 kb=kb_si
661 else
662 mp=mp_cgs
663 kb=kb_cgs
664 end if
665 if(eq_state_units) then
666 a = 1d0 + 4d0 * he_abundance
668 b = 2.d0+3.d0*he_abundance
669 else
670 b = 1d0 + h_ion_fr + he_abundance*(he_ion_fr*(he_ion_fr2 + 1d0)+1d0)
671 end if
672 rr = 1d0
673 else
674 a = 1d0
675 b = 1d0
676 rr = (1d0 + h_ion_fr + he_abundance*(he_ion_fr*(he_ion_fr2 + 1d0)+1d0))/(1d0 + 4d0 * he_abundance)
677 end if
678 if(unit_density/=1.d0) then
680 else
681 ! unit of numberdensity is independent by default
683 end if
684 if(unit_velocity/=1.d0) then
687 else if(unit_pressure/=1.d0) then
690 else
691 ! unit of temperature is independent by default
694 end if
695 if(unit_time/=1.d0) then
697 else
698 ! unit of length is independent by default
700 end if
701
702 !> Units for radiative flux and opacity
705
706 end subroutine rhd_physical_units
707
708 !> Returns logical argument flag where values are ok
709 subroutine rhd_check_w(primitive, ixI^L, ixO^L, w, flag)
711 use mod_dust, only: dust_check_w
712
713 logical, intent(in) :: primitive
714 integer, intent(in) :: ixi^l, ixo^l
715 double precision, intent(in) :: w(ixi^s, nw)
716 logical, intent(inout) :: flag(ixi^s,1:nw)
717 double precision :: tmp(ixi^s)
718
719 flag=.false.
720
721 if (rhd_energy) then
722 if (primitive) then
723 where(w(ixo^s, e_) < small_pressure) flag(ixo^s,e_) = .true.
724 else
725 tmp(ixo^s) = (rhd_gamma - 1.0d0)*(w(ixo^s, e_) - &
726 rhd_kin_en(w, ixi^l, ixo^l))
727 where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_) = .true.
728 endif
729 end if
730
731 where(w(ixo^s, rho_) < small_density) flag(ixo^s,rho_) = .true.
732
733 where(w(ixo^s, r_e) < small_r_e) flag(ixo^s,r_e) = .true.
734
735 if(rhd_dust) call dust_check_w(ixi^l,ixo^l,w,flag)
736
737 end subroutine rhd_check_w
738
739 !> Transform primitive variables into conservative ones
740 subroutine rhd_to_conserved(ixI^L, ixO^L, w, x)
742 use mod_dust, only: dust_to_conserved
743 integer, intent(in) :: ixi^l, ixo^l
744 double precision, intent(inout) :: w(ixi^s, nw)
745 double precision, intent(in) :: x(ixi^s, 1:ndim)
746 double precision :: invgam
747 integer :: idir, itr
748
749 !!if (fix_small_values) then
750 !! call rhd_handle_small_values(.true., w, x, ixI^L, ixO^L, 'rhd_to_conserved')
751 !!end if
752
753 if (rhd_energy) then
754 invgam = 1.d0/(rhd_gamma - 1.0d0)
755 ! Calculate total energy from pressure and kinetic energy
756 w(ixo^s, e_) = w(ixo^s, e_) * invgam + &
757 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * w(ixo^s, rho_)
758 end if
759
760 ! Convert velocity to momentum
761 do idir = 1, ndir
762 w(ixo^s, mom(idir)) = w(ixo^s, rho_) * w(ixo^s, mom(idir))
763 end do
764
765 if (rhd_dust) then
766 call dust_to_conserved(ixi^l, ixo^l, w, x)
767 end if
768
769 end subroutine rhd_to_conserved
770
771 !> Transform conservative variables into primitive ones
772 subroutine rhd_to_primitive(ixI^L, ixO^L, w, x)
774 use mod_dust, only: dust_to_primitive
775 integer, intent(in) :: ixi^l, ixo^l
776 double precision, intent(inout) :: w(ixi^s, nw)
777 double precision, intent(in) :: x(ixi^s, 1:ndim)
778 integer :: itr, idir
779 double precision :: inv_rho(ixo^s)
780
781 if (fix_small_values) then
782 call rhd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'rhd_to_primitive')
783 end if
784
785 inv_rho = 1.0d0 / w(ixo^s, rho_)
786
787 if (rhd_energy) then
788 ! Compute pressure
789 w(ixo^s, e_) = (rhd_gamma - 1.0d0) * (w(ixo^s, e_) - &
790 rhd_kin_en(w, ixi^l, ixo^l, inv_rho))
791 end if
792
793 ! Convert momentum to velocity
794 do idir = 1, ndir
795 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir)) * inv_rho
796 end do
797
798 ! Convert dust momentum to dust velocity
799 if (rhd_dust) then
800 call dust_to_primitive(ixi^l, ixo^l, w, x)
801 end if
802
803 end subroutine rhd_to_primitive
804
805 !> Transform internal energy to total energy
806 subroutine rhd_ei_to_e(ixI^L,ixO^L,w,x)
808 integer, intent(in) :: ixi^l, ixo^l
809 double precision, intent(inout) :: w(ixi^s, nw)
810 double precision, intent(in) :: x(ixi^s, 1:ndim)
811
812 ! Calculate total energy from internal and kinetic energy
813 w(ixo^s,e_)=w(ixo^s,e_)&
814 +rhd_kin_en(w,ixi^l,ixo^l)
815
816 end subroutine rhd_ei_to_e
817
818 !> Transform total energy to internal energy
819 subroutine rhd_e_to_ei(ixI^L,ixO^L,w,x)
821 integer, intent(in) :: ixi^l, ixo^l
822 double precision, intent(inout) :: w(ixi^s, nw)
823 double precision, intent(in) :: x(ixi^s, 1:ndim)
824
825 ! Calculate ei = e - ek
826 w(ixo^s,e_)=w(ixo^s,e_)&
827 -rhd_kin_en(w,ixi^l,ixo^l)
828
829 end subroutine rhd_e_to_ei
830
831 subroutine e_to_rhos(ixI^L, ixO^L, w, x)
833
834 integer, intent(in) :: ixi^l, ixo^l
835 double precision :: w(ixi^s, nw)
836 double precision, intent(in) :: x(ixi^s, 1:ndim)
837
838 if (rhd_energy) then
839 w(ixo^s, e_) = (rhd_gamma - 1.0d0) * w(ixo^s, rho_)**(1.0d0 - rhd_gamma) * &
840 (w(ixo^s, e_) - rhd_kin_en(w, ixi^l, ixo^l))
841 else
842 call mpistop("energy from entropy can not be used with -eos = iso !")
843 end if
844 end subroutine e_to_rhos
845
846 subroutine rhos_to_e(ixI^L, ixO^L, w, x)
848
849 integer, intent(in) :: ixi^l, ixo^l
850 double precision :: w(ixi^s, nw)
851 double precision, intent(in) :: x(ixi^s, 1:ndim)
852
853 if (rhd_energy) then
854 w(ixo^s, e_) = w(ixo^s, rho_)**(rhd_gamma - 1.0d0) * w(ixo^s, e_) &
855 / (rhd_gamma - 1.0d0) + rhd_kin_en(w, ixi^l, ixo^l)
856 else
857 call mpistop("entropy from energy can not be used with -eos = iso !")
858 end if
859 end subroutine rhos_to_e
860
861 !> Calculate v_i = m_i / rho within ixO^L
862 subroutine rhd_get_v(w, x, ixI^L, ixO^L, idim, v)
864 integer, intent(in) :: ixi^l, ixo^l, idim
865 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
866 double precision, intent(out) :: v(ixi^s)
867
868 v(ixo^s) = w(ixo^s, mom(idim)) / w(ixo^s, rho_)
869 end subroutine rhd_get_v
870
871 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
872 subroutine rhd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
876
877 integer, intent(in) :: ixi^l, ixo^l, idim
878 ! w in primitive form
879 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
880 double precision, intent(inout) :: cmax(ixi^s)
881
882 if(rhd_energy) then
883 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(rhd_gamma*w(ixo^s,p_)/w(ixo^s,rho_))
884 else
885 if (.not. associated(usr_set_pthermal)) then
886 select case (rhd_pressure)
887 case ('Trad')
888 cmax(ixo^s) = (w(ixo^s,r_e)*unit_pressure/const_rad_a)**0.25d0&
889 /unit_temperature*w(ixo^s, rho_)
890 case ('adiabatic')
891 cmax(ixo^s) = rhd_adiab * w(ixo^s, rho_)**rhd_gamma
892 case ('Tcond') !> Thermal conduction?!
893 cmax(ixo^s) = (rhd_gamma-1.d0)*w(ixo^s,r_e)
894 case default
895 call mpistop('rhd_pressure unknown, use Trad or adiabatic')
896 end select
897 else
898 call usr_set_pthermal(w,x,ixi^l,ixo^l,cmax)
899 end if
900 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(rhd_gamma*cmax(ixo^s)/w(ixo^s,rho_))
901 end if
902
903 if (rhd_dust) then
904 call dust_get_cmax_prim(w, x, ixi^l, ixo^l, idim, cmax)
905 end if
906 end subroutine rhd_get_cmax
907
908 subroutine rhd_get_a2max(w,x,ixI^L,ixO^L,a2max)
910
911 integer, intent(in) :: ixi^l, ixo^l
912 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
913 double precision, intent(inout) :: a2max(ndim)
914 double precision :: a2(ixi^s,ndim,nw)
915 integer :: gxo^l,hxo^l,jxo^l,kxo^l,i,j
916
917 a2=zero
918 do i = 1,ndim
919 !> 4th order
920 hxo^l=ixo^l-kr(i,^d);
921 gxo^l=hxo^l-kr(i,^d);
922 jxo^l=ixo^l+kr(i,^d);
923 kxo^l=jxo^l+kr(i,^d);
924 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
925 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
926 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/dxlevel(i)**2
927 end do
928 end subroutine rhd_get_a2max
929
930 !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
931 subroutine rhd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
933 integer, intent(in) :: ixi^l,ixo^l
934 double precision, intent(in) :: x(ixi^s,1:ndim)
935 double precision, intent(inout) :: w(ixi^s,1:nw)
936 double precision, intent(out) :: tco_local, tmax_local
937
938 double precision, parameter :: trac_delta=0.25d0
939 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
940 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
941 integer :: jxo^l,hxo^l
942 integer :: jxp^l,hxp^l,ixp^l
943 logical :: lrlt(ixi^s)
944
945 {^ifoned
946 call rhd_get_temperature_from_etot(w,x,ixi^l,ixi^l,te)
947
948 tco_local=zero
949 tmax_local=maxval(te(ixo^s))
950 select case(rhd_trac_type)
951 case(0)
952 w(ixi^s,tcoff_)=3.d5/unit_temperature
953 case(1)
954 hxo^l=ixo^l-1;
955 jxo^l=ixo^l+1;
956 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
957 lrlt=.false.
958 where(lts(ixo^s) > trac_delta)
959 lrlt(ixo^s)=.true.
960 end where
961 if(any(lrlt(ixo^s))) then
962 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
963 end if
964 case(2)
965 !> iijima et al. 2021, LTRAC method
966 ltrc=1.5d0
967 ltrp=2.5d0
968 ixp^l=ixo^l^ladd1;
969 hxo^l=ixo^l-1;
970 jxo^l=ixo^l+1;
971 hxp^l=ixp^l-1;
972 jxp^l=ixp^l+1;
973 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
974 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
975 w(ixo^s,tcoff_)=te(ixo^s)*&
976 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
977 case default
978 call mpistop("mrhd_trac_type not allowed for 1D simulation")
979 end select
980 }
981 end subroutine rhd_get_tcutoff
982
983 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
984 subroutine rhd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
986 use mod_dust, only: dust_get_cmax
987 use mod_variables
988
989 integer, intent(in) :: ixi^l, ixo^l, idim
990 ! conservative left and right status
991 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
992 ! primitive left and right status
993 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
994 double precision, intent(in) :: x(ixi^s, 1:ndim)
995 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
996 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
997 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
998
999 double precision :: wmean(ixi^s,nw)
1000 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1001 integer :: ix^d
1002
1003 select case(boundspeed)
1004 case (1)
1005 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1006 ! Methods for Fluid Dynamics" by Toro.
1007
1008 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1009 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1010 tmp3(ixo^s)=1.d0/(dsqrt(wlp(ixo^s,rho_))+dsqrt(wrp(ixo^s,rho_)))
1011 umean(ixo^s)=(wlp(ixo^s,mom(idim))*tmp1(ixo^s)+wrp(ixo^s,mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1012
1013 if(rhd_energy) then
1014 csoundl(ixo^s)=rhd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1015 csoundr(ixo^s)=rhd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1016 else
1017 select case (rhd_pressure)
1018 case ('Trad')
1019 csoundl(ixo^s)=rhd_gamma*kbmpmua4*wlp(ixo^s,r_e)**(1.d0/4)
1020 csoundr(ixo^s)=rhd_gamma*kbmpmua4*wrp(ixo^s,r_e)**(1.d0/4)
1021 case ('adiabatic')
1022 csoundl(ixo^s)=rhd_gamma*rhd_adiab*wlp(ixo^s,rho_)**(rhd_gamma-one)
1023 csoundr(ixo^s)=rhd_gamma*rhd_adiab*wrp(ixo^s,rho_)**(rhd_gamma-one)
1024 end select
1025 end if
1026
1027 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1028 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1029 (wrp(ixo^s,mom(idim))-wlp(ixo^s,mom(idim)))**2
1030
1031 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1032 if(present(cmin)) then
1033 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1034 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1035 if(h_correction) then
1036 {do ix^db=ixomin^db,ixomax^db\}
1037 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1038 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1039 {end do\}
1040 end if
1041 else
1042 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1043 end if
1044
1045 if (rhd_dust) then
1046 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1047 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1048 end if
1049
1050 case (2)
1051 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1052 tmp1(ixo^s)=wmean(ixo^s,mom(idim))/wmean(ixo^s,rho_)
1053 call rhd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1054 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1055
1056 if(present(cmin)) then
1057 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1058 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
1059 if(h_correction) then
1060 {do ix^db=ixomin^db,ixomax^db\}
1061 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1062 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1063 {end do\}
1064 end if
1065 else
1066 cmax(ixo^s,1)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1067 end if
1068
1069 if (rhd_dust) then
1070 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1071 end if
1072 case (3)
1073 ! Miyoshi 2005 JCP 208, 315 equation (67)
1074 if(rhd_energy) then
1075 csoundl(ixo^s)=rhd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1076 csoundr(ixo^s)=rhd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1077 else
1078 select case (rhd_pressure)
1079 case ('Trad')
1080 csoundl(ixo^s)=rhd_gamma*kbmpmua4*wlp(ixo^s,r_e)**(1.d0/4)
1081 csoundr(ixo^s)=rhd_gamma*kbmpmua4*wrp(ixo^s,r_e)**(1.d0/4)
1082 case ('adiabatic')
1083 csoundl(ixo^s)=rhd_gamma*rhd_adiab*wlp(ixo^s,rho_)**(rhd_gamma-one)
1084 csoundr(ixo^s)=rhd_gamma*rhd_adiab*wrp(ixo^s,rho_)**(rhd_gamma-one)
1085 end select
1086 end if
1087 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1088 if(present(cmin)) then
1089 cmin(ixo^s,1)=min(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))-csoundl(ixo^s)
1090 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1091 if(h_correction) then
1092 {do ix^db=ixomin^db,ixomax^db\}
1093 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1094 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1095 {end do\}
1096 end if
1097 else
1098 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1099 end if
1100 if (rhd_dust) then
1101 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1102 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1103 end if
1104 end select
1105
1106 end subroutine rhd_get_cbounds
1107
1108 !> Calculate the square of the thermal sound speed csound2 within ixO^L.
1109 !> csound2=gamma*p/rho
1110 subroutine rhd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1112 integer, intent(in) :: ixi^l, ixo^l
1113 double precision, intent(in) :: w(ixi^s,nw)
1114 double precision, intent(in) :: x(ixi^s,1:ndim)
1115 double precision, intent(out) :: csound2(ixi^s)
1116
1117 call rhd_get_ptot(w,x,ixi^l,ixo^l,csound2)
1118 csound2(ixo^s)=max(rhd_gamma,4.d0/3.d0)*csound2(ixo^s)/w(ixo^s,rho_)
1119 !> Turner & Stone 2001
1120
1121 end subroutine rhd_get_csound2
1122
1123 !> Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L
1124 subroutine rhd_get_pthermal(w, x, ixI^L, ixO^L, pth)
1128
1129 integer, intent(in) :: ixi^l, ixo^l
1130 double precision, intent(in) :: w(ixi^s, 1:nw)
1131 double precision, intent(in) :: x(ixi^s, 1:ndim)
1132 double precision, intent(out):: pth(ixi^s)
1133 integer :: iw, ix^d
1134
1135 if (rhd_energy) then
1136 pth(ixi^s) = (rhd_gamma - 1.d0) * (w(ixi^s, e_) - &
1137 0.5d0 * sum(w(ixi^s, mom(:))**2, dim=ndim+1) / w(ixi^s, rho_))
1138 else
1139 if (.not. associated(usr_set_pthermal)) then
1140 select case (rhd_pressure)
1141 case ('Trad')
1142 pth(ixi^s) = (w(ixi^s,r_e)*unit_pressure/const_rad_a)**0.25d0&
1143 /unit_temperature*w(ixi^s, rho_)
1144 case ('adiabatic')
1145 pth(ixi^s) = rhd_adiab * w(ixi^s, rho_)**rhd_gamma
1146 case ('Tcond') !> Thermal conduction?!
1147 pth(ixi^s) = (rhd_gamma-1.d0)*w(ixi^s,r_e)
1148 case default
1149 call mpistop('rhd_pressure unknown, use Trad or adiabatic')
1150 end select
1151 else
1152 call usr_set_pthermal(w,x,ixi^l,ixo^l,pth)
1153 end if
1154 end if
1155
1156 if (fix_small_values) then
1157 {do ix^db= ixo^lim^db\}
1158 if(pth(ix^d)<small_pressure) then
1159 pth(ix^d)=small_pressure
1160 endif
1161 {enddo^d&\}
1162 else if (check_small_values) then
1163 {do ix^db= ixo^lim^db\}
1164 if(pth(ix^d)<small_pressure) then
1165 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1166 " encountered when call rhd_get_pthermal"
1167 write(*,*) "Iteration: ", it, " Time: ", global_time
1168 write(*,*) "Location: ", x(ix^d,:)
1169 write(*,*) "Cell number: ", ix^d
1170 do iw=1,nw
1171 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1172 end do
1173 ! use erroneous arithmetic operation to crash the run
1174 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1175 write(*,*) "Saving status at the previous time step"
1176 crash=.true.
1177 end if
1178 {enddo^d&\}
1179 end if
1180
1181 end subroutine rhd_get_pthermal
1182
1183 !> Calculate radiation pressure within ixO^L
1184 subroutine rhd_get_pradiation(w, x, ixI^L, ixO^L, prad)
1186 use mod_fld
1187 use mod_afld
1188
1189 integer, intent(in) :: ixi^l, ixo^l
1190 double precision, intent(in) :: w(ixi^s, 1:nw)
1191 double precision, intent(in) :: x(ixi^s, 1:ndim)
1192 double precision, intent(out):: prad(ixo^s, 1:ndim, 1:ndim)
1193
1194 select case (rhd_radiation_formalism)
1195 case('fld')
1196 call fld_get_radpress(w, x, ixi^l, ixo^l, prad, nghostcells)
1197 case('afld')
1198 call afld_get_radpress(w, x, ixi^l, ixo^l, prad, nghostcells)
1199 case default
1200 call mpistop('Radiation formalism unknown')
1201 end select
1202 end subroutine rhd_get_pradiation
1203
1204 !> calculates the sum of the gas pressure and max Prad tensor element
1205 subroutine rhd_get_ptot(w, x, ixI^L, ixO^L, ptot)
1207
1208 integer, intent(in) :: ixi^l, ixo^l
1209 double precision, intent(in) :: w(ixi^s, 1:nw)
1210 double precision, intent(in) :: x(ixi^s, 1:ndim)
1211 double precision :: pth(ixi^s)
1212 double precision :: prad_tensor(ixo^s, 1:ndim, 1:ndim)
1213 double precision :: prad_max(ixo^s)
1214 double precision, intent(out):: ptot(ixi^s)
1215 integer :: ix^d
1216
1217 call rhd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1218 call rhd_get_pradiation(w, x, ixi^l, ixo^l, prad_tensor)
1219
1220 {do ix^d = ixomin^d,ixomax^d\}
1221 prad_max(ix^d) = maxval(prad_tensor(ix^d,:,:))
1222 {enddo\}
1223
1224 !> filter cmax
1225 if (radio_acoustic_filter) then
1226 call rhd_radio_acoustic_filter(x, ixi^l, ixo^l, prad_max)
1227 endif
1228
1229 ptot(ixo^s) = pth(ixo^s) + prad_max(ixo^s)
1230
1231 end subroutine rhd_get_ptot
1232
1233 !> Filter peaks in cmax due to radiation energy density, used for debugging
1234 subroutine rhd_radio_acoustic_filter(x, ixI^L, ixO^L, prad_max)
1236
1237 integer, intent(in) :: ixi^l, ixo^l
1238 double precision, intent(in) :: x(ixi^s, 1:ndim)
1239 double precision, intent(inout) :: prad_max(ixo^s)
1240
1241 double precision :: tmp_prad(ixi^s)
1242 integer :: ix^d, filter, idim
1243
1244 if (size_ra_filter .lt. 1) call mpistop("ra filter of size < 1 makes no sense")
1245 if (size_ra_filter .gt. nghostcells) call mpistop("ra filter of size < nghostcells makes no sense")
1246
1247 tmp_prad(ixi^s) = zero
1248 tmp_prad(ixo^s) = prad_max(ixo^s)
1249
1250 do filter = 1,size_ra_filter
1251 do idim = 1,ndim
1252 ! {do ix^D = ixOmin^D+filter,ixOmax^D-filter\}
1253 {do ix^d = ixomin^d,ixomax^d\}
1254 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d+filter*kr(idim,^d)))
1255 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d-filter*kr(idim,^d)))
1256 {enddo\}
1257 enddo
1258 enddo
1259 end subroutine rhd_radio_acoustic_filter
1260
1261 !> Calculate temperature=p/rho when in e_ the total energy is stored
1262 subroutine rhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1264 integer, intent(in) :: ixi^l, ixo^l
1265 double precision, intent(in) :: w(ixi^s, 1:nw)
1266 double precision, intent(in) :: x(ixi^s, 1:ndim)
1267 double precision, intent(out):: res(ixi^s)
1268
1269 double precision :: r(ixi^s)
1270
1271 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1272 call rhd_get_pthermal(w, x, ixi^l, ixo^l, res)
1273 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1274 end subroutine rhd_get_temperature_from_etot
1275
1276
1277 !> Calculate temperature=p/rho when in e_ the internal energy is stored
1278 subroutine rhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1280 integer, intent(in) :: ixi^l, ixo^l
1281 double precision, intent(in) :: w(ixi^s, 1:nw)
1282 double precision, intent(in) :: x(ixi^s, 1:ndim)
1283 double precision, intent(out):: res(ixi^s)
1284 double precision :: r(ixi^s)
1285
1286 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1287 res(ixo^s) = (rhd_gamma - 1.0d0) * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1288 end subroutine rhd_get_temperature_from_eint
1289
1290 !> Calculates gas temperature
1291 subroutine rhd_get_tgas(w, x, ixI^L, ixO^L, tgas)
1293
1294 integer, intent(in) :: ixi^l, ixo^l
1295 double precision, intent(in) :: w(ixi^s, 1:nw)
1296 double precision, intent(in) :: x(ixi^s, 1:ndim)
1297 double precision :: pth(ixi^s)
1298 double precision, intent(out):: tgas(ixi^s)
1299
1300 call rhd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1301 tgas(ixi^s) = pth(ixi^s)/w(ixi^s,rho_)
1302
1303 end subroutine rhd_get_tgas
1304
1305 !> Calculates radiation temperature
1306 subroutine rhd_get_trad(w, x, ixI^L, ixO^L, trad)
1308 use mod_constants
1309
1310 integer, intent(in) :: ixi^l, ixo^l
1311 double precision, intent(in) :: w(ixi^s, 1:nw)
1312 double precision, intent(in) :: x(ixi^s, 1:ndim)
1313 double precision, intent(out):: trad(ixi^s)
1314
1315 trad(ixi^s) = (w(ixi^s,r_e)*unit_pressure&
1316 /const_rad_a)**(1.d0/4.d0)/unit_temperature
1317
1318 end subroutine rhd_get_trad
1319
1320
1321 !these are very similar to the subroutines without 1, used in mod_thermal_conductivity
1322 !but no check on whether energy variable is present
1323 subroutine rhd_ei_to_e1(ixI^L,ixO^L,w,x)
1325 integer, intent(in) :: ixi^l, ixo^l
1326 double precision, intent(inout) :: w(ixi^s, nw)
1327 double precision, intent(in) :: x(ixi^s, 1:ndim)
1328
1329 ! Calculate total energy from internal and kinetic energy
1330 w(ixo^s,e_)=w(ixo^s,e_)&
1331 +rhd_kin_en(w,ixi^l,ixo^l)
1332
1333 end subroutine rhd_ei_to_e1
1334
1335 !> Transform total energy to internal energy
1336 !but no check on whether energy variable is present
1337 subroutine rhd_e_to_ei1(ixI^L,ixO^L,w,x)
1339 integer, intent(in) :: ixi^l, ixo^l
1340 double precision, intent(inout) :: w(ixi^s, nw)
1341 double precision, intent(in) :: x(ixi^s, 1:ndim)
1342
1343 ! Calculate ei = e - ek
1344 w(ixo^s,e_)=w(ixo^s,e_)&
1345 -rhd_kin_en(w,ixi^l,ixo^l)
1346
1347 end subroutine rhd_e_to_ei1
1348
1349 ! Calculate flux f_idim[iw]
1350 subroutine rhd_get_flux_cons(w, x, ixI^L, ixO^L, idim, f)
1352 use mod_dust, only: dust_get_flux
1353
1354 integer, intent(in) :: ixi^l, ixo^l, idim
1355 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1356 double precision, intent(out) :: f(ixi^s, nwflux)
1357 double precision :: pth(ixi^s), v(ixi^s),frame_vel(ixi^s)
1358 integer :: idir, itr
1359
1360 call rhd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1361 call rhd_get_v(w, x, ixi^l, ixo^l, idim, v)
1362
1363 f(ixo^s, rho_) = v(ixo^s) * w(ixo^s, rho_)
1364
1365 ! Momentum flux is v_i*m_i, +p in direction idim
1366 do idir = 1, ndir
1367 f(ixo^s, mom(idir)) = v(ixo^s) * w(ixo^s, mom(idir))
1368 end do
1369
1370 f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1371
1372 if(rhd_energy) then
1373 ! Energy flux is v_i*(e + p)
1374 f(ixo^s, e_) = v(ixo^s) * (w(ixo^s, e_) + pth(ixo^s))
1375 end if
1376
1377 if (rhd_radiation_advection) then
1378 f(ixo^s, r_e) = v(ixo^s) * w(ixo^s,r_e)
1379 else
1380 f(ixo^s, r_e) = zero
1381 endif
1382
1383 do itr = 1, rhd_n_tracer
1384 f(ixo^s, tracer(itr)) = v(ixo^s) * w(ixo^s, tracer(itr))
1385 end do
1386
1387 ! Dust fluxes
1388 if (rhd_dust) then
1389 call dust_get_flux(w, x, ixi^l, ixo^l, idim, f)
1390 end if
1391
1392 end subroutine rhd_get_flux_cons
1393
1394 ! Calculate flux f_idim[iw]
1395 subroutine rhd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1397 use mod_dust, only: dust_get_flux_prim
1398 use mod_viscosity, only: visc_get_flux_prim ! viscInDiv
1399
1400 integer, intent(in) :: ixi^l, ixo^l, idim
1401 ! conservative w
1402 double precision, intent(in) :: wc(ixi^s, 1:nw)
1403 ! primitive w
1404 double precision, intent(in) :: w(ixi^s, 1:nw)
1405 double precision, intent(in) :: x(ixi^s, 1:ndim)
1406 double precision, intent(out) :: f(ixi^s, nwflux)
1407 double precision :: pth(ixi^s),frame_vel(ixi^s)
1408 integer :: idir, itr
1409
1410 if (rhd_energy) then
1411 pth(ixo^s) = w(ixo^s,p_)
1412 else
1413 call rhd_get_pthermal(wc, x, ixi^l, ixo^l, pth)
1414 end if
1415
1416 f(ixo^s, rho_) = w(ixo^s,mom(idim)) * w(ixo^s, rho_)
1417
1418 ! Momentum flux is v_i*m_i, +p in direction idim
1419 do idir = 1, ndir
1420 f(ixo^s, mom(idir)) = w(ixo^s,mom(idim)) * wc(ixo^s, mom(idir))
1421 end do
1422
1423 f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1424
1425 if(rhd_energy) then
1426 ! Energy flux is v_i*(e + p)
1427 f(ixo^s, e_) = w(ixo^s,mom(idim)) * (wc(ixo^s, e_) + w(ixo^s,p_))
1428 end if
1429
1430 if (rhd_radiation_advection) then
1431 f(ixo^s, r_e) = w(ixo^s,mom(idim)) * wc(ixo^s,r_e)
1432 else
1433 f(ixo^s, r_e) = zero
1434 endif
1435
1436 do itr = 1, rhd_n_tracer
1437 f(ixo^s, tracer(itr)) = w(ixo^s,mom(idim)) * w(ixo^s, tracer(itr))
1438 end do
1439
1440 ! Dust fluxes
1441 if (rhd_dust) then
1442 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1443 end if
1444
1445 ! Viscosity fluxes - viscInDiv
1446 if (rhd_viscosity) then
1447 call visc_get_flux_prim(w, x, ixi^l, ixo^l, idim, f, rhd_energy)
1448 endif
1449
1450 end subroutine rhd_get_flux
1451
1452 !> Add geometrical source terms to w
1453 !>
1454 !> Notice that the expressions of the geometrical terms depend only on ndir,
1455 !> not ndim. Eg, they are the same in 2.5D and in 3D, for any geometry.
1456 !>
1457 subroutine rhd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1460 use mod_viscosity, only: visc_add_source_geom ! viscInDiv
1463 use mod_geometry
1464 integer, intent(in) :: ixi^l, ixo^l
1465 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:ndim)
1466 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw), w(ixi^s, 1:nw)
1467 ! to change and to set as a parameter in the parfile once the possibility to
1468 ! solve the equations in an angular momentum conserving form has been
1469 ! implemented (change tvdlf.t eg)
1470 double precision :: pth(ixi^s), source(ixi^s), minrho
1471 integer :: iw,idir, h1x^l{^nooned, h2x^l}
1472 integer :: mr_,mphi_ ! Polar var. names
1473 integer :: irho, ifluid, n_fluids
1474 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1475
1476 if (rhd_dust) then
1477 n_fluids = 1 + dust_n_species
1478 else
1479 n_fluids = 1
1480 end if
1481
1482 select case (coordinate)
1483
1485 !the user provides the functions of exp_factor and del_exp_factor
1486 if(associated(usr_set_surface)) call usr_set_surface(ixi^l,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
1487 if(rhd_energy) then
1488 source(ixo^s) =wprim(ixo^s, p_)
1489 else
1490 if(.not. associated(usr_set_pthermal)) then
1491 source(ixo^s)=rhd_adiab * wprim(ixo^s, rho_)**rhd_gamma
1492 else
1493 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1494 end if
1495 end if
1496 source(ixo^s) = source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1497 w(ixo^s,mom(1)) = w(ixo^s,mom(1)) + qdt*source(ixo^s)
1498
1499 case (cylindrical)
1500 if ((rhd_radiation_formalism .eq. 'fld') .or. (rhd_radiation_formalism .eq. 'afld')) then
1501 call mpistop("Diffusion term not implemented yet with cylkindrical geometries")
1502 end if
1503
1504 do ifluid = 0, n_fluids-1
1505 ! s[mr]=(pthermal+mphi**2/rho)/radius
1506 if (ifluid == 0) then
1507 ! gas
1508 irho = rho_
1509 mr_ = mom(r_)
1510 if(phi_>0) mphi_ = mom(phi_)
1511 if(rhd_energy) then
1512 source(ixo^s) =wprim(ixo^s, p_)
1513 else
1514 if(.not. associated(usr_set_pthermal)) then
1515 source(ixo^s)=rhd_adiab * wprim(ixo^s, rho_)**rhd_gamma
1516 else
1517 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1518 end if
1519 end if
1520 minrho = 0.0d0
1521 else
1522 ! dust : no pressure
1523 irho = dust_rho(ifluid)
1524 mr_ = dust_mom(r_, ifluid)
1525 if(phi_>0) mphi_ = dust_mom(phi_, ifluid)
1526 source(ixi^s) = zero
1527 minrho = 0.0d0
1528 end if
1529 if (phi_ > 0) then
1530 source(ixo^s) = source(ixo^s) + wprim(ixo^s, mphi_)**2 * wprim(ixo^s, irho)
1531 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1532 ! s[mphi]=(-vphi*vr*rho)/radius
1533 source(ixo^s) = -wprim(ixo^s, mphi_) * wprim(ixo^s, mr_) * wprim(ixo^s, irho)
1534 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * source(ixo^s) / x(ixo^s, r_)
1535 else
1536 ! s[mr]=2pthermal/radius
1537 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1538 end if
1539 end do
1540 case (spherical)
1541
1542 if ((rhd_radiation_formalism .eq. 'fld') .or. (rhd_radiation_formalism .eq. 'afld')) then
1543 call mpistop("Diffusion term not implemented yet with spherical geometries")
1544 end if
1545
1546 if (rhd_dust) then
1547 call mpistop("Dust geom source terms not implemented yet with spherical geometries")
1548 end if
1549 mr_ = mom(r_)
1550 if(phi_>0) mphi_ = mom(phi_)
1551 h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
1552 ! s[mr]=((vtheta**2+vphi**2)*rho+2*p)/r
1553 if(rhd_energy) then
1554 pth(ixo^s) =wprim(ixo^s, p_)
1555 else
1556 if(.not. associated(usr_set_pthermal)) then
1557 pth(ixo^s)=rhd_adiab * wprim(ixo^s, rho_)**rhd_gamma
1558 else
1559 call usr_set_pthermal(wct,x,ixi^l,ixo^l,pth)
1560 end if
1561 end if
1562 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1563 *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
1564 /block%dvolume(ixo^s)
1565 do idir = 2, ndir
1566 source(ixo^s) = source(ixo^s) + wprim(ixo^s, mom(idir))**2 * wprim(ixo^s, rho_)
1567 end do
1568 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, 1)
1569
1570 {^nooned
1571 ! s[mtheta]=-(vr*vtheta*rho)/r+cot(theta)*(vphi**2*rho+p)/r
1572 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1573 * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
1574 / block%dvolume(ixo^s)
1575 if (ndir == 3) then
1576 source(ixo^s) = source(ixo^s) + (wprim(ixo^s, mom(3))**2 * wprim(ixo^s, rho_)) / tan(x(ixo^s, 2))
1577 end if
1578 source(ixo^s) = source(ixo^s) - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)
1579 w(ixo^s, mom(2)) = w(ixo^s, mom(2)) + qdt * source(ixo^s) / x(ixo^s, 1)
1580
1581 if (ndir == 3) then
1582 ! s[mphi]=-(vphi*vr/rho)/r-cot(theta)*(vtheta*vphi/rho)/r
1583 source(ixo^s) = -(wprim(ixo^s, mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)&
1584 - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mom(3))) * wprim(ixo^s, rho_) / tan(x(ixo^s, 2))
1585 w(ixo^s, mom(3)) = w(ixo^s, mom(3)) + qdt * source(ixo^s) / x(ixo^s, 1)
1586 end if
1587 }
1588 end select
1589
1590 if (rhd_viscosity) call visc_add_source_geom(qdt,ixi^l,ixo^l,wprim,w,x)
1591
1592 if (rhd_rotating_frame) then
1593 if (rhd_dust) then
1594 call mpistop("Rotating frame not implemented yet with dust")
1595 else
1596 call rotating_frame_add_source(qdt,dtfactor,ixi^l,ixo^l,wprim,w,x)
1597 end if
1598 end if
1599
1600 end subroutine rhd_add_source_geom
1601
1602 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
1603 subroutine rhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1608 use mod_usr_methods, only: usr_gravity
1610
1611 integer, intent(in) :: ixi^l, ixo^l
1612 double precision, intent(in) :: qdt,dtfactor
1613 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw),x(ixi^s, 1:ndim)
1614 double precision, intent(inout) :: w(ixi^s, 1:nw)
1615 logical, intent(in) :: qsourcesplit
1616 logical, intent(inout) :: active
1617
1618 double precision :: gravity_field(ixi^s, 1:ndim)
1619 integer :: idust, idim
1620
1621 if(rhd_dust) then
1622 call dust_add_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1623 end if
1624
1625 if(rhd_radiative_cooling) then
1626 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1627 qsourcesplit,active, rc_fl)
1628 end if
1629
1630 if(rhd_viscosity) then
1631 call viscosity_add_source(qdt,ixi^l,ixo^l,wct,w,x,&
1632 rhd_energy,qsourcesplit,active)
1633 end if
1634
1635 if (rhd_gravity) then
1636 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1637 rhd_energy,.false.,qsourcesplit,active)
1638
1639 if (rhd_dust .and. qsourcesplit .eqv. grav_split) then
1640 active = .true.
1641
1642 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1643 do idust = 1, dust_n_species
1644 do idim = 1, ndim
1645 w(ixo^s, dust_mom(idim, idust)) = w(ixo^s, dust_mom(idim, idust)) &
1646 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s, dust_rho(idust))
1647 end do
1648 end do
1649 end if
1650 end if
1651
1652 !> This is where the radiation force and heating/cooling are added/
1653 call rhd_add_radiation_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1654
1655 if(rhd_partial_ionization) then
1656 if(.not.qsourcesplit) then
1657 active = .true.
1658 call rhd_update_temperature(ixi^l,ixo^l,wct,w,x)
1659 end if
1660 end if
1661
1662 end subroutine rhd_add_source
1663
1664 subroutine rhd_add_radiation_source(qdt,ixI^L,ixO^L,wCT,w,x,qsourcesplit,active)
1665 use mod_constants
1667 use mod_usr_methods
1668 use mod_fld!, only: fld_get_diffcoef_central, get_fld_rad_force, get_fld_energy_interact
1669 use mod_afld!, only: afld_get_diffcoef_central, get_afld_rad_force, get_afld_energy_interact
1670
1671 integer, intent(in) :: ixi^l, ixo^l
1672 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
1673 double precision, intent(in) :: wct(ixi^s,1:nw)
1674 double precision, intent(inout) :: w(ixi^s,1:nw)
1675 logical, intent(in) :: qsourcesplit
1676 logical, intent(inout) :: active
1677 double precision :: cmax(ixi^s)
1678
1679 select case(rhd_radiation_formalism)
1680 case('fld')
1681 if(fld_diff_scheme .eq. 'mg') call fld_get_diffcoef_central(w, wct, x, ixi^l, ixo^l)
1682 !> radiation force
1683 if(rhd_radiation_force) call get_fld_rad_force(qdt,ixi^l,ixo^l,wct,w,x,&
1684 rhd_energy,qsourcesplit,active)
1685 call rhd_handle_small_values(.true., w, x, ixi^l, ixo^l, 'fld_e_interact')
1686 case('afld')
1687 if (fld_diff_scheme .eq. 'mg') call afld_get_diffcoef_central(w, wct, x, ixi^l, ixo^l)
1688 !> radiation force
1689 if (rhd_radiation_force) call get_afld_rad_force(qdt,ixi^l,ixo^l,wct,w,x,&
1690 rhd_energy,qsourcesplit,active)
1691 call rhd_handle_small_values(.true., w, x, ixi^l, ixo^l, 'fld_e_interact')
1692 !> photon tiring, heating and cooling
1693 if (rhd_energy) then
1694 if (rhd_energy_interact) call get_afld_energy_interact(qdt,ixi^l,ixo^l,wct,w,x,&
1695 rhd_energy,qsourcesplit,active)
1696 endif
1697 case default
1698 call mpistop('Radiation formalism unknown')
1699 end select
1700
1701 ! ! !> NOT necessary for calculation, just want to know the grid-dependent-timestep
1702 ! call rhd_get_cmax(w, x, ixI^L, ixO^L, 2, cmax)
1703 ! w(ixI^S,i_test) = cmax(ixI^S)
1704 end subroutine rhd_add_radiation_source
1705
1706 subroutine rhd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
1708 use mod_dust, only: dust_get_dt
1711 use mod_gravity, only: gravity_get_dt
1712 use mod_fld, only: fld_radforce_get_dt
1714
1715 integer, intent(in) :: ixi^l, ixo^l
1716 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
1717 double precision, intent(in) :: w(ixi^s, 1:nw)
1718 double precision, intent(inout) :: dtnew
1719
1720 dtnew = bigdouble
1721
1722 if (.not. dt_c) then
1723
1724 if(rhd_dust) then
1725 call dust_get_dt(w, ixi^l, ixo^l, dtnew, dx^d, x)
1726 end if
1727
1728 if(rhd_radiation_force) then
1729 select case(rhd_radiation_formalism)
1730 case('fld')
1731 call fld_radforce_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1732 case('afld')
1733 call afld_radforce_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1734 case default
1735 call mpistop('Radiation formalism unknown')
1736 end select
1737 endif
1738
1739 if(rhd_radiative_cooling) then
1740 call cooling_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x,rc_fl)
1741 end if
1742
1743 if(rhd_viscosity) then
1744 call viscosity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1745 end if
1746
1747 if(rhd_gravity) then
1748 call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1749 end if
1750 else
1751 {^ifoned dtnew = dx1*unit_velocity/const_c}
1752 {^nooned dtnew = min(dx^d*unit_velocity/const_c)}
1753 endif
1754
1755 end subroutine rhd_get_dt
1756
1757 function rhd_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
1758 use mod_global_parameters, only: nw, ndim
1759 integer, intent(in) :: ixi^l, ixo^l
1760 double precision, intent(in) :: w(ixi^s, nw)
1761 double precision :: ke(ixo^s)
1762 double precision, intent(in), optional :: inv_rho(ixo^s)
1763
1764 if (present(inv_rho)) then
1765 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * inv_rho
1766 else
1767 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) / w(ixo^s, rho_)
1768 end if
1769 end function rhd_kin_en
1770
1771 function rhd_inv_rho(w, ixI^L, ixO^L) result(inv_rho)
1772 use mod_global_parameters, only: nw, ndim
1773 integer, intent(in) :: ixi^l, ixo^l
1774 double precision, intent(in) :: w(ixi^s, nw)
1775 double precision :: inv_rho(ixo^s)
1776
1777 ! Can make this more robust
1778 inv_rho = 1.0d0 / w(ixo^s, rho_)
1779 end function rhd_inv_rho
1780
1781 subroutine rhd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1782 ! handles hydro (density,pressure,velocity) bootstrapping
1783 ! any negative dust density is flagged as well (and throws an error)
1784 ! small_values_method=replace also for dust
1788 logical, intent(in) :: primitive
1789 integer, intent(in) :: ixi^l,ixo^l
1790 double precision, intent(inout) :: w(ixi^s,1:nw)
1791 double precision, intent(in) :: x(ixi^s,1:ndim)
1792 character(len=*), intent(in) :: subname
1793
1794 integer :: n,idir
1795 logical :: flag(ixi^s,1:nw)
1796
1797 call rhd_check_w(primitive, ixi^l, ixo^l, w, flag)
1798
1799 if (any(flag)) then
1800 select case (small_values_method)
1801 case ("replace")
1802 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
1803 do idir = 1, ndir
1804 if(small_values_fix_iw(mom(idir))) then
1805 where(flag(ixo^s,rho_)) w(ixo^s, mom(idir)) = 0.0d0
1806 end if
1807 end do
1808
1809 if (small_values_fix_iw(r_e)) then
1810 where(flag(ixo^s,r_e)) w(ixo^s,r_e) = small_r_e
1811 end if
1812
1813 if(rhd_energy)then
1814 if(small_values_fix_iw(e_)) then
1815 if(primitive) then
1816 where(flag(ixo^s,rho_)) w(ixo^s, p_) = small_pressure
1817 else
1818 where(flag(ixo^s,rho_)) w(ixo^s, e_) = small_e + rhd_kin_en(w,ixi^l,ixo^l)
1819 endif
1820 end if
1821 endif
1822
1823 if(rhd_energy) then
1824 if(primitive) then
1825 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
1826 else
1827 where(flag(ixo^s,e_))
1828 ! Add kinetic energy
1829 w(ixo^s,e_) = small_e + rhd_kin_en(w,ixi^l,ixo^l)
1830 end where
1831 end if
1832 end if
1833
1834 if(rhd_dust)then
1835 do n=1,dust_n_species
1836 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1837 do idir = 1, ndir
1838 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1839 enddo
1840 enddo
1841 endif
1842 case ("average")
1843 if(primitive)then
1844 ! averaging for all primitive fields, including dust
1845 call small_values_average(ixi^l, ixo^l, w, x, flag)
1846 else
1847 ! do averaging of density
1848 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
1849 if(rhd_energy) then
1850 ! do averaging of pressure
1851 w(ixi^s,p_)=(rhd_gamma-1.d0)*(w(ixi^s,e_) &
1852 -0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_))
1853 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
1854 w(ixi^s,e_)=w(ixi^s,p_)/(rhd_gamma-1.d0) &
1855 +0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_)
1856 end if
1857 ! do averaging of radiation energy
1858 call small_values_average(ixi^l, ixo^l, w, x, flag, r_e)
1859 if(rhd_dust)then
1860 do n=1,dust_n_species
1861 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1862 do idir = 1, ndir
1863 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1864 enddo
1865 enddo
1866 endif
1867 endif
1868 case default
1869 if(.not.primitive) then
1870 !convert w to primitive
1871 ! Calculate pressure = (gamma-1) * (e-ek)
1872 if(rhd_energy) then
1873 w(ixo^s,p_)=(rhd_gamma-1.d0)*(w(ixo^s,e_)-rhd_kin_en(w,ixi^l,ixo^l))
1874 end if
1875 ! Convert gas momentum to velocity
1876 do idir = 1, ndir
1877 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir))/w(ixo^s,rho_)
1878 end do
1879 end if
1880 ! NOTE: dust entries may still have conserved values here
1881 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1882 end select
1883 end if
1884 end subroutine rhd_handle_small_values
1885
1886 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1889 integer, intent(in) :: ixi^l, ixo^l
1890 double precision, intent(in) :: w(ixi^s,1:nw)
1891 double precision, intent(in) :: x(ixi^s,1:ndim)
1892 double precision, intent(out):: rfactor(ixi^s)
1893
1894 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1895
1896 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1897 ! assume the first and second ionization of Helium have the same degree
1898 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
1899
1900 end subroutine rfactor_from_temperature_ionization
1901
1902 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1904 integer, intent(in) :: ixi^l, ixo^l
1905 double precision, intent(in) :: w(ixi^s,1:nw)
1906 double precision, intent(in) :: x(ixi^s,1:ndim)
1907 double precision, intent(out):: rfactor(ixi^s)
1908
1909 rfactor(ixo^s)=rr
1910
1911 end subroutine rfactor_from_constant_ionization
1912
1913 subroutine rhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1916
1917 integer, intent(in) :: ixi^l, ixo^l
1918 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1919 double precision, intent(inout) :: w(ixi^s,1:nw)
1920
1921 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1922
1923 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1924
1925 call rhd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1926
1927 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1928 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1929
1930 end subroutine rhd_update_temperature
1931
1932end module mod_rhd_phys
Calculate w(iw)=w(iw)+qdt*SOURCE[wCT,qtC,x] within ixO for all indices iw=iwmin......
Module for including anisotropic flux limited diffusion (AFLD)-approximation in Radiation-hydrodynami...
Definition mod_afld.t:8
subroutine afld_get_diffcoef_central(w, wct, x, ixil, ixol)
Calculates cell-centered diffusion coefficient to be used in multigrid.
Definition mod_afld.t:684
subroutine, public get_afld_rad_force(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 This subroutine handles th...
Definition mod_afld.t:141
subroutine, public afld_init(he_abundance, rhd_radiation_diffusion, afld_gamma)
Initialising FLD-module: Read opacities Initialise Multigrid adimensionalise kappa Add extra variable...
Definition mod_afld.t:93
subroutine, public afld_radforce_get_dt(w, ixil, ixol, dtnew, dxd, x)
Definition mod_afld.t:217
subroutine, public afld_get_radpress(w, x, ixil, ixol, rad_pressure, nth)
Calculate Radiation Pressure Returns Radiation Pressure as tensor.
Definition mod_afld.t:518
subroutine, public get_afld_energy_interact(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 This subroutine handles th...
Definition mod_afld.t:244
Module with basic data types used in amrvac.
integer, parameter std_len
Default length for strings.
subroutine, public mpistop(message)
Exit MPI-AMRVAC with an error message.
Module for physical and numeric constants.
double precision, parameter bigdouble
A very large real number.
double precision, parameter const_rad_a
Module for including dust species, which interact with the gas through a drag force.
Definition mod_dust.t:3
subroutine, public dust_add_source(qdt, ixil, ixol, wct, w, x, qsourcesplit, active)
w[iw]= w[iw]+qdt*S[wCT, x] where S is the source based on wCT within ixO
Definition mod_dust.t:530
subroutine, public dust_to_primitive(ixil, ixol, w, x)
Definition mod_dust.t:228
subroutine, public dust_get_dt(w, ixil, ixol, dtnew, dxd, x)
Get dt related to dust and gas stopping time (Laibe 2011)
Definition mod_dust.t:898
subroutine, public dust_get_flux(w, x, ixil, ixol, idim, f)
Definition mod_dust.t:252
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:1011
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:1035
subroutine, public dust_init(g_rho, g_mom, g_energy)
Definition mod_dust.t:95
Module for flux conservation near refinement boundaries.
Nicolas Moens Module for including flux limited diffusion (FLD)-approximation in Radiation-hydrodynam...
Definition mod_fld.t:9
double precision, public fld_mu
mean particle mass
Definition mod_fld.t:19
subroutine, public get_fld_rad_force(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 This subroutine handles th...
Definition mod_fld.t:146
subroutine, public fld_get_radpress(w, x, ixil, ixol, rad_pressure, nth)
Calculate Radiation Pressure Returns Radiation Pressure as tensor.
Definition mod_fld.t:456
subroutine, public fld_init(he_abundance, radiation_diffusion, energy_interact, r_gamma)
Initialising FLD-module: Read opacities Initialise Multigrid adimensionalise kappa Add extra variable...
Definition mod_fld.t:94
subroutine fld_get_diffcoef_central(w, wct, x, ixil, ixol)
Calculates cell-centered diffusion coefficient to be used in multigrid.
Definition mod_fld.t:719
character(len=8) fld_diff_scheme
Which method to solve diffusion part.
Definition mod_fld.t:35
subroutine, public fld_radforce_get_dt(w, ixil, ixol, dtnew, dxd, x)
Definition mod_fld.t:181
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.
integer, parameter bc_noinflow
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.
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.
integer, dimension(:, :), allocatable typeboundary
Array indicating the type of boundary condition per variable and per physical boundary.
double precision unit_numberdensity
Physical scaling factor for number density.
character(len=std_len) convert_type
Which format to use when converting.
double precision unit_pressure
Physical scaling factor for pressure.
integer, parameter ndim
Number of spatial dimensions for grid variables.
double precision unit_length
Physical scaling factor for length.
logical use_particles
Use particles module or not.
character(len=std_len), dimension(:), allocatable par_files
Which par files are used as input.
integer mype
The rank of the current MPI task.
double precision, dimension(:), allocatable, parameter d
integer ndir
Number of spatial dimensions (components) for vector variables.
integer, parameter bc_periodic
integer, parameter bc_special
boundary condition types
double precision unit_velocity
Physical scaling factor for velocity.
double precision unit_temperature
Physical scaling factor for temperature.
double precision unit_radflux
Physical scaling factor for radiation flux.
logical si_unit
Use SI units (.true.) or use cgs units (.false.)
double precision, dimension(:,:), allocatable dx
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.
double precision, dimension(^nd) dxlevel
store unstretched cell size of current level
logical use_multigrid
Use multigrid (only available in 2D and 3D)
integer r_
Indices for cylindrical coordinates FOR TESTS, negative value when not used:
integer boundspeed
bound (left/min and right.max) speed of Riemann fan
integer, parameter unitconvert
logical check_small_values
check and optionally fix unphysical small values (density, gas pressure)
Module for including gravity in (magneto)hydrodynamics simulations.
Definition mod_gravity.t:2
logical grav_split
source split or not
Definition mod_gravity.t:6
subroutine gravity_get_dt(w, ixil, ixol, dtnew, dxd, x)
Definition mod_gravity.t:87
subroutine gravity_init()
Initialize the module.
Definition mod_gravity.t:26
subroutine gravity_add_source(qdt, ixil, ixol, wct, wctprim, w, x, energy, rhov, qsourcesplit, active)
w[iw]=w[iw]+qdt*S[wCT,qtC,x] where S is the source based on wCT within ixO
Definition mod_gravity.t:43
module ionization degree - get ionization degree for given temperature
subroutine ionization_degree_from_temperature(ixil, ixol, te, iz_h, iz_he)
Module to couple the octree-mg library to AMRVAC. This file uses the VACPP preprocessor,...
type(mg_t) mg
Data structure containing the multigrid tree.
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 for HD and MHD
subroutine radiative_cooling_init_params(phys_gamma, he_abund)
Radiative cooling initialization.
subroutine cooling_get_dt(w, ixil, ixol, dtnew, dxd, x, fl)
subroutine radiative_cooling_init(fl, read_params)
subroutine radiative_cooling_add_source(qdt, ixil, ixol, wct, wctprim, w, x, qsourcesplit, active, fl)
Radiation-Hydrodynamics physics module Module aims at solving the Hydrodynamic equations toghether wi...
logical, public, protected rhd_radiative_cooling
Whether radiative cooling is added.
subroutine, public rhd_to_primitive(ixil, ixol, w, x)
Transform conservative variables into primitive ones.
integer, public, protected rhd_n_tracer
Number of tracer species.
integer, public, protected p_
Index of the gas pressure (-1 if not present) should equal e_.
type(tc_fluid), allocatable, public tc_fl
logical, public, protected rhd_dust
Whether dust is added.
integer, public, protected rhd_trac_type
integer, public, protected te_
Indices of temperature.
logical, public, protected rhd_rotating_frame
Whether rotating frame is activated.
logical, public, protected rhd_viscosity
Whether viscosity is added.
double precision, public kbmpmua4
kb/(m_p mu)* 1/a_rad**4,
subroutine, public rhd_check_params
logical, public, protected rhd_trac
Whether TRAC method is used.
subroutine, public rhd_phys_init()
Initialize the module.
character(len=8), public rhd_radiation_formalism
Formalism to treat radiation.
double precision function, dimension(ixo^s), public rhd_kin_en(w, ixil, ixol, inv_rho)
integer, public, protected rho_
Index of the density (in the w array)
integer, public, protected r_e
Index of the radiation energy.
logical, public, protected rhd_radiation_diffusion
Treat radiation energy diffusion.
integer, dimension(:), allocatable, public, protected mom
Indices of the momentum density.
character(len=8), public rhd_pressure
In the case of no rhd_energy, how to compute pressure.
subroutine, public rhd_get_trad(w, x, ixil, ixol, trad)
Calculates radiation temperature.
subroutine, public rhd_check_w(primitive, ixil, ixol, w, flag)
Returns logical argument flag where values are ok.
logical, public, protected rhd_partial_ionization
Whether plasma is partially ionized.
logical, public, protected rhd_radiation_force
Treat radiation fld_Rad_force.
logical, public, protected rhd_energy
Whether an energy equation is used.
double precision, public rhd_adiab
The adiabatic constant.
subroutine, public rhd_get_tgas(w, x, ixil, ixol, tgas)
Calculates gas temperature.
double precision, public, protected small_r_e
The smallest allowed radiation energy.
subroutine, public rhd_get_ptot(w, x, ixil, ixol, ptot)
calculates the sum of the gas pressure and max Prad tensor element
double precision, public, protected h_ion_fr
Ionization fraction of H H_ion_fr = H+/(H+ + H)
integer, dimension(:), allocatable, public, protected tracer
Indices of the tracers.
subroutine, public rhd_to_conserved(ixil, ixol, w, x)
Transform primitive variables into conservative ones.
integer, public, protected e_
Index of the energy density (-1 if not present)
double precision, public, protected he_abundance
Helium abundance over Hydrogen.
subroutine, public rhd_set_mg_bounds
Set the boundaries for the diffusion of E.
subroutine, public rhd_get_csound2(w, x, ixil, ixol, csound2)
Calculate the square of the thermal sound speed csound2 within ixO^L. csound2=gamma*p/rho.
double precision, public, protected he_ion_fr2
Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)
logical, public, protected rhd_thermal_conduction
Whether thermal conduction is added.
double precision, public, protected rr
subroutine, public rhd_get_pthermal(w, x, ixil, ixol, pth)
Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.
logical, public, protected rhd_radiation_advection
Treat radiation advection.
logical, public, protected rhd_particles
Whether particles module is added.
double precision, public, protected he_ion_fr
Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)
type(te_fluid), allocatable, public te_fl_rhd
logical, public, protected rhd_energy_interact
Treat radiation-gas energy interaction.
logical, public, protected eq_state_units
type(rc_fluid), allocatable, public rc_fl
logical, public, protected rhd_gravity
Whether gravity is added.
integer, public, protected tcoff_
Index of the cutoff temperature for the TRAC method.
subroutine, public rhd_get_pradiation(w, x, ixil, ixol, prad)
Calculate radiation pressure within ixO^L.
double precision, public rhd_gamma
The adiabatic index.
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 1) in amrvac.par in sts_list set the following parameters which have...
subroutine, public add_sts_method(sts_getdt, sts_set_sources, startvar, nflux, startwbc, nwbc, evolve_b)
subroutine which added programatically a term to be calculated using STS Params: sts_getdt function c...
subroutine, public set_conversion_methods_to_head(sts_before_first_cycle, sts_after_last_cycle)
Set the hooks called before the first cycle and after the last cycle in the STS update This method sh...
subroutine, public set_error_handling_to_head(sts_error_handling)
Set the hook of error handling in the STS update. This method is called before updating the BC....
subroutine, public sts_init()
Initialize sts module.
Thermal conduction for HD and MHD 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)
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)
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(special_mg_bc), pointer usr_special_mg_bc
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 visc_get_flux_prim(w, x, ixil, ixol, idim, f, energy)
subroutine viscosity_add_source(qdt, ixil, ixol, wct, w, x, energy, qsourcesplit, active)
subroutine viscosity_init(phys_wider_stencil)
Initialize the module.
subroutine viscosity_get_dt(w, ixil, ixol, dtnew, dxd, x)
subroutine visc_add_source_geom(qdt, ixil, ixol, wct, w, x)