MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
Functions/Subroutines | Variables
mod_rhd_phys Module Reference

Radiation-Hydrodynamics physics module Module aims at solving the Hydrodynamic equations toghether with the zeroth moment of the radiative transfer equation. A closure is provided by the flux limited diffusion (FLD)-approximation in the mod_fld.t module. See [1]Moens, N., Sundqvist, J. O., El Mellah, I., Poniatowski, L., Teunissen, J., and Keppens, R., “Radiation-hydrodynamics with MPI-AMRVAC . Flux-limited diffusion”, Astronomy and Astrophysics, vol. 657, 2022. doi:10.1051/0004-6361/202141023. For more information. Another possible closure in the works is the anisotropic flux limited diffusion approximation (AFLD) described in mod_afld.t. More...

Functions/Subroutines

subroutine rhd_write_info (fh)
 Write this module's parameters to a snapsoht. More...
 
subroutine, public rhd_phys_init ()
 Initialize the module. More...
 
subroutine rhd_te_images ()
 
subroutine rhd_sts_set_source_tc_rhd (ixIL, ixOL, w, x, wres, fix_conserve_at_step, my_dt, igrid, nflux)
 
double precision function rhd_get_tc_dt_rhd (w, ixIL, ixOL, dxD, x)
 
subroutine rhd_tc_handle_small_e (w, x, ixIL, ixOL, step)
 
subroutine tc_params_read_rhd (fl)
 
subroutine rhd_get_rho (w, x, ixIL, ixOL, rho)
 
subroutine rc_params_read (fl)
 
subroutine, public rhd_check_params
 
subroutine, public rhd_set_mg_bounds
 Set the boundaries for the diffusion of E. More...
 
subroutine rhd_physical_units
 
subroutine, public rhd_check_w (primitive, ixIL, ixOL, w, flag)
 Returns logical argument flag where values are ok. More...
 
subroutine, public rhd_to_conserved (ixIL, ixOL, w, x)
 Transform primitive variables into conservative ones. More...
 
subroutine, public rhd_to_primitive (ixIL, ixOL, w, x)
 Transform conservative variables into primitive ones. More...
 
subroutine rhd_ei_to_e (ixIL, ixOL, w, x)
 Transform internal energy to total energy. More...
 
subroutine rhd_e_to_ei (ixIL, ixOL, w, x)
 Transform total energy to internal energy. More...
 
subroutine e_to_rhos (ixIL, ixOL, w, x)
 
subroutine rhos_to_e (ixIL, ixOL, w, x)
 
subroutine rhd_get_v (w, x, ixIL, ixOL, idim, v)
 Calculate v_i = m_i / rho within ixO^L. More...
 
subroutine rhd_get_cmax (w, x, ixIL, ixOL, idim, cmax)
 Calculate cmax_idim = csound + abs(v_idim) within ixO^L. More...
 
subroutine rhd_get_a2max (w, x, ixIL, ixOL, a2max)
 
subroutine rhd_get_tcutoff (ixIL, ixOL, w, x, tco_local, Tmax_local)
 get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22) More...
 
subroutine rhd_get_cbounds (wLC, wRC, wLp, wRp, x, ixIL, ixOL, idim, Hspeed, cmax, cmin)
 Calculate cmax_idim = csound + abs(v_idim) within ixO^L. More...
 
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. More...
 
subroutine, public rhd_get_pthermal (w, x, ixIL, ixOL, pth)
 Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L. More...
 
subroutine, public rhd_get_pradiation (w, x, ixIL, ixOL, prad)
 Calculate radiation pressure within ixO^L. More...
 
subroutine, public rhd_get_ptot (w, x, ixIL, ixOL, ptot)
 calculates the sum of the gas pressure and max Prad tensor element More...
 
subroutine rhd_radio_acoustic_filter (x, ixIL, ixOL, prad_max)
 Filter peaks in cmax due to radiation energy density, used for debugging. More...
 
subroutine rhd_get_temperature_from_etot (w, x, ixIL, ixOL, res)
 Calculate temperature=p/rho when in e_ the total energy is stored. More...
 
subroutine rhd_get_temperature_from_eint (w, x, ixIL, ixOL, res)
 Calculate temperature=p/rho when in e_ the internal energy is stored. More...
 
subroutine, public rhd_get_tgas (w, x, ixIL, ixOL, tgas)
 Calculates gas temperature. More...
 
