MPI-AMRVAC 3.2
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 if(rhd_trac) then
340 tcoff_ = var_set_wextra()
341 iw_tcoff=tcoff_
342 else
343 tcoff_ = -1
344 end if
345
346 ! choose Rfactor in ideal gas law
348 rhd_get_rfactor=>rfactor_from_temperature_ionization
349 phys_update_temperature => rhd_update_temperature
350 else if(associated(usr_rfactor)) then
351 rhd_get_rfactor=>usr_rfactor
352 else
353 rhd_get_rfactor=>rfactor_from_constant_ionization
354 end if
355
356 ! initialize thermal conduction module
357 if (rhd_thermal_conduction) then
358 if (.not. rhd_energy) &
359 call mpistop("thermal conduction needs rhd_energy=T")
360
361 call sts_init()
363
364 allocate(tc_fl)
365 call tc_get_hd_params(tc_fl,tc_params_read_rhd)
366 tc_fl%get_temperature_from_conserved => rhd_get_temperature_from_etot
367 call add_sts_method(rhd_get_tc_dt_rhd,rhd_sts_set_source_tc_rhd,e_,1,e_,1,.false.)
368 call set_conversion_methods_to_head(rhd_e_to_ei, rhd_ei_to_e)
369 call set_error_handling_to_head(rhd_tc_handle_small_e)
370 tc_fl%get_temperature_from_eint => rhd_get_temperature_from_eint
371 tc_fl%get_rho => rhd_get_rho
372 tc_fl%e_ = e_
373 end if
374
375 ! Initialize radiative cooling module
376 if (rhd_radiative_cooling) then
377 if (.not. rhd_energy) &
378 call mpistop("radiative cooling needs rhd_energy=T")
380 allocate(rc_fl)
381 call radiative_cooling_init(rc_fl,rc_params_read)
382 rc_fl%get_rho => rhd_get_rho
383 rc_fl%get_pthermal => rhd_get_pthermal
384 rc_fl%get_var_Rfactor => rhd_get_rfactor
385 rc_fl%e_ = e_
386 rc_fl%Tcoff_ = tcoff_
387 end if
388 allocate(te_fl_rhd)
389 te_fl_rhd%get_rho=> rhd_get_rho
390 te_fl_rhd%get_pthermal=> rhd_get_pthermal
391 te_fl_rhd%get_var_Rfactor => rhd_get_rfactor
392{^ifthreed
393 phys_te_images => rhd_te_images
394}
395 ! Initialize viscosity module
396 if (rhd_viscosity) call viscosity_init(phys_wider_stencil)
397
398 ! Initialize gravity module
399 if (rhd_gravity) call gravity_init()
400
401 ! Initialize rotating_frame module
403
404 ! Initialize particles module
405 if (rhd_particles) then
406 call particles_init()
407 end if
408
409 ! Check whether custom flux types have been defined
410 if (.not. allocated(flux_type)) then
411 allocate(flux_type(ndir, nw))
412 flux_type = flux_default
413 else if (any(shape(flux_type) /= [ndir, nw])) then
414 call mpistop("phys_check error: flux_type has wrong shape")
415 end if
416
417 nvector = 1 ! No. vector vars
418 allocate(iw_vector(nvector))
419 iw_vector(1) = mom(1) - 1
420
421 !> Usefull constante
422 kbmpmua4 = unit_pressure**(-3.d0/4.d0)*unit_density*const_kb/(const_mp*fld_mu)*const_rad_a**(-1.d0/4.d0)
423 ! initialize ionization degree table
425 end subroutine rhd_phys_init
426
427{^ifthreed
428 subroutine rhd_te_images()
431 select case(convert_type)
432 case('EIvtiCCmpi','EIvtuCCmpi')
434 case('ESvtiCCmpi','ESvtuCCmpi')
436 case('SIvtiCCmpi','SIvtuCCmpi')
438 case default
439 call mpistop("Error in synthesize emission: Unknown convert_type")
440 end select
441 end subroutine rhd_te_images
442}
443!!start th cond
444 ! wrappers for STS functions in thermal_conductivity module
445 ! which take as argument the tc_fluid (defined in the physics module)
446 subroutine rhd_sts_set_source_tc_rhd(ixI^L,ixO^L,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux)
450 integer, intent(in) :: ixi^l, ixo^l, igrid, nflux
451 double precision, intent(in) :: x(ixi^s,1:ndim)
452 double precision, intent(inout) :: wres(ixi^s,1:nw), w(ixi^s,1:nw)
453 double precision, intent(in) :: my_dt
454 logical, intent(in) :: fix_conserve_at_step
455 call sts_set_source_tc_hd(ixi^l,ixo^l,w,x,wres,fix_conserve_at_step,my_dt,igrid,nflux,tc_fl)
456 end subroutine rhd_sts_set_source_tc_rhd
457
458
459 function rhd_get_tc_dt_rhd(w,ixI^L,ixO^L,dx^D,x) result(dtnew)
460 !Check diffusion time limit dt < dx_i**2/((gamma-1)*tc_k_para_i/rho)
461 !where tc_k_para_i=tc_k_para*B_i**2/B**2
462 !and T=p/rho
465
466 integer, intent(in) :: ixi^l, ixo^l
467 double precision, intent(in) :: dx^d, x(ixi^s,1:ndim)
468 double precision, intent(in) :: w(ixi^s,1:nw)
469 double precision :: dtnew
470
471 dtnew=get_tc_dt_hd(w,ixi^l,ixo^l,dx^d,x,tc_fl)
472 end function rhd_get_tc_dt_rhd
473
474
475 subroutine rhd_tc_handle_small_e(w, x, ixI^L, ixO^L, step)
476 ! move this in a different routine as in mhd if needed in more places
479
480 integer, intent(in) :: ixi^l,ixo^l
481 double precision, intent(inout) :: w(ixi^s,1:nw)
482 double precision, intent(in) :: x(ixi^s,1:ndim)
483 integer, intent(in) :: step
484
485 integer :: idir
486 logical :: flag(ixi^s,1:nw)
487 character(len=140) :: error_msg
488
489 flag=.false.
490 where(w(ixo^s,e_)<small_e) flag(ixo^s,e_)=.true.
491 if(any(flag(ixo^s,e_))) then
492 select case (small_values_method)
493 case ("replace")
494 where(flag(ixo^s,e_)) w(ixo^s,e_)=small_e
495 case ("average")
496 call small_values_average(ixi^l, ixo^l, w, x, flag, e_)
497 case default
498 ! small values error shows primitive variables
499 w(ixo^s,e_)=w(ixo^s,e_)*(rhd_gamma - 1.0d0)
500 do idir = 1, ndir
501 w(ixo^s, iw_mom(idir)) = w(ixo^s, iw_mom(idir))/w(ixo^s,rho_)
502 end do
503 write(error_msg,"(a,i3)") "Thermal conduction step ", step
504 call small_values_error(w, x, ixi^l, ixo^l, flag, error_msg)
505 end select
506 end if
507 end subroutine rhd_tc_handle_small_e
508
509 ! fill in tc_fluid fields from namelist
510 subroutine tc_params_read_rhd(fl)
513 type(tc_fluid), intent(inout) :: fl
514 integer :: n
515 logical :: tc_saturate=.false.
516 double precision :: tc_k_para=0d0
517
518 namelist /tc_list/ tc_saturate, tc_k_para
519
520 do n = 1, size(par_files)
521 open(unitpar, file=trim(par_files(n)), status="old")
522 read(unitpar, tc_list, end=111)
523111 close(unitpar)
524 end do
525 fl%tc_saturate = tc_saturate
526 fl%tc_k_para = tc_k_para
527
528 end subroutine tc_params_read_rhd
529
530 subroutine rhd_get_rho(w,x,ixI^L,ixO^L,rho)
532 integer, intent(in) :: ixi^l, ixo^l
533 double precision, intent(in) :: w(ixi^s,1:nw),x(ixi^s,1:ndim)
534 double precision, intent(out) :: rho(ixi^s)
535
536 rho(ixo^s) = w(ixo^s,rho_)
537
538 end subroutine rhd_get_rho
539
540!!end th cond
541!!rad cool
542 subroutine rc_params_read(fl)
544 use mod_constants, only: bigdouble
545 use mod_basic_types, only: std_len
546 type(rc_fluid), intent(inout) :: fl
547 integer :: n
548 ! list parameters
549 integer :: ncool = 4000
550
551 !> Name of cooling curve
552 character(len=std_len) :: coolcurve='JCcorona'
553
554 !> Fixed temperature not lower than tlow
555 logical :: tfix=.false.
556
557 !> Lower limit of temperature
558 double precision :: tlow=bigdouble
559
560 !> Add cooling source in a split way (.true.) or un-split way (.false.)
561 logical :: rc_split=.false.
562
563
564 namelist /rc_list/ coolcurve, ncool, tlow, tfix, rc_split
565
566 do n = 1, size(par_files)
567 open(unitpar, file=trim(par_files(n)), status="old")
568 read(unitpar, rc_list, end=111)
569111 close(unitpar)
570 end do
571
572 fl%ncool=ncool
573 fl%coolcurve=coolcurve
574 fl%tlow=tlow
575 fl%Tfix=tfix
576 fl%rc_split=rc_split
577 end subroutine rc_params_read
578!! end rad cool
579
582 use mod_dust, only: dust_check_params
583
584 if (.not. rhd_energy .and. rhd_pressure == 'adiabatic') then
585 if (rhd_gamma <= 0.0d0) call mpistop ("Error: rhd_gamma <= 0")
586 if (rhd_adiab < 0.0d0) call mpistop ("Error: rhd_adiab < 0")
588 elseif (rhd_pressure == 'Tcond') then
589 small_pressure = smalldouble
590 else
591 if (rhd_gamma <= 0.0d0 .or. rhd_gamma == 1.0d0) &
592 call mpistop ("Error: rhd_gamma <= 0 or rhd_gamma == 1.0")
593 small_e = small_pressure/(rhd_gamma - 1.0d0)
594 end if
595
597
598 if (rhd_dust) call dust_check_params()
599
600 ! if (rhd_radiation_diffusion .and. (.not. use_imex_scheme) ) &
601 if (rhd_radiation_diffusion .and. (.not. use_imex_scheme) .and. ((rhd_radiation_formalism .eq. 'fld') .or. (rhd_radiation_formalism .eq. 'afld')) ) &
602 call mpistop("Use an IMEX scheme when doing FLD")
603
605
606 end subroutine rhd_check_params
607
608 !> Set the boundaries for the diffusion of E
613
614 integer :: ib
615
616 ! Set boundary conditions for the multigrid solver
617 do ib = 1, 2*ndim
618 select case (typeboundary(r_e, ib))
619 case (bc_symm)
620 ! d/dx u = 0
621 mg%bc(ib, mg_iphi)%bc_type = mg_bc_neumann
622 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
623 case (bc_asymm)
624 ! u = 0
625 mg%bc(ib, mg_iphi)%bc_type = mg_bc_dirichlet
626 mg%bc(ib, mg_iphi)%bc_value = 0.0_dp
627 case (bc_cont)
628 ! d/dx u = 0
629 ! mg%bc(iB, mg_iphi)%bc_type = mg_bc_continuous
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_periodic)
633 ! Nothing to do here
634 case (bc_noinflow)
635 call usr_special_mg_bc(ib)
636 case (bc_special)
637 call usr_special_mg_bc(ib)
638 case default
639 call mpistop("divE_multigrid warning: unknown b.c. ")
640 end select
641 end do
642 end subroutine rhd_set_mg_bounds
643
644 subroutine rhd_physical_units
646 double precision :: mp,kb
647 double precision :: a,b
648 ! Derive scaling units
649 if(si_unit) then
650 mp=mp_si
651 kb=kb_si
652 else
653 mp=mp_cgs
654 kb=kb_cgs
655 end if
656 if(eq_state_units) then
657 a=1d0+4d0*he_abundance
659 b=1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0)
660 else
661 b=2d0+3d0*he_abundance
662 end if
663 rr=1d0
664 else
665 a=1d0
666 b=1d0
667 rr=(1d0+h_ion_fr+he_abundance*(he_ion_fr*(he_ion_fr2+1d0)+1d0))/(1d0+4d0*he_abundance)
668 end if
669 if(unit_density/=1.d0 .or. unit_numberdensity/=1.d0) then
670 if(unit_density/=1.d0) then
672 else if(unit_numberdensity/=1.d0) then
674 end if
675 if(unit_temperature/=1.d0) then
678 if(unit_length/=1.d0) then
680 else if(unit_time/=1.d0) then
682 end if
683 else if(unit_pressure/=1.d0) then
686 if(unit_length/=1.d0) then
688 else if(unit_time/=1.d0) then
690 end if
691 else if(unit_velocity/=1.d0) then
694 if(unit_length/=1.d0) then
696 else if(unit_time/=1.d0) then
698 end if
699 else if(unit_time/=1.d0) then
703 end if
704 else if(unit_temperature/=1.d0) then
705 ! units of temperature and velocity are dependent
706 if(unit_pressure/=1.d0) then
710 if(unit_length/=1.d0) then
712 else if(unit_time/=1.d0) then
714 end if
715 end if
716 else if(unit_pressure/=1.d0) then
717 if(unit_velocity/=1.d0) then
721 if(unit_length/=1.d0) then
723 else if(unit_time/=1.d0) then
725 end if
726 else if(unit_time/=0.d0) then
731 end if
732 end if
733
734 !> Units for radiative flux and opacity
737
738 end subroutine rhd_physical_units
739
740 !> Returns logical argument flag where values are ok
741 subroutine rhd_check_w(primitive, ixI^L, ixO^L, w, flag)
743 use mod_dust, only: dust_check_w
744
745 logical, intent(in) :: primitive
746 integer, intent(in) :: ixi^l, ixo^l
747 double precision, intent(in) :: w(ixi^s, nw)
748 logical, intent(inout) :: flag(ixi^s,1:nw)
749 double precision :: tmp(ixi^s)
750
751 flag=.false.
752
753 if (rhd_energy) then
754 if (primitive) then
755 where(w(ixo^s, e_) < small_pressure) flag(ixo^s,e_) = .true.
756 else
757 tmp(ixo^s) = (rhd_gamma - 1.0d0)*(w(ixo^s, e_) - &
758 rhd_kin_en(w, ixi^l, ixo^l))
759 where(tmp(ixo^s) < small_pressure) flag(ixo^s,e_) = .true.
760 endif
761 end if
762
763 where(w(ixo^s, rho_) < small_density) flag(ixo^s,rho_) = .true.
764
765 where(w(ixo^s, r_e) < small_r_e) flag(ixo^s,r_e) = .true.
766
767 if(rhd_dust) call dust_check_w(ixi^l,ixo^l,w,flag)
768
769 end subroutine rhd_check_w
770
771 !> Transform primitive variables into conservative ones
772 subroutine rhd_to_conserved(ixI^L, ixO^L, w, x)
774 use mod_dust, only: dust_to_conserved
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 double precision :: invgam
779 integer :: idir, itr
780
781 !!if (fix_small_values) then
782 !! call rhd_handle_small_values(.true., w, x, ixI^L, ixO^L, 'rhd_to_conserved')
783 !!end if
784
785 if (rhd_energy) then
786 invgam = 1.d0/(rhd_gamma - 1.0d0)
787 ! Calculate total energy from pressure and kinetic energy
788 w(ixo^s, e_) = w(ixo^s, e_) * invgam + &
789 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * w(ixo^s, rho_)
790 end if
791
792 ! Convert velocity to momentum
793 do idir = 1, ndir
794 w(ixo^s, mom(idir)) = w(ixo^s, rho_) * w(ixo^s, mom(idir))
795 end do
796
797 if (rhd_dust) then
798 call dust_to_conserved(ixi^l, ixo^l, w, x)
799 end if
800
801 end subroutine rhd_to_conserved
802
803 !> Transform conservative variables into primitive ones
804 subroutine rhd_to_primitive(ixI^L, ixO^L, w, x)
806 use mod_dust, only: dust_to_primitive
807 integer, intent(in) :: ixi^l, ixo^l
808 double precision, intent(inout) :: w(ixi^s, nw)
809 double precision, intent(in) :: x(ixi^s, 1:ndim)
810 integer :: itr, idir
811 double precision :: inv_rho(ixo^s)
812
813 if (fix_small_values) then
814 call rhd_handle_small_values(.false., w, x, ixi^l, ixo^l, 'rhd_to_primitive')
815 end if
816
817 inv_rho = 1.0d0 / w(ixo^s, rho_)
818
819 if (rhd_energy) then
820 ! Compute pressure
821 w(ixo^s, e_) = (rhd_gamma - 1.0d0) * (w(ixo^s, e_) - &
822 rhd_kin_en(w, ixi^l, ixo^l, inv_rho))
823 end if
824
825 ! Convert momentum to velocity
826 do idir = 1, ndir
827 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir)) * inv_rho
828 end do
829
830 ! Convert dust momentum to dust velocity
831 if (rhd_dust) then
832 call dust_to_primitive(ixi^l, ixo^l, w, x)
833 end if
834
835 end subroutine rhd_to_primitive
836
837 !> Transform internal energy to total energy
838 subroutine rhd_ei_to_e(ixI^L,ixO^L,w,x)
840 integer, intent(in) :: ixi^l, ixo^l
841 double precision, intent(inout) :: w(ixi^s, nw)
842 double precision, intent(in) :: x(ixi^s, 1:ndim)
843
844 ! Calculate total energy from internal and kinetic energy
845 w(ixo^s,e_)=w(ixo^s,e_)&
846 +rhd_kin_en(w,ixi^l,ixo^l)
847
848 end subroutine rhd_ei_to_e
849
850 !> Transform total energy to internal energy
851 subroutine rhd_e_to_ei(ixI^L,ixO^L,w,x)
853 integer, intent(in) :: ixi^l, ixo^l
854 double precision, intent(inout) :: w(ixi^s, nw)
855 double precision, intent(in) :: x(ixi^s, 1:ndim)
856
857 ! Calculate ei = e - ek
858 w(ixo^s,e_)=w(ixo^s,e_)&
859 -rhd_kin_en(w,ixi^l,ixo^l)
860
861 end subroutine rhd_e_to_ei
862
863 subroutine e_to_rhos(ixI^L, ixO^L, w, x)
865
866 integer, intent(in) :: ixi^l, ixo^l
867 double precision :: w(ixi^s, nw)
868 double precision, intent(in) :: x(ixi^s, 1:ndim)
869
870 if (rhd_energy) then
871 w(ixo^s, e_) = (rhd_gamma - 1.0d0) * w(ixo^s, rho_)**(1.0d0 - rhd_gamma) * &
872 (w(ixo^s, e_) - rhd_kin_en(w, ixi^l, ixo^l))
873 else
874 call mpistop("energy from entropy can not be used with -eos = iso !")
875 end if
876 end subroutine e_to_rhos
877
878 subroutine rhos_to_e(ixI^L, ixO^L, w, x)
880
881 integer, intent(in) :: ixi^l, ixo^l
882 double precision :: w(ixi^s, nw)
883 double precision, intent(in) :: x(ixi^s, 1:ndim)
884
885 if (rhd_energy) then
886 w(ixo^s, e_) = w(ixo^s, rho_)**(rhd_gamma - 1.0d0) * w(ixo^s, e_) &
887 / (rhd_gamma - 1.0d0) + rhd_kin_en(w, ixi^l, ixo^l)
888 else
889 call mpistop("entropy from energy can not be used with -eos = iso !")
890 end if
891 end subroutine rhos_to_e
892
893 !> Calculate v_i = m_i / rho within ixO^L
894 subroutine rhd_get_v(w, x, ixI^L, ixO^L, idim, v)
896 integer, intent(in) :: ixi^l, ixo^l, idim
897 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
898 double precision, intent(out) :: v(ixi^s)
899
900 v(ixo^s) = w(ixo^s, mom(idim)) / w(ixo^s, rho_)
901 end subroutine rhd_get_v
902
903 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
904 subroutine rhd_get_cmax(w, x, ixI^L, ixO^L, idim, cmax)
908
909 integer, intent(in) :: ixi^l, ixo^l, idim
910 ! w in primitive form
911 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s, 1:ndim)
912 double precision, intent(inout) :: cmax(ixi^s)
913
914 if(rhd_energy) then
915 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(rhd_gamma*w(ixo^s,p_)/w(ixo^s,rho_))
916 else
917 if (.not. associated(usr_set_pthermal)) then
918 select case (rhd_pressure)
919 case ('Trad')
920 cmax(ixo^s) = (w(ixo^s,r_e)*unit_pressure/const_rad_a)**0.25d0&
921 /unit_temperature*w(ixo^s, rho_)
922 case ('adiabatic')
923 cmax(ixo^s) = rhd_adiab * w(ixo^s, rho_)**rhd_gamma
924 case ('Tcond') !> Thermal conduction?!
925 cmax(ixo^s) = (rhd_gamma-1.d0)*w(ixo^s,r_e)
926 case default
927 call mpistop('rhd_pressure unknown, use Trad or adiabatic')
928 end select
929 else
930 call usr_set_pthermal(w,x,ixi^l,ixo^l,cmax)
931 end if
932 cmax(ixo^s)=dabs(w(ixo^s,mom(idim)))+dsqrt(rhd_gamma*cmax(ixo^s)/w(ixo^s,rho_))
933 end if
934
935 if (rhd_dust) then
936 call dust_get_cmax_prim(w, x, ixi^l, ixo^l, idim, cmax)
937 end if
938 end subroutine rhd_get_cmax
939
940 subroutine rhd_get_a2max(w,x,ixI^L,ixO^L,a2max)
942
943 integer, intent(in) :: ixi^l, ixo^l
944 double precision, intent(in) :: w(ixi^s, nw), x(ixi^s,1:ndim)
945 double precision, intent(inout) :: a2max(ndim)
946 double precision :: a2(ixi^s,ndim,nw)
947 integer :: gxo^l,hxo^l,jxo^l,kxo^l,i,j
948
949 a2=zero
950 do i = 1,ndim
951 !> 4th order
952 hxo^l=ixo^l-kr(i,^d);
953 gxo^l=hxo^l-kr(i,^d);
954 jxo^l=ixo^l+kr(i,^d);
955 kxo^l=jxo^l+kr(i,^d);
956 a2(ixo^s,i,1:nw)=dabs(-w(kxo^s,1:nw)+16.d0*w(jxo^s,1:nw)&
957 -30.d0*w(ixo^s,1:nw)+16.d0*w(hxo^s,1:nw)-w(gxo^s,1:nw))
958 a2max(i)=maxval(a2(ixo^s,i,1:nw))/12.d0/dxlevel(i)**2
959 end do
960 end subroutine rhd_get_a2max
961
962 !> get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)
963 subroutine rhd_get_tcutoff(ixI^L,ixO^L,w,x,tco_local,Tmax_local)
965 integer, intent(in) :: ixi^l,ixo^l
966 double precision, intent(in) :: x(ixi^s,1:ndim)
967 double precision, intent(inout) :: w(ixi^s,1:nw)
968 double precision, intent(out) :: tco_local, tmax_local
969
970 double precision, parameter :: trac_delta=0.25d0
971 double precision :: tmp1(ixi^s),te(ixi^s),lts(ixi^s)
972 double precision :: ltr(ixi^s),ltrc,ltrp,tcoff(ixi^s)
973 integer :: jxo^l,hxo^l
974 integer :: jxp^l,hxp^l,ixp^l
975 logical :: lrlt(ixi^s)
976
977 {^ifoned
978 call rhd_get_temperature_from_etot(w,x,ixi^l,ixi^l,te)
979
980 tco_local=zero
981 tmax_local=maxval(te(ixo^s))
982 select case(rhd_trac_type)
983 case(0)
984 w(ixi^s,tcoff_)=3.d5/unit_temperature
985 case(1)
986 hxo^l=ixo^l-1;
987 jxo^l=ixo^l+1;
988 lts(ixo^s)=0.5d0*dabs(te(jxo^s)-te(hxo^s))/te(ixo^s)
989 lrlt=.false.
990 where(lts(ixo^s) > trac_delta)
991 lrlt(ixo^s)=.true.
992 end where
993 if(any(lrlt(ixo^s))) then
994 tco_local=maxval(te(ixo^s), mask=lrlt(ixo^s))
995 end if
996 case(2)
997 !> iijima et al. 2021, LTRAC method
998 ltrc=1.5d0
999 ltrp=2.5d0
1000 ixp^l=ixo^l^ladd1;
1001 hxo^l=ixo^l-1;
1002 jxo^l=ixo^l+1;
1003 hxp^l=ixp^l-1;
1004 jxp^l=ixp^l+1;
1005 lts(ixp^s)=0.5d0*abs(te(jxp^s)-te(hxp^s))/te(ixp^s)
1006 ltr(ixp^s)=max(one, (exp(lts(ixp^s))/ltrc)**ltrp)
1007 w(ixo^s,tcoff_)=te(ixo^s)*&
1008 (0.25*(ltr(jxo^s)+two*ltr(ixo^s)+ltr(hxo^s)))**0.4d0
1009 case default
1010 call mpistop("mrhd_trac_type not allowed for 1D simulation")
1011 end select
1012 }
1013 end subroutine rhd_get_tcutoff
1014
1015 !> Calculate cmax_idim = csound + abs(v_idim) within ixO^L
1016 subroutine rhd_get_cbounds(wLC, wRC, wLp, wRp, x, ixI^L, ixO^L, idim,Hspeed,cmax, cmin)
1018 use mod_dust, only: dust_get_cmax
1019 use mod_variables
1020
1021 integer, intent(in) :: ixi^l, ixo^l, idim
1022 ! conservative left and right status
1023 double precision, intent(in) :: wlc(ixi^s, nw), wrc(ixi^s, nw)
1024 ! primitive left and right status
1025 double precision, intent(in) :: wlp(ixi^s, nw), wrp(ixi^s, nw)
1026 double precision, intent(in) :: x(ixi^s, 1:ndim)
1027 double precision, intent(inout) :: cmax(ixi^s,1:number_species)
1028 double precision, intent(inout), optional :: cmin(ixi^s,1:number_species)
1029 double precision, intent(in) :: hspeed(ixi^s,1:number_species)
1030
1031 double precision :: wmean(ixi^s,nw)
1032 double precision, dimension(ixI^S) :: umean, dmean, csoundl, csoundr, tmp1,tmp2,tmp3
1033 integer :: ix^d
1034
1035 select case(boundspeed)
1036 case (1)
1037 ! This implements formula (10.52) from "Riemann Solvers and Numerical
1038 ! Methods for Fluid Dynamics" by Toro.
1039
1040 tmp1(ixo^s)=dsqrt(wlp(ixo^s,rho_))
1041 tmp2(ixo^s)=dsqrt(wrp(ixo^s,rho_))
1042 tmp3(ixo^s)=1.d0/(dsqrt(wlp(ixo^s,rho_))+dsqrt(wrp(ixo^s,rho_)))
1043 umean(ixo^s)=(wlp(ixo^s,mom(idim))*tmp1(ixo^s)+wrp(ixo^s,mom(idim))*tmp2(ixo^s))*tmp3(ixo^s)
1044
1045 if(rhd_energy) then
1046 csoundl(ixo^s)=rhd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1047 csoundr(ixo^s)=rhd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1048 else
1049 select case (rhd_pressure)
1050 case ('Trad')
1051 csoundl(ixo^s)=rhd_gamma*kbmpmua4*wlp(ixo^s,r_e)**(1.d0/4)
1052 csoundr(ixo^s)=rhd_gamma*kbmpmua4*wrp(ixo^s,r_e)**(1.d0/4)
1053 case ('adiabatic')
1054 csoundl(ixo^s)=rhd_gamma*rhd_adiab*wlp(ixo^s,rho_)**(rhd_gamma-one)
1055 csoundr(ixo^s)=rhd_gamma*rhd_adiab*wrp(ixo^s,rho_)**(rhd_gamma-one)
1056 end select
1057 end if
1058
1059 dmean(ixo^s) = (tmp1(ixo^s)*csoundl(ixo^s)+tmp2(ixo^s)*csoundr(ixo^s)) * &
1060 tmp3(ixo^s) + 0.5d0*tmp1(ixo^s)*tmp2(ixo^s)*tmp3(ixo^s)**2 * &
1061 (wrp(ixo^s,mom(idim))-wlp(ixo^s,mom(idim)))**2
1062
1063 dmean(ixo^s)=dsqrt(dmean(ixo^s))
1064 if(present(cmin)) then
1065 cmin(ixo^s,1)=umean(ixo^s)-dmean(ixo^s)
1066 cmax(ixo^s,1)=umean(ixo^s)+dmean(ixo^s)
1067 if(h_correction) then
1068 {do ix^db=ixomin^db,ixomax^db\}
1069 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1070 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1071 {end do\}
1072 end if
1073 else
1074 cmax(ixo^s,1)=dabs(umean(ixo^s))+dmean(ixo^s)
1075 end if
1076
1077 if (rhd_dust) then
1078 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1079 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1080 end if
1081
1082 case (2)
1083 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1084 tmp1(ixo^s)=wmean(ixo^s,mom(idim))/wmean(ixo^s,rho_)
1085 call rhd_get_csound2(wmean,x,ixi^l,ixo^l,csoundr)
1086 csoundr(ixo^s) = dsqrt(csoundr(ixo^s))
1087
1088 if(present(cmin)) then
1089 cmax(ixo^s,1)=max(tmp1(ixo^s)+csoundr(ixo^s),zero)
1090 cmin(ixo^s,1)=min(tmp1(ixo^s)-csoundr(ixo^s),zero)
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)=dabs(tmp1(ixo^s))+csoundr(ixo^s)
1099 end if
1100
1101 if (rhd_dust) then
1102 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1103 end if
1104 case (3)
1105 ! Miyoshi 2005 JCP 208, 315 equation (67)
1106 if(rhd_energy) then
1107 csoundl(ixo^s)=rhd_gamma*wlp(ixo^s,p_)/wlp(ixo^s,rho_)
1108 csoundr(ixo^s)=rhd_gamma*wrp(ixo^s,p_)/wrp(ixo^s,rho_)
1109 else
1110 select case (rhd_pressure)
1111 case ('Trad')
1112 csoundl(ixo^s)=rhd_gamma*kbmpmua4*wlp(ixo^s,r_e)**(1.d0/4)
1113 csoundr(ixo^s)=rhd_gamma*kbmpmua4*wrp(ixo^s,r_e)**(1.d0/4)
1114 case ('adiabatic')
1115 csoundl(ixo^s)=rhd_gamma*rhd_adiab*wlp(ixo^s,rho_)**(rhd_gamma-one)
1116 csoundr(ixo^s)=rhd_gamma*rhd_adiab*wrp(ixo^s,rho_)**(rhd_gamma-one)
1117 end select
1118 end if
1119 csoundl(ixo^s)=max(dsqrt(csoundl(ixo^s)),dsqrt(csoundr(ixo^s)))
1120 if(present(cmin)) then
1121 cmin(ixo^s,1)=min(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))-csoundl(ixo^s)
1122 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1123 if(h_correction) then
1124 {do ix^db=ixomin^db,ixomax^db\}
1125 cmin(ix^d,1)=sign(one,cmin(ix^d,1))*max(abs(cmin(ix^d,1)),hspeed(ix^d,1))
1126 cmax(ix^d,1)=sign(one,cmax(ix^d,1))*max(abs(cmax(ix^d,1)),hspeed(ix^d,1))
1127 {end do\}
1128 end if
1129 else
1130 cmax(ixo^s,1)=max(wlp(ixo^s,mom(idim)),wrp(ixo^s,mom(idim)))+csoundl(ixo^s)
1131 end if
1132 if (rhd_dust) then
1133 wmean(ixo^s,1:nwflux)=0.5d0*(wlc(ixo^s,1:nwflux)+wrc(ixo^s,1:nwflux))
1134 call dust_get_cmax(wmean, x, ixi^l, ixo^l, idim, cmax, cmin)
1135 end if
1136 end select
1137
1138 end subroutine rhd_get_cbounds
1139
1140 !> Calculate the square of the thermal sound speed csound2 within ixO^L.
1141 !> csound2=gamma*p/rho
1142 subroutine rhd_get_csound2(w,x,ixI^L,ixO^L,csound2)
1144 integer, intent(in) :: ixi^l, ixo^l
1145 double precision, intent(in) :: w(ixi^s,nw)
1146 double precision, intent(in) :: x(ixi^s,1:ndim)
1147 double precision, intent(out) :: csound2(ixi^s)
1148
1149 call rhd_get_ptot(w,x,ixi^l,ixo^l,csound2)
1150 csound2(ixo^s)=max(rhd_gamma,4.d0/3.d0)*csound2(ixo^s)/w(ixo^s,rho_)
1151 !> Turner & Stone 2001
1152
1153 end subroutine rhd_get_csound2
1154
1155 !> Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L
1156 subroutine rhd_get_pthermal(w, x, ixI^L, ixO^L, pth)
1160
1161 integer, intent(in) :: ixi^l, ixo^l
1162 double precision, intent(in) :: w(ixi^s, 1:nw)
1163 double precision, intent(in) :: x(ixi^s, 1:ndim)
1164 double precision, intent(out):: pth(ixi^s)
1165 integer :: iw, ix^d
1166
1167 if (rhd_energy) then
1168 pth(ixi^s) = (rhd_gamma - 1.d0) * (w(ixi^s, e_) - &
1169 0.5d0 * sum(w(ixi^s, mom(:))**2, dim=ndim+1) / w(ixi^s, rho_))
1170 else
1171 if (.not. associated(usr_set_pthermal)) then
1172 select case (rhd_pressure)
1173 case ('Trad')
1174 pth(ixi^s) = (w(ixi^s,r_e)*unit_pressure/const_rad_a)**0.25d0&
1175 /unit_temperature*w(ixi^s, rho_)
1176 case ('adiabatic')
1177 pth(ixi^s) = rhd_adiab * w(ixi^s, rho_)**rhd_gamma
1178 case ('Tcond') !> Thermal conduction?!
1179 pth(ixi^s) = (rhd_gamma-1.d0)*w(ixi^s,r_e)
1180 case default
1181 call mpistop('rhd_pressure unknown, use Trad or adiabatic')
1182 end select
1183 else
1184 call usr_set_pthermal(w,x,ixi^l,ixo^l,pth)
1185 end if
1186 end if
1187
1188 if (fix_small_values) then
1189 {do ix^db= ixo^lim^db\}
1190 if(pth(ix^d)<small_pressure) then
1191 pth(ix^d)=small_pressure
1192 endif
1193 {enddo^d&\}
1194 else if (check_small_values) then
1195 {do ix^db= ixo^lim^db\}
1196 if(pth(ix^d)<small_pressure) then
1197 write(*,*) "Error: small value of gas pressure",pth(ix^d),&
1198 " encountered when call rhd_get_pthermal"
1199 write(*,*) "Iteration: ", it, " Time: ", global_time
1200 write(*,*) "Location: ", x(ix^d,:)
1201 write(*,*) "Cell number: ", ix^d
1202 do iw=1,nw
1203 write(*,*) trim(cons_wnames(iw)),": ",w(ix^d,iw)
1204 end do
1205 ! use erroneous arithmetic operation to crash the run
1206 if(trace_small_values) write(*,*) dsqrt(pth(ix^d)-bigdouble)
1207 write(*,*) "Saving status at the previous time step"
1208 crash=.true.
1209 end if
1210 {enddo^d&\}
1211 end if
1212
1213 end subroutine rhd_get_pthermal
1214
1215 !> Calculate radiation pressure within ixO^L
1216 subroutine rhd_get_pradiation(w, x, ixI^L, ixO^L, prad)
1218 use mod_fld
1219 use mod_afld
1220
1221 integer, intent(in) :: ixi^l, ixo^l
1222 double precision, intent(in) :: w(ixi^s, 1:nw)
1223 double precision, intent(in) :: x(ixi^s, 1:ndim)
1224 double precision, intent(out):: prad(ixo^s, 1:ndim, 1:ndim)
1225
1226 select case (rhd_radiation_formalism)
1227 case('fld')
1228 call fld_get_radpress(w, x, ixi^l, ixo^l, prad, nghostcells)
1229 case('afld')
1230 call afld_get_radpress(w, x, ixi^l, ixo^l, prad, nghostcells)
1231 case default
1232 call mpistop('Radiation formalism unknown')
1233 end select
1234 end subroutine rhd_get_pradiation
1235
1236 !> calculates the sum of the gas pressure and max Prad tensor element
1237 subroutine rhd_get_ptot(w, x, ixI^L, ixO^L, ptot)
1239
1240 integer, intent(in) :: ixi^l, ixo^l
1241 double precision, intent(in) :: w(ixi^s, 1:nw)
1242 double precision, intent(in) :: x(ixi^s, 1:ndim)
1243 double precision :: pth(ixi^s)
1244 double precision :: prad_tensor(ixo^s, 1:ndim, 1:ndim)
1245 double precision :: prad_max(ixo^s)
1246 double precision, intent(out):: ptot(ixi^s)
1247 integer :: ix^d
1248
1249 call rhd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1250 call rhd_get_pradiation(w, x, ixi^l, ixo^l, prad_tensor)
1251
1252 {do ix^d = ixomin^d,ixomax^d\}
1253 prad_max(ix^d) = maxval(prad_tensor(ix^d,:,:))
1254 {enddo\}
1255
1256 !> filter cmax
1257 if (radio_acoustic_filter) then
1258 call rhd_radio_acoustic_filter(x, ixi^l, ixo^l, prad_max)
1259 endif
1260
1261 ptot(ixo^s) = pth(ixo^s) + prad_max(ixo^s)
1262
1263 end subroutine rhd_get_ptot
1264
1265 !> Filter peaks in cmax due to radiation energy density, used for debugging
1266 subroutine rhd_radio_acoustic_filter(x, ixI^L, ixO^L, prad_max)
1268
1269 integer, intent(in) :: ixi^l, ixo^l
1270 double precision, intent(in) :: x(ixi^s, 1:ndim)
1271 double precision, intent(inout) :: prad_max(ixo^s)
1272
1273 double precision :: tmp_prad(ixi^s)
1274 integer :: ix^d, filter, idim
1275
1276 if (size_ra_filter .lt. 1) call mpistop("ra filter of size < 1 makes no sense")
1277 if (size_ra_filter .gt. nghostcells) call mpistop("ra filter of size < nghostcells makes no sense")
1278
1279 tmp_prad(ixi^s) = zero
1280 tmp_prad(ixo^s) = prad_max(ixo^s)
1281
1282 do filter = 1,size_ra_filter
1283 do idim = 1,ndim
1284 ! {do ix^D = ixOmin^D+filter,ixOmax^D-filter\}
1285 {do ix^d = ixomin^d,ixomax^d\}
1286 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d+filter*kr(idim,^d)))
1287 prad_max(ix^d) = min(tmp_prad(ix^d),tmp_prad(ix^d-filter*kr(idim,^d)))
1288 {enddo\}
1289 enddo
1290 enddo
1291 end subroutine rhd_radio_acoustic_filter
1292
1293 !> Calculate temperature=p/rho when in e_ the total energy is stored
1294 subroutine rhd_get_temperature_from_etot(w, x, ixI^L, ixO^L, res)
1296 integer, intent(in) :: ixi^l, ixo^l
1297 double precision, intent(in) :: w(ixi^s, 1:nw)
1298 double precision, intent(in) :: x(ixi^s, 1:ndim)
1299 double precision, intent(out):: res(ixi^s)
1300
1301 double precision :: r(ixi^s)
1302
1303 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1304 call rhd_get_pthermal(w, x, ixi^l, ixo^l, res)
1305 res(ixo^s)=res(ixo^s)/(r(ixo^s)*w(ixo^s,rho_))
1306 end subroutine rhd_get_temperature_from_etot
1307
1308
1309 !> Calculate temperature=p/rho when in e_ the internal energy is stored
1310 subroutine rhd_get_temperature_from_eint(w, x, ixI^L, ixO^L, res)
1312 integer, intent(in) :: ixi^l, ixo^l
1313 double precision, intent(in) :: w(ixi^s, 1:nw)
1314 double precision, intent(in) :: x(ixi^s, 1:ndim)
1315 double precision, intent(out):: res(ixi^s)
1316 double precision :: r(ixi^s)
1317
1318 call rhd_get_rfactor(w,x,ixi^l,ixo^l,r)
1319 res(ixo^s) = (rhd_gamma - 1.0d0) * w(ixo^s, e_)/(w(ixo^s,rho_)*r(ixo^s))
1320 end subroutine rhd_get_temperature_from_eint
1321
1322 !> Calculates gas temperature
1323 subroutine rhd_get_tgas(w, x, ixI^L, ixO^L, tgas)
1325
1326 integer, intent(in) :: ixi^l, ixo^l
1327 double precision, intent(in) :: w(ixi^s, 1:nw)
1328 double precision, intent(in) :: x(ixi^s, 1:ndim)
1329 double precision :: pth(ixi^s)
1330 double precision, intent(out):: tgas(ixi^s)
1331
1332 call rhd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1333 tgas(ixi^s) = pth(ixi^s)/w(ixi^s,rho_)
1334
1335 end subroutine rhd_get_tgas
1336
1337 !> Calculates radiation temperature
1338 subroutine rhd_get_trad(w, x, ixI^L, ixO^L, trad)
1340 use mod_constants
1341
1342 integer, intent(in) :: ixi^l, ixo^l
1343 double precision, intent(in) :: w(ixi^s, 1:nw)
1344 double precision, intent(in) :: x(ixi^s, 1:ndim)
1345 double precision, intent(out):: trad(ixi^s)
1346
1347 trad(ixi^s) = (w(ixi^s,r_e)*unit_pressure&
1348 /const_rad_a)**(1.d0/4.d0)/unit_temperature
1349
1350 end subroutine rhd_get_trad
1351
1352
1353 !these are very similar to the subroutines without 1, used in mod_thermal_conductivity
1354 !but no check on whether energy variable is present
1355 subroutine rhd_ei_to_e1(ixI^L,ixO^L,w,x)
1357 integer, intent(in) :: ixi^l, ixo^l
1358 double precision, intent(inout) :: w(ixi^s, nw)
1359 double precision, intent(in) :: x(ixi^s, 1:ndim)
1360
1361 ! Calculate total energy from internal and kinetic energy
1362 w(ixo^s,e_)=w(ixo^s,e_)&
1363 +rhd_kin_en(w,ixi^l,ixo^l)
1364
1365 end subroutine rhd_ei_to_e1
1366
1367 !> Transform total energy to internal energy
1368 !but no check on whether energy variable is present
1369 subroutine rhd_e_to_ei1(ixI^L,ixO^L,w,x)
1371 integer, intent(in) :: ixi^l, ixo^l
1372 double precision, intent(inout) :: w(ixi^s, nw)
1373 double precision, intent(in) :: x(ixi^s, 1:ndim)
1374
1375 ! Calculate ei = e - ek
1376 w(ixo^s,e_)=w(ixo^s,e_)&
1377 -rhd_kin_en(w,ixi^l,ixo^l)
1378
1379 end subroutine rhd_e_to_ei1
1380
1381 ! Calculate flux f_idim[iw]
1382 subroutine rhd_get_flux_cons(w, x, ixI^L, ixO^L, idim, f)
1384 use mod_dust, only: dust_get_flux
1385
1386 integer, intent(in) :: ixi^l, ixo^l, idim
1387 double precision, intent(in) :: w(ixi^s, 1:nw), x(ixi^s, 1:ndim)
1388 double precision, intent(out) :: f(ixi^s, nwflux)
1389 double precision :: pth(ixi^s), v(ixi^s),frame_vel(ixi^s)
1390 integer :: idir, itr
1391
1392 call rhd_get_pthermal(w, x, ixi^l, ixo^l, pth)
1393 call rhd_get_v(w, x, ixi^l, ixo^l, idim, v)
1394
1395 f(ixo^s, rho_) = v(ixo^s) * w(ixo^s, rho_)
1396
1397 ! Momentum flux is v_i*m_i, +p in direction idim
1398 do idir = 1, ndir
1399 f(ixo^s, mom(idir)) = v(ixo^s) * w(ixo^s, mom(idir))
1400 end do
1401
1402 f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1403
1404 if(rhd_energy) then
1405 ! Energy flux is v_i*(e + p)
1406 f(ixo^s, e_) = v(ixo^s) * (w(ixo^s, e_) + pth(ixo^s))
1407 end if
1408
1409 if (rhd_radiation_advection) then
1410 f(ixo^s, r_e) = v(ixo^s) * w(ixo^s,r_e)
1411 else
1412 f(ixo^s, r_e) = zero
1413 endif
1414
1415 do itr = 1, rhd_n_tracer
1416 f(ixo^s, tracer(itr)) = v(ixo^s) * w(ixo^s, tracer(itr))
1417 end do
1418
1419 ! Dust fluxes
1420 if (rhd_dust) then
1421 call dust_get_flux(w, x, ixi^l, ixo^l, idim, f)
1422 end if
1423
1424 end subroutine rhd_get_flux_cons
1425
1426 ! Calculate flux f_idim[iw]
1427 subroutine rhd_get_flux(wC, w, x, ixI^L, ixO^L, idim, f)
1429 use mod_dust, only: dust_get_flux_prim
1430
1431 integer, intent(in) :: ixi^l, ixo^l, idim
1432 ! conservative w
1433 double precision, intent(in) :: wc(ixi^s, 1:nw)
1434 ! primitive w
1435 double precision, intent(in) :: w(ixi^s, 1:nw)
1436 double precision, intent(in) :: x(ixi^s, 1:ndim)
1437 double precision, intent(out) :: f(ixi^s, nwflux)
1438 double precision :: pth(ixi^s),frame_vel(ixi^s)
1439 integer :: idir, itr
1440
1441 if (rhd_energy) then
1442 pth(ixo^s) = w(ixo^s,p_)
1443 else
1444 call rhd_get_pthermal(wc, x, ixi^l, ixo^l, pth)
1445 end if
1446
1447 f(ixo^s, rho_) = w(ixo^s,mom(idim)) * w(ixo^s, rho_)
1448
1449 ! Momentum flux is v_i*m_i, +p in direction idim
1450 do idir = 1, ndir
1451 f(ixo^s, mom(idir)) = w(ixo^s,mom(idim)) * wc(ixo^s, mom(idir))
1452 end do
1453
1454 f(ixo^s, mom(idim)) = f(ixo^s, mom(idim)) + pth(ixo^s)
1455
1456 if(rhd_energy) then
1457 ! Energy flux is v_i*(e + p)
1458 f(ixo^s, e_) = w(ixo^s,mom(idim)) * (wc(ixo^s, e_) + w(ixo^s,p_))
1459 end if
1460
1461 if (rhd_radiation_advection) then
1462 f(ixo^s, r_e) = w(ixo^s,mom(idim)) * wc(ixo^s,r_e)
1463 else
1464 f(ixo^s, r_e) = zero
1465 endif
1466
1467 do itr = 1, rhd_n_tracer
1468 f(ixo^s, tracer(itr)) = w(ixo^s,mom(idim)) * w(ixo^s, tracer(itr))
1469 end do
1470
1471 ! Dust fluxes
1472 if (rhd_dust) then
1473 call dust_get_flux_prim(w, x, ixi^l, ixo^l, idim, f)
1474 end if
1475
1476 end subroutine rhd_get_flux
1477
1478 !> Add geometrical source terms to w
1479 !>
1480 !> Notice that the expressions of the geometrical terms depend only on ndir,
1481 !> not ndim. Eg, they are the same in 2.5D and in 3D, for any geometry.
1482 !>
1483 subroutine rhd_add_source_geom(qdt, dtfactor, ixI^L, ixO^L, wCT, wprim, w, x)
1488 use mod_geometry
1489 integer, intent(in) :: ixi^l, ixo^l
1490 double precision, intent(in) :: qdt, dtfactor, x(ixi^s, 1:ndim)
1491 double precision, intent(inout) :: wct(ixi^s, 1:nw), wprim(ixi^s,1:nw), w(ixi^s, 1:nw)
1492 ! to change and to set as a parameter in the parfile once the possibility to
1493 ! solve the equations in an angular momentum conserving form has been
1494 ! implemented (change tvdlf.t eg)
1495 double precision :: pth(ixi^s), source(ixi^s), minrho
1496 integer :: iw,idir, h1x^l{^nooned, h2x^l}
1497 integer :: mr_,mphi_ ! Polar var. names
1498 integer :: irho, ifluid, n_fluids
1499 double precision :: exp_factor(ixi^s), del_exp_factor(ixi^s), exp_factor_primitive(ixi^s)
1500
1501 if (rhd_dust) then
1502 n_fluids = 1 + dust_n_species
1503 else
1504 n_fluids = 1
1505 end if
1506
1507 select case (coordinate)
1508
1510 !the user provides the functions of exp_factor and del_exp_factor
1511 if(associated(usr_set_surface)) call usr_set_surface(ixi^l,x,block%dx,exp_factor,del_exp_factor,exp_factor_primitive)
1512 if(rhd_energy) then
1513 source(ixo^s) =wprim(ixo^s, p_)
1514 else
1515 if(.not. associated(usr_set_pthermal)) then
1516 source(ixo^s)=rhd_adiab * wprim(ixo^s, rho_)**rhd_gamma
1517 else
1518 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1519 end if
1520 end if
1521 source(ixo^s) = source(ixo^s)*del_exp_factor(ixo^s)/exp_factor(ixo^s)
1522 w(ixo^s,mom(1)) = w(ixo^s,mom(1)) + qdt*source(ixo^s)
1523
1524 case (cylindrical)
1525 if ((rhd_radiation_formalism .eq. 'fld') .or. (rhd_radiation_formalism .eq. 'afld')) then
1526 call mpistop("Diffusion term not implemented yet with cylkindrical geometries")
1527 end if
1528
1529 do ifluid = 0, n_fluids-1
1530 ! s[mr]=(pthermal+mphi**2/rho)/radius
1531 if (ifluid == 0) then
1532 ! gas
1533 irho = rho_
1534 mr_ = mom(r_)
1535 if(phi_>0) mphi_ = mom(phi_)
1536 if(rhd_energy) then
1537 source(ixo^s) =wprim(ixo^s, p_)
1538 else
1539 if(.not. associated(usr_set_pthermal)) then
1540 source(ixo^s)=rhd_adiab * wprim(ixo^s, rho_)**rhd_gamma
1541 else
1542 call usr_set_pthermal(wct,x,ixi^l,ixo^l,source)
1543 end if
1544 end if
1545 minrho = 0.0d0
1546 else
1547 ! dust : no pressure
1548 irho = dust_rho(ifluid)
1549 mr_ = dust_mom(r_, ifluid)
1550 if(phi_>0) mphi_ = dust_mom(phi_, ifluid)
1551 source(ixi^s) = zero
1552 minrho = 0.0d0
1553 end if
1554 if (phi_ > 0) then
1555 source(ixo^s) = source(ixo^s) + wprim(ixo^s, mphi_)**2 * wprim(ixo^s, irho)
1556 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1557 ! s[mphi]=(-vphi*vr*rho)/radius
1558 source(ixo^s) = -wprim(ixo^s, mphi_) * wprim(ixo^s, mr_) * wprim(ixo^s, irho)
1559 w(ixo^s, mphi_) = w(ixo^s, mphi_) + qdt * source(ixo^s) / x(ixo^s, r_)
1560 else
1561 ! s[mr]=2pthermal/radius
1562 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, r_)
1563 end if
1564 end do
1565 case (spherical)
1566
1567 if ((rhd_radiation_formalism .eq. 'fld') .or. (rhd_radiation_formalism .eq. 'afld')) then
1568 call mpistop("Diffusion term not implemented yet with spherical geometries")
1569 end if
1570
1571 if (rhd_dust) then
1572 call mpistop("Dust geom source terms not implemented yet with spherical geometries")
1573 end if
1574 mr_ = mom(r_)
1575 if(phi_>0) mphi_ = mom(phi_)
1576 h1x^l=ixo^l-kr(1,^d); {^nooned h2x^l=ixo^l-kr(2,^d);}
1577 ! s[mr]=((vtheta**2+vphi**2)*rho+2*p)/r
1578 if(rhd_energy) then
1579 pth(ixo^s) =wprim(ixo^s, p_)
1580 else
1581 if(.not. associated(usr_set_pthermal)) then
1582 pth(ixo^s)=rhd_adiab * wprim(ixo^s, rho_)**rhd_gamma
1583 else
1584 call usr_set_pthermal(wct,x,ixi^l,ixo^l,pth)
1585 end if
1586 end if
1587 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1588 *(block%surfaceC(ixo^s, 1) - block%surfaceC(h1x^s, 1)) &
1589 /block%dvolume(ixo^s)
1590 do idir = 2, ndir
1591 source(ixo^s) = source(ixo^s) + wprim(ixo^s, mom(idir))**2 * wprim(ixo^s, rho_)
1592 end do
1593 w(ixo^s, mr_) = w(ixo^s, mr_) + qdt * source(ixo^s) / x(ixo^s, 1)
1594
1595 {^nooned
1596 ! s[mtheta]=-(vr*vtheta*rho)/r+cot(theta)*(vphi**2*rho+p)/r
1597 source(ixo^s) = pth(ixo^s) * x(ixo^s, 1) &
1598 * (block%surfaceC(ixo^s, 2) - block%surfaceC(h2x^s, 2)) &
1599 / block%dvolume(ixo^s)
1600 if (ndir == 3) then
1601 source(ixo^s) = source(ixo^s) + (wprim(ixo^s, mom(3))**2 * wprim(ixo^s, rho_)) / tan(x(ixo^s, 2))
1602 end if
1603 source(ixo^s) = source(ixo^s) - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)
1604 w(ixo^s, mom(2)) = w(ixo^s, mom(2)) + qdt * source(ixo^s) / x(ixo^s, 1)
1605
1606 if (ndir == 3) then
1607 ! s[mphi]=-(vphi*vr/rho)/r-cot(theta)*(vtheta*vphi/rho)/r
1608 source(ixo^s) = -(wprim(ixo^s, mom(3)) * wprim(ixo^s, mr_)) * wprim(ixo^s, rho_)&
1609 - (wprim(ixo^s, mom(2)) * wprim(ixo^s, mom(3))) * wprim(ixo^s, rho_) / tan(x(ixo^s, 2))
1610 w(ixo^s, mom(3)) = w(ixo^s, mom(3)) + qdt * source(ixo^s) / x(ixo^s, 1)
1611 end if
1612 }
1613 end select
1614
1615 if (rhd_rotating_frame) then
1616 if (rhd_dust) then
1617 call mpistop("Rotating frame not implemented yet with dust")
1618 else
1619 call rotating_frame_add_source(qdt,dtfactor,ixi^l,ixo^l,wprim,w,x)
1620 end if
1621 end if
1622
1623 end subroutine rhd_add_source_geom
1624
1625 ! w[iw]= w[iw]+qdt*S[wCT, qtC, x] where S is the source based on wCT within ixO
1626 subroutine rhd_add_source(qdt,dtfactor,ixI^L,ixO^L,wCT,wCTprim,w,x,qsourcesplit,active)
1631 use mod_usr_methods, only: usr_gravity
1633
1634 integer, intent(in) :: ixi^l, ixo^l
1635 double precision, intent(in) :: qdt,dtfactor
1636 double precision, intent(in) :: wct(ixi^s, 1:nw),wctprim(ixi^s,1:nw),x(ixi^s, 1:ndim)
1637 double precision, intent(inout) :: w(ixi^s, 1:nw)
1638 logical, intent(in) :: qsourcesplit
1639 logical, intent(inout) :: active
1640
1641 double precision :: gravity_field(ixi^s, 1:ndim)
1642 integer :: idust, idim
1643
1644 if(rhd_dust) then
1645 call dust_add_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1646 end if
1647
1648 if(rhd_radiative_cooling) then
1649 call radiative_cooling_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1650 qsourcesplit,active, rc_fl)
1651 end if
1652
1653 if(rhd_viscosity) then
1654 call viscosity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1655 rhd_energy,qsourcesplit,active)
1656 end if
1657
1658 if (rhd_gravity) then
1659 call gravity_add_source(qdt,ixi^l,ixo^l,wct,wctprim,w,x,&
1660 rhd_energy,qsourcesplit,active)
1661
1662 if (rhd_dust .and. qsourcesplit .eqv. grav_split) then
1663 active = .true.
1664
1665 call usr_gravity(ixi^l, ixo^l, wct, x, gravity_field)
1666 do idust = 1, dust_n_species
1667 do idim = 1, ndim
1668 w(ixo^s, dust_mom(idim, idust)) = w(ixo^s, dust_mom(idim, idust)) &
1669 + qdt * gravity_field(ixo^s, idim) * wct(ixo^s, dust_rho(idust))
1670 end do
1671 end do
1672 end if
1673 end if
1674
1675 !> This is where the radiation force and heating/cooling are added/
1676 call rhd_add_radiation_source(qdt,ixi^l,ixo^l,wct,w,x,qsourcesplit,active)
1677
1678 if(rhd_partial_ionization) then
1679 if(.not.qsourcesplit) then
1680 active = .true.
1681 call rhd_update_temperature(ixi^l,ixo^l,wct,w,x)
1682 end if
1683 end if
1684
1685 end subroutine rhd_add_source
1686
1687 subroutine rhd_add_radiation_source(qdt,ixI^L,ixO^L,wCT,w,x,qsourcesplit,active)
1688 use mod_constants
1690 use mod_usr_methods
1691 use mod_fld!, only: fld_get_diffcoef_central, get_fld_rad_force, get_fld_energy_interact
1692 use mod_afld!, only: afld_get_diffcoef_central, get_afld_rad_force, get_afld_energy_interact
1693
1694 integer, intent(in) :: ixi^l, ixo^l
1695 double precision, intent(in) :: qdt, x(ixi^s,1:ndim)
1696 double precision, intent(in) :: wct(ixi^s,1:nw)
1697 double precision, intent(inout) :: w(ixi^s,1:nw)
1698 logical, intent(in) :: qsourcesplit
1699 logical, intent(inout) :: active
1700 double precision :: cmax(ixi^s)
1701
1702 select case(rhd_radiation_formalism)
1703 case('fld')
1704 if(fld_diff_scheme .eq. 'mg') call fld_get_diffcoef_central(w, wct, x, ixi^l, ixo^l)
1705 !> radiation force
1706 if(rhd_radiation_force) call get_fld_rad_force(qdt,ixi^l,ixo^l,wct,w,x,&
1707 rhd_energy,qsourcesplit,active)
1708 call rhd_handle_small_values(.true., w, x, ixi^l, ixo^l, 'fld_e_interact')
1709 case('afld')
1710 if (fld_diff_scheme .eq. 'mg') call afld_get_diffcoef_central(w, wct, x, ixi^l, ixo^l)
1711 !> radiation force
1712 if (rhd_radiation_force) call get_afld_rad_force(qdt,ixi^l,ixo^l,wct,w,x,&
1713 rhd_energy,qsourcesplit,active)
1714 call rhd_handle_small_values(.true., w, x, ixi^l, ixo^l, 'fld_e_interact')
1715 !> photon tiring, heating and cooling
1716 if (rhd_energy) then
1717 if (rhd_energy_interact) call get_afld_energy_interact(qdt,ixi^l,ixo^l,wct,w,x,&
1718 rhd_energy,qsourcesplit,active)
1719 endif
1720 case default
1721 call mpistop('Radiation formalism unknown')
1722 end select
1723
1724 ! ! !> NOT necessary for calculation, just want to know the grid-dependent-timestep
1725 ! call rhd_get_cmax(w, x, ixI^L, ixO^L, 2, cmax)
1726 ! w(ixI^S,i_test) = cmax(ixI^S)
1727 end subroutine rhd_add_radiation_source
1728
1729 subroutine rhd_get_dt(w, ixI^L, ixO^L, dtnew, dx^D, x)
1731 use mod_dust, only: dust_get_dt
1733 use mod_gravity, only: gravity_get_dt
1734 use mod_fld, only: fld_radforce_get_dt
1736
1737 integer, intent(in) :: ixi^l, ixo^l
1738 double precision, intent(in) :: dx^d, x(ixi^s, 1:^nd)
1739 double precision, intent(in) :: w(ixi^s, 1:nw)
1740 double precision, intent(inout) :: dtnew
1741
1742 dtnew = bigdouble
1743
1744 if (.not. dt_c) then
1745
1746 if(rhd_dust) then
1747 call dust_get_dt(w, ixi^l, ixo^l, dtnew, dx^d, x)
1748 end if
1749
1750 if(rhd_radiation_force) then
1751 select case(rhd_radiation_formalism)
1752 case('fld')
1753 call fld_radforce_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1754 case('afld')
1755 call afld_radforce_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1756 case default
1757 call mpistop('Radiation formalism unknown')
1758 end select
1759 endif
1760
1761 if(rhd_viscosity) then
1762 call viscosity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1763 end if
1764
1765 if(rhd_gravity) then
1766 call gravity_get_dt(w,ixi^l,ixo^l,dtnew,dx^d,x)
1767 end if
1768 else
1769 {^ifoned dtnew = dx1*unit_velocity/const_c}
1770 {^nooned dtnew = min(dx^d*unit_velocity/const_c)}
1771 endif
1772
1773 end subroutine rhd_get_dt
1774
1775 function rhd_kin_en(w, ixI^L, ixO^L, inv_rho) result(ke)
1776 use mod_global_parameters, only: nw, ndim
1777 integer, intent(in) :: ixi^l, ixo^l
1778 double precision, intent(in) :: w(ixi^s, nw)
1779 double precision :: ke(ixo^s)
1780 double precision, intent(in), optional :: inv_rho(ixo^s)
1781
1782 if (present(inv_rho)) then
1783 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) * inv_rho
1784 else
1785 ke = 0.5d0 * sum(w(ixo^s, mom(:))**2, dim=ndim+1) / w(ixo^s, rho_)
1786 end if
1787 end function rhd_kin_en
1788
1789 function rhd_inv_rho(w, ixI^L, ixO^L) result(inv_rho)
1790 use mod_global_parameters, only: nw, ndim
1791 integer, intent(in) :: ixi^l, ixo^l
1792 double precision, intent(in) :: w(ixi^s, nw)
1793 double precision :: inv_rho(ixo^s)
1794
1795 ! Can make this more robust
1796 inv_rho = 1.0d0 / w(ixo^s, rho_)
1797 end function rhd_inv_rho
1798
1799 subroutine rhd_handle_small_values(primitive, w, x, ixI^L, ixO^L, subname)
1800 ! handles hydro (density,pressure,velocity) bootstrapping
1801 ! any negative dust density is flagged as well (and throws an error)
1802 ! small_values_method=replace also for dust
1806 logical, intent(in) :: primitive
1807 integer, intent(in) :: ixi^l,ixo^l
1808 double precision, intent(inout) :: w(ixi^s,1:nw)
1809 double precision, intent(in) :: x(ixi^s,1:ndim)
1810 character(len=*), intent(in) :: subname
1811
1812 integer :: n,idir
1813 logical :: flag(ixi^s,1:nw)
1814
1815 call rhd_check_w(primitive, ixi^l, ixo^l, w, flag)
1816
1817 if (any(flag)) then
1818 select case (small_values_method)
1819 case ("replace")
1820 where(flag(ixo^s,rho_)) w(ixo^s,rho_) = small_density
1821 do idir = 1, ndir
1822 if(small_values_fix_iw(mom(idir))) then
1823 where(flag(ixo^s,rho_)) w(ixo^s, mom(idir)) = 0.0d0
1824 end if
1825 end do
1826
1827 if (small_values_fix_iw(r_e)) then
1828 where(flag(ixo^s,r_e)) w(ixo^s,r_e) = small_r_e
1829 end if
1830
1831 if(rhd_energy)then
1832 if(small_values_fix_iw(e_)) then
1833 if(primitive) then
1834 where(flag(ixo^s,rho_)) w(ixo^s, p_) = small_pressure
1835 else
1836 where(flag(ixo^s,rho_)) w(ixo^s, e_) = small_e + rhd_kin_en(w,ixi^l,ixo^l)
1837 endif
1838 end if
1839 endif
1840
1841 if(rhd_energy) then
1842 if(primitive) then
1843 where(flag(ixo^s,e_)) w(ixo^s,p_) = small_pressure
1844 else
1845 where(flag(ixo^s,e_))
1846 ! Add kinetic energy
1847 w(ixo^s,e_) = small_e + rhd_kin_en(w,ixi^l,ixo^l)
1848 end where
1849 end if
1850 end if
1851
1852 if(rhd_dust)then
1853 do n=1,dust_n_species
1854 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1855 do idir = 1, ndir
1856 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1857 enddo
1858 enddo
1859 endif
1860 case ("average")
1861 if(primitive)then
1862 ! averaging for all primitive fields, including dust
1863 call small_values_average(ixi^l, ixo^l, w, x, flag)
1864 else
1865 ! do averaging of density
1866 call small_values_average(ixi^l, ixo^l, w, x, flag, rho_)
1867 if(rhd_energy) then
1868 ! do averaging of pressure
1869 w(ixi^s,p_)=(rhd_gamma-1.d0)*(w(ixi^s,e_) &
1870 -0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_))
1871 call small_values_average(ixi^l, ixo^l, w, x, flag, p_)
1872 w(ixi^s,e_)=w(ixi^s,p_)/(rhd_gamma-1.d0) &
1873 +0.5d0*sum(w(ixi^s, mom(:))**2, dim=ndim+1)/w(ixi^s,rho_)
1874 end if
1875 ! do averaging of radiation energy
1876 call small_values_average(ixi^l, ixo^l, w, x, flag, r_e)
1877 if(rhd_dust)then
1878 do n=1,dust_n_species
1879 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_rho(n)) = 0.0d0
1880 do idir = 1, ndir
1881 where(flag(ixo^s,dust_rho(n))) w(ixo^s,dust_mom(idir,n)) = 0.0d0
1882 enddo
1883 enddo
1884 endif
1885 endif
1886 case default
1887 if(.not.primitive) then
1888 !convert w to primitive
1889 ! Calculate pressure = (gamma-1) * (e-ek)
1890 if(rhd_energy) then
1891 w(ixo^s,p_)=(rhd_gamma-1.d0)*(w(ixo^s,e_)-rhd_kin_en(w,ixi^l,ixo^l))
1892 end if
1893 ! Convert gas momentum to velocity
1894 do idir = 1, ndir
1895 w(ixo^s, mom(idir)) = w(ixo^s, mom(idir))/w(ixo^s,rho_)
1896 end do
1897 end if
1898 ! NOTE: dust entries may still have conserved values here
1899 call small_values_error(w, x, ixi^l, ixo^l, flag, subname)
1900 end select
1901 end if
1902 end subroutine rhd_handle_small_values
1903
1904 subroutine rfactor_from_temperature_ionization(w,x,ixI^L,ixO^L,Rfactor)
1907 integer, intent(in) :: ixi^l, ixo^l
1908 double precision, intent(in) :: w(ixi^s,1:nw)
1909 double precision, intent(in) :: x(ixi^s,1:ndim)
1910 double precision, intent(out):: rfactor(ixi^s)
1911
1912 double precision :: iz_h(ixo^s),iz_he(ixo^s)
1913
1914 call ionization_degree_from_temperature(ixi^l,ixo^l,w(ixi^s,te_),iz_h,iz_he)
1915 ! assume the first and second ionization of Helium have the same degree
1916 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
1917
1918 end subroutine rfactor_from_temperature_ionization
1919
1920 subroutine rfactor_from_constant_ionization(w,x,ixI^L,ixO^L,Rfactor)
1922 integer, intent(in) :: ixi^l, ixo^l
1923 double precision, intent(in) :: w(ixi^s,1:nw)
1924 double precision, intent(in) :: x(ixi^s,1:ndim)
1925 double precision, intent(out):: rfactor(ixi^s)
1926
1927 rfactor(ixo^s)=rr
1928
1929 end subroutine rfactor_from_constant_ionization
1930
1931 subroutine rhd_update_temperature(ixI^L,ixO^L,wCT,w,x)
1934
1935 integer, intent(in) :: ixi^l, ixo^l
1936 double precision, intent(in) :: wct(ixi^s,1:nw), x(ixi^s,1:ndim)
1937 double precision, intent(inout) :: w(ixi^s,1:nw)
1938
1939 double precision :: iz_h(ixo^s),iz_he(ixo^s), pth(ixi^s)
1940
1941 call ionization_degree_from_temperature(ixi^l,ixo^l,wct(ixi^s,te_),iz_h,iz_he)
1942
1943 call rhd_get_pthermal(w,x,ixi^l,ixo^l,pth)
1944
1945 w(ixo^s,te_)=(2.d0+3.d0*he_abundance)*pth(ixo^s)/(w(ixo^s,rho_)*(1.d0+iz_h(ixo^s)+&
1946 he_abundance*(iz_he(ixo^s)*(iz_he(ixo^s)+1.d0)+1.d0)))
1947
1948 end subroutine rhd_update_temperature
1949
1950end 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
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.
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: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
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
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)
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 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)
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 viscosity_init(phys_wider_stencil)
Initialize the module.
subroutine viscosity_get_dt(w, ixil, ixol, dtnew, dxd, x)
procedure(sub_add_source), pointer, public viscosity_add_source