subroutine, public rhd_get_trad (w, x, ixIL, ixOL, trad)
 Calculates radiation temperature. More...
 
subroutine rhd_ei_to_e1 (ixIL, ixOL, w, x)
 
subroutine rhd_e_to_ei1 (ixIL, ixOL, w, x)
 Transform total energy to internal energy. More...
 
subroutine rhd_get_flux_cons (w, x, ixIL, ixOL, idim, f)
 
subroutine rhd_get_flux (wC, w, x, ixIL, ixOL, idim, f)
 
subroutine rhd_add_source_geom (qdt, dtfactor, ixIL, ixOL, wCT, w, x)
 Add geometrical source terms to w. More...
 
subroutine rhd_add_source (qdt, dtfactor, ixIL, ixOL, wCT, wCTprim, w, x, qsourcesplit, active)
 
subroutine rhd_add_radiation_source (qdt, ixIL, ixOL, wCT, w, x, qsourcesplit, active)
 
subroutine rhd_get_dt (w, ixIL, ixOL, dtnew, dxD, x)
 
double precision function, dimension(ixo^s), public rhd_kin_en (w, ixIL, ixOL, inv_rho)
 
double precision function, dimension(ixo^s) rhd_inv_rho (w, ixIL, ixOL)
 
subroutine rhd_handle_small_values (primitive, w, x, ixIL, ixOL, subname)
 
subroutine rfactor_from_temperature_ionization (w, x, ixIL, ixOL, Rfactor)
 
subroutine rfactor_from_constant_ionization (w, x, ixIL, ixOL, Rfactor)
 
subroutine rhd_update_temperature (ixIL, ixOL, wCT, w, x)
 

Variables

logical, public, protected rhd_energy = .true.
 Whether an energy equation is used. More...
 
logical, public, protected rhd_thermal_conduction = .false.
 Whether thermal conduction is added. More...
 
type(tc_fluid), allocatable, public tc_fl
 
type(te_fluid), allocatable, public te_fl_rhd
 
logical, public, protected rhd_radiative_cooling = .false.
 Whether radiative cooling is added. More...
 
type(rc_fluid), allocatable, public rc_fl
 
logical, public, protected rhd_dust = .false.
 Whether dust is added. More...
 
logical, public, protected rhd_viscosity = .false.
 Whether viscosity is added. More...
 
logical, public, protected rhd_gravity = .false.
 Whether gravity is added. More...
 
logical, public, protected rhd_particles = .false.
 Whether particles module is added. More...
 
logical, public, protected rhd_rotating_frame = .false.
 Whether rotating frame is activated. More...
 
integer, public, protected rhd_n_tracer = 0
 Number of tracer species. More...
 
integer, public, protected rho_
 Index of the density (in the w array) More...
 
integer, dimension(:), allocatable, public, protected mom
 Indices of the momentum density. More...
 
integer, dimension(:), allocatable, public, protected tracer
 Indices of the tracers. More...
 
integer, public, protected e_
 Index of the energy density (-1 if not present) More...
 
integer, public, protected p_
 Index of the gas pressure (-1 if not present) should equal e_. More...
 
integer, public, protected r_e
 Index of the radiation energy. More...
 
integer, public, protected te_
 Indices of temperature. More...
 
integer, public, protected tcoff_
 Index of the cutoff temperature for the TRAC method. More...
 
double precision, public rhd_gamma = 5.d0/3.0d0
 The adiabatic index. More...
 
double precision, public rhd_adiab = 1.0d0
 The adiabatic constant. More...
 
double precision, public, protected small_r_e = 0.d0
 The smallest allowed radiation energy. More...
 
logical, public, protected rhd_trac = .false.
 Whether TRAC method is used. More...
 
integer, public, protected rhd_trac_type = 1
 
logical, public, protected rhd_force_diagonal = .false.
 Allows overruling default corner filling (for debug mode, since otherwise corner primitives fail) More...
 
double precision, public, protected he_abundance =0.1d0
 Helium abundance over Hydrogen. More...
 
character(len=8), public rhd_radiation_formalism = 'fld'
 Formalism to treat radiation. More...
 
character(len=8), public rhd_pressure = 'Trad'
 In the case of no rhd_energy, how to compute pressure. More...
 
logical, public, protected rhd_radiation_force = .true.
 Treat radiation fld_Rad_force. More...
 
logical, public, protected rhd_energy_interact = .true.
 Treat radiation-gas energy interaction. More...
 
logical, public, protected rhd_radiation_diffusion = .true.
 Treat radiation energy diffusion. More...
 
logical, public, protected rhd_radiation_advection = .true.
 Treat radiation advection. More...
 
logical, public, protected rhd_partial_ionization = .false.
 Whether plasma is partially ionized. More...
 
double precision, public kbmpmua4
 kb/(m_p mu)* 1/a_rad**4, More...
 
double precision, public, protected h_ion_fr =1d0
 Ionization fraction of H H_ion_fr = H+/(H+ + H) More...
 
double precision, public, protected he_ion_fr =1d0
 Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He) More...
 
double precision, public, protected he_ion_fr2 =1d0
 Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+) More...
 
double precision, public, protected rr =1d0
 
logical, public, protected eq_state_units = .true.
 

Detailed Description

Radiation-Hydrodynamics physics module Module aims at solving the Hydrodynamic equations toghether with the zeroth moment of the radiative transfer equation. A closure is provided by the flux limited diffusion (FLD)-approximation in the mod_fld.t module. See [1]Moens, N., Sundqvist, J. O., El Mellah, I., Poniatowski, L., Teunissen, J., and Keppens, R., “Radiation-hydrodynamics with MPI-AMRVAC . Flux-limited diffusion”, Astronomy and Astrophysics, vol. 657, 2022. doi:10.1051/0004-6361/202141023. For more information. Another possible closure in the works is the anisotropic flux limited diffusion approximation (AFLD) described in mod_afld.t.

Function/Subroutine Documentation

◆ e_to_rhos()

subroutine mod_rhd_phys::e_to_rhos ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Definition at line 843 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rc_params_read()

subroutine mod_rhd_phys::rc_params_read ( type(rc_fluid), intent(inout)  fl)

Name of cooling curve

Name of cooling method

Fixed temperature not lower than tlow

Lower limit of temperature

Add cooling source in a split way (.true.) or un-split way (.false.)

Definition at line 557 of file mod_rhd_phys.t.

◆ rfactor_from_constant_ionization()

subroutine mod_rhd_phys::rfactor_from_constant_ionization ( double precision, dimension(ixi^s,1:nw), intent(in)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  Rfactor 
)

Definition at line 1899 of file mod_rhd_phys.t.

◆ rfactor_from_temperature_ionization()

subroutine mod_rhd_phys::rfactor_from_temperature_ionization ( double precision, dimension(ixi^s,1:nw), intent(in)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  Rfactor 
)

Definition at line 1883 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_add_radiation_source()

subroutine mod_rhd_phys::rhd_add_radiation_source ( double precision, intent(in)  qdt,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s,1:nw), intent(in)  wCT,
double precision, dimension(ixi^s,1:nw), intent(inout)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
logical, intent(in)  qsourcesplit,
logical, intent(inout)  active 
)

radiation force

photon tiring, heating and cooling

radiation force

photon tiring, heating and cooling

Definition at line 1643 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_add_source()

subroutine mod_rhd_phys::rhd_add_source ( double precision, intent(in)  qdt,
double precision, intent(in)  dtfactor,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, 1:nw), intent(in)  wCT,
double precision, dimension(ixi^s,1:nw), intent(in)  wCTprim,
double precision, dimension(ixi^s, 1:nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
logical, intent(in)  qsourcesplit,
logical, intent(inout)  active 
)

This is where the radiation force and heating/cooling are added/

Definition at line 1582 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_add_source_geom()

subroutine mod_rhd_phys::rhd_add_source_geom ( double precision, intent(in)  qdt,
double precision, intent(in)  dtfactor,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, 1:nw), intent(inout)  wCT,
double precision, dimension(ixi^s, 1:nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Add geometrical source terms to w.

Notice that the expressions of the geometrical terms depend only on ndir, not ndim. Eg, they are the same in 2.5D and in 3D, for any geometry.

Definition at line 1454 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_check_params()

subroutine, public mod_rhd_phys::rhd_check_params

Definition at line 601 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_check_w()

subroutine, public mod_rhd_phys::rhd_check_w ( logical, intent(in)  primitive,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw), intent(in)  w,
logical, dimension(ixi^s,1:nw), intent(inout)  flag 
)

Returns logical argument flag where values are ok.

Definition at line 721 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_e_to_ei()

subroutine mod_rhd_phys::rhd_e_to_ei ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Transform total energy to internal energy.

Definition at line 831 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_e_to_ei1()

subroutine mod_rhd_phys::rhd_e_to_ei1 ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Transform total energy to internal energy.

Definition at line 1334 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_ei_to_e()

subroutine mod_rhd_phys::rhd_ei_to_e ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Transform internal energy to total energy.

Definition at line 818 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_ei_to_e1()

subroutine mod_rhd_phys::rhd_ei_to_e1 ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Definition at line 1320 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_a2max()

subroutine mod_rhd_phys::rhd_get_a2max ( double precision, dimension(ixi^s, nw), intent(in)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ndim), intent(inout)  a2max 
)

4th order

Definition at line 905 of file mod_rhd_phys.t.

◆ rhd_get_cbounds()

subroutine mod_rhd_phys::rhd_get_cbounds ( double precision, dimension(ixi^s, nw), intent(in)  wLC,
double precision, dimension(ixi^s, nw), intent(in)  wRC,
double precision, dimension(ixi^s, nw), intent(in)  wLp,
double precision, dimension(ixi^s, nw), intent(in)  wRp,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
integer, intent(in)  idim,
double precision, dimension(ixi^s,1:number_species), intent(in)  Hspeed,
double precision, dimension(ixi^s,1:number_species), intent(inout)  cmax,
double precision, dimension(ixi^s,1:number_species), intent(inout), optional  cmin 
)

Calculate cmax_idim = csound + abs(v_idim) within ixO^L.

Definition at line 981 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_cmax()

subroutine mod_rhd_phys::rhd_get_cmax ( double precision, dimension(ixi^s, nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
integer, intent(in)  idim,
double precision, dimension(ixi^s), intent(inout)  cmax 
)

Calculate cmax_idim = csound + abs(v_idim) within ixO^L.

Definition at line 884 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_csound2()

subroutine, public mod_rhd_phys::rhd_get_csound2 ( double precision, dimension(ixi^s,nw), intent(in)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  csound2 
)

Calculate the square of the thermal sound speed csound2 within ixO^L. csound2=gamma*p/rho.

Turner & Stone 2001

Definition at line 1107 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_dt()

subroutine mod_rhd_phys::rhd_get_dt ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, intent(inout)  dtnew,
double precision, intent(in)  dx,
double precision, intent(in)  D,
double precision, dimension(ixi^s, 1:^nd), intent(in)  x 
)

Definition at line 1703 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_flux()

subroutine mod_rhd_phys::rhd_get_flux ( double precision, dimension(ixi^s, 1:nw), intent(in)  wC,
double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
integer, intent(in)  idim,
double precision, dimension(ixi^s, nwflux), intent(out)  f 
)

Definition at line 1392 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_flux_cons()

subroutine mod_rhd_phys::rhd_get_flux_cons ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
integer, intent(in)  idim,
double precision, dimension(ixi^s, nwflux), intent(out)  f 
)

Definition at line 1347 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_pradiation()

subroutine, public mod_rhd_phys::rhd_get_pradiation ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixo^s, 1:ndim, 1:ndim), intent(out)  prad 
)

Calculate radiation pressure within ixO^L.

Definition at line 1181 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_pthermal()

subroutine, public mod_rhd_phys::rhd_get_pthermal ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  pth 
)

Calculate thermal pressure=(gamma-1)*(e-0.5*m**2/rho) within ixO^L.

Thermal conduction?!

Definition at line 1121 of file mod_rhd_phys.t.

◆ rhd_get_ptot()

subroutine, public mod_rhd_phys::rhd_get_ptot ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  ptot 
)

calculates the sum of the gas pressure and max Prad tensor element

filter cmax

Definition at line 1202 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_rho()

subroutine mod_rhd_phys::rhd_get_rho ( double precision, dimension(ixi^s,1:nw), intent(in)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  rho 
)

Definition at line 545 of file mod_rhd_phys.t.

◆ rhd_get_tc_dt_rhd()

double precision function mod_rhd_phys::rhd_get_tc_dt_rhd ( double precision, dimension(ixi^s,1:nw), intent(in)  w,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, intent(in)  dx,
double precision, intent(in)  D,
double precision, dimension(ixi^s,1:ndim), intent(in)  x 
)

Definition at line 474 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_tcutoff()

subroutine mod_rhd_phys::rhd_get_tcutoff ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s,1:nw), intent(inout)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
double precision, intent(out)  tco_local,
double precision, intent(out)  Tmax_local 
)

get adaptive cutoff temperature for TRAC (Johnston 2019 ApJL, 873, L22)

iijima et al. 2021, LTRAC method

Definition at line 928 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_temperature_from_eint()

subroutine mod_rhd_phys::rhd_get_temperature_from_eint ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  res 
)

Calculate temperature=p/rho when in e_ the internal energy is stored.

Definition at line 1275 of file mod_rhd_phys.t.

◆ rhd_get_temperature_from_etot()

subroutine mod_rhd_phys::rhd_get_temperature_from_etot ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  res 
)

Calculate temperature=p/rho when in e_ the total energy is stored.

Definition at line 1259 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_tgas()

subroutine, public mod_rhd_phys::rhd_get_tgas ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  tgas 
)

Calculates gas temperature.

Definition at line 1288 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_get_trad()

subroutine, public mod_rhd_phys::rhd_get_trad ( double precision, dimension(ixi^s, 1:nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s), intent(out)  trad 
)

Calculates radiation temperature.

Definition at line 1303 of file mod_rhd_phys.t.

◆ rhd_get_v()

subroutine mod_rhd_phys::rhd_get_v ( double precision, dimension(ixi^s, nw), intent(in)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
integer, intent(in)  idim,
double precision, dimension(ixi^s), intent(out)  v 
)

Calculate v_i = m_i / rho within ixO^L.

Definition at line 874 of file mod_rhd_phys.t.

◆ rhd_handle_small_values()

subroutine mod_rhd_phys::rhd_handle_small_values ( logical, intent(in)  primitive,
double precision, dimension(ixi^s,1:nw), intent(inout)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
character(len=*), intent(in)  subname 
)

Definition at line 1778 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_inv_rho()

double precision function, dimension(ixo^s) mod_rhd_phys::rhd_inv_rho ( double precision, dimension(ixi^s, nw), intent(in)  w,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L 
)

Definition at line 1768 of file mod_rhd_phys.t.

◆ rhd_kin_en()

double precision function, dimension(ixo^s), public mod_rhd_phys::rhd_kin_en ( double precision, dimension(ixi^s, nw), intent(in)  w,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixo^s), intent(in), optional  inv_rho 
)

Definition at line 1754 of file mod_rhd_phys.t.

◆ rhd_phys_init()

subroutine, public mod_rhd_phys::rhd_phys_init

Initialize the module.

set radiation energy

Initiate radiation-closure module

Usefull constante

Definition at line 204 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_physical_units()

subroutine mod_rhd_phys::rhd_physical_units

Units for radiative flux and opacity

Definition at line 665 of file mod_rhd_phys.t.

◆ rhd_radio_acoustic_filter()

subroutine mod_rhd_phys::rhd_radio_acoustic_filter ( double precision, dimension(ixi^s, 1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixo^s), intent(inout)  prad_max 
)

Filter peaks in cmax due to radiation energy density, used for debugging.

Definition at line 1231 of file mod_rhd_phys.t.

◆ rhd_set_mg_bounds()

subroutine, public mod_rhd_phys::rhd_set_mg_bounds

Set the boundaries for the diffusion of E.

Definition at line 630 of file mod_rhd_phys.t.

◆ rhd_sts_set_source_tc_rhd()

subroutine mod_rhd_phys::rhd_sts_set_source_tc_rhd ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s,1:nw), intent(inout)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
double precision, dimension(ixi^s,1:nw), intent(inout)  wres,
logical, intent(in)  fix_conserve_at_step,
double precision, intent(in)  my_dt,
integer, intent(in)  igrid,
integer, intent(in)  nflux 
)

Definition at line 461 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_tc_handle_small_e()

subroutine mod_rhd_phys::rhd_tc_handle_small_e ( double precision, dimension(ixi^s,1:nw), intent(inout)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x,
integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
integer, intent(in)  step 
)

Definition at line 490 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_te_images()

subroutine mod_rhd_phys::rhd_te_images

Definition at line 443 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_to_conserved()

subroutine, public mod_rhd_phys::rhd_to_conserved ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Transform primitive variables into conservative ones.

Definition at line 752 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_to_primitive()

subroutine, public mod_rhd_phys::rhd_to_primitive ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw), intent(inout)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Transform conservative variables into primitive ones.

Definition at line 784 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_update_temperature()

subroutine mod_rhd_phys::rhd_update_temperature ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s,1:nw), intent(in)  wCT,
double precision, dimension(ixi^s,1:nw), intent(inout)  w,
double precision, dimension(ixi^s,1:ndim), intent(in)  x 
)

Definition at line 1910 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ rhd_write_info()

subroutine mod_rhd_phys::rhd_write_info ( integer, intent(in)  fh)

Write this module's parameters to a snapsoht.

Definition at line 186 of file mod_rhd_phys.t.

◆ rhos_to_e()

subroutine mod_rhd_phys::rhos_to_e ( integer, intent(in)  ixI,
integer, intent(in)  L,
integer, intent(in)  ixO,
  L,
double precision, dimension(ixi^s, nw)  w,
double precision, dimension(ixi^s, 1:ndim), intent(in)  x 
)

Definition at line 858 of file mod_rhd_phys.t.

Here is the call graph for this function:

◆ tc_params_read_rhd()

subroutine mod_rhd_phys::tc_params_read_rhd ( type(tc_fluid), intent(inout)  fl)

Definition at line 525 of file mod_rhd_phys.t.

Variable Documentation

◆ e_

integer, public, protected mod_rhd_phys::e_

Index of the energy density (-1 if not present)

Definition at line 60 of file mod_rhd_phys.t.

◆ eq_state_units

logical, public, protected mod_rhd_phys::eq_state_units = .true.

Definition at line 143 of file mod_rhd_phys.t.

◆ h_ion_fr

double precision, public, protected mod_rhd_phys::h_ion_fr =1d0

Ionization fraction of H H_ion_fr = H+/(H+ + H)

Definition at line 129 of file mod_rhd_phys.t.

◆ he_abundance

double precision, public, protected mod_rhd_phys::he_abundance =0.1d0

Helium abundance over Hydrogen.

Definition at line 94 of file mod_rhd_phys.t.

◆ he_ion_fr

double precision, public, protected mod_rhd_phys::he_ion_fr =1d0

Ionization fraction of He He_ion_fr = (He2+ + He+)/(He2+ + He+ + He)

Definition at line 132 of file mod_rhd_phys.t.

◆ he_ion_fr2

double precision, public, protected mod_rhd_phys::he_ion_fr2 =1d0

Ratio of number He2+ / number He+ + He2+ He_ion_fr2 = He2+/(He2+ + He+)

Definition at line 135 of file mod_rhd_phys.t.

◆ kbmpmua4

double precision, public mod_rhd_phys::kbmpmua4

kb/(m_p mu)* 1/a_rad**4,

Definition at line 122 of file mod_rhd_phys.t.

◆ mom

integer, dimension(:), allocatable, public, protected mod_rhd_phys::mom

Indices of the momentum density.

Definition at line 54 of file mod_rhd_phys.t.

◆ p_

integer, public, protected mod_rhd_phys::p_

Index of the gas pressure (-1 if not present) should equal e_.

Definition at line 63 of file mod_rhd_phys.t.

◆ r_e

integer, public, protected mod_rhd_phys::r_e

Index of the radiation energy.

Definition at line 66 of file mod_rhd_phys.t.

◆ rc_fl

type(rc_fluid), allocatable, public mod_rhd_phys::rc_fl

Definition at line 30 of file mod_rhd_phys.t.

◆ rhd_adiab

double precision, public mod_rhd_phys::rhd_adiab = 1.0d0

The adiabatic constant.

Definition at line 78 of file mod_rhd_phys.t.

◆ rhd_dust

logical, public, protected mod_rhd_phys::rhd_dust = .false.

Whether dust is added.

Definition at line 33 of file mod_rhd_phys.t.

◆ rhd_energy

logical, public, protected mod_rhd_phys::rhd_energy = .true.

Whether an energy equation is used.

Definition at line 21 of file mod_rhd_phys.t.

◆ rhd_energy_interact

logical, public, protected mod_rhd_phys::rhd_energy_interact = .true.

Treat radiation-gas energy interaction.

Definition at line 106 of file mod_rhd_phys.t.

◆ rhd_force_diagonal

logical, public, protected mod_rhd_phys::rhd_force_diagonal = .false.

Allows overruling default corner filling (for debug mode, since otherwise corner primitives fail)

Definition at line 91 of file mod_rhd_phys.t.

◆ rhd_gamma

double precision, public mod_rhd_phys::rhd_gamma = 5.d0/3.0d0

The adiabatic index.

Definition at line 75 of file mod_rhd_phys.t.

◆ rhd_gravity

logical, public, protected mod_rhd_phys::rhd_gravity = .false.

Whether gravity is added.

Definition at line 39 of file mod_rhd_phys.t.

◆ rhd_n_tracer

integer, public, protected mod_rhd_phys::rhd_n_tracer = 0

Number of tracer species.

Definition at line 48 of file mod_rhd_phys.t.

◆ rhd_partial_ionization

logical, public, protected mod_rhd_phys::rhd_partial_ionization = .false.

Whether plasma is partially ionized.

Definition at line 115 of file mod_rhd_phys.t.

◆ rhd_particles

logical, public, protected mod_rhd_phys::rhd_particles = .false.

Whether particles module is added.

Definition at line 42 of file mod_rhd_phys.t.

◆ rhd_pressure

character(len=8), public mod_rhd_phys::rhd_pressure = 'Trad'

In the case of no rhd_energy, how to compute pressure.

Definition at line 100 of file mod_rhd_phys.t.

◆ rhd_radiation_advection

logical, public, protected mod_rhd_phys::rhd_radiation_advection = .true.

Treat radiation advection.

Definition at line 112 of file mod_rhd_phys.t.

◆ rhd_radiation_diffusion

logical, public, protected mod_rhd_phys::rhd_radiation_diffusion = .true.

Treat radiation energy diffusion.

Definition at line 109 of file mod_rhd_phys.t.

◆ rhd_radiation_force

logical, public, protected mod_rhd_phys::rhd_radiation_force = .true.

Treat radiation fld_Rad_force.

Definition at line 103 of file mod_rhd_phys.t.

◆ rhd_radiation_formalism

character(len=8), public mod_rhd_phys::rhd_radiation_formalism = 'fld'

Formalism to treat radiation.

Definition at line 97 of file mod_rhd_phys.t.

◆ rhd_radiative_cooling

logical, public, protected mod_rhd_phys::rhd_radiative_cooling = .false.

Whether radiative cooling is added.

Definition at line 29 of file mod_rhd_phys.t.

◆ rhd_rotating_frame

logical, public, protected mod_rhd_phys::rhd_rotating_frame = .false.

Whether rotating frame is activated.

Definition at line 45 of file mod_rhd_phys.t.

◆ rhd_thermal_conduction

logical, public, protected mod_rhd_phys::rhd_thermal_conduction = .false.

Whether thermal conduction is added.

Definition at line 24 of file mod_rhd_phys.t.

◆ rhd_trac

logical, public, protected mod_rhd_phys::rhd_trac = .false.

Whether TRAC method is used.

Definition at line 87 of file mod_rhd_phys.t.

◆ rhd_trac_type

integer, public, protected mod_rhd_phys::rhd_trac_type = 1

Definition at line 88 of file mod_rhd_phys.t.

◆ rhd_viscosity

logical, public, protected mod_rhd_phys::rhd_viscosity = .false.

Whether viscosity is added.

Definition at line 36 of file mod_rhd_phys.t.

◆ rho_

integer, public, protected mod_rhd_phys::rho_

Index of the density (in the w array)

Definition at line 51 of file mod_rhd_phys.t.

◆ rr

double precision, public, protected mod_rhd_phys::rr =1d0

Definition at line 139 of file mod_rhd_phys.t.

◆ small_r_e

double precision, public, protected mod_rhd_phys::small_r_e = 0.d0

The smallest allowed radiation energy.

Definition at line 84 of file mod_rhd_phys.t.

◆ tc_fl

type(tc_fluid), allocatable, public mod_rhd_phys::tc_fl

Definition at line 25 of file mod_rhd_phys.t.

◆ tcoff_

integer, public, protected mod_rhd_phys::tcoff_

Index of the cutoff temperature for the TRAC method.

Definition at line 72 of file mod_rhd_phys.t.

◆ te_

integer, public, protected mod_rhd_phys::te_

Indices of temperature.

Definition at line 69 of file mod_rhd_phys.t.

◆ te_fl_rhd

type(te_fluid), allocatable, public mod_rhd_phys::te_fl_rhd

Definition at line 26 of file mod_rhd_phys.t.

◆ tracer

integer, dimension(:), allocatable, public, protected mod_rhd_phys::tracer

Indices of the tracers.

Definition at line 57 of file mod_rhd_phys.t.