MPI-AMRVAC  3.1
The MPI - Adaptive Mesh Refinement - Versatile Advection Code
Modules | Functions/Subroutines | Variables
mod_fld.t File Reference

Go to the source code of this file.

Modules

module  mod_fld
 Nicolas Moens Module for including flux limited diffusion (FLD)-approximation in Radiation-hydrodynamics simulations using mod_rhd Based on Turner and stone 2001. 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.
 

Functions/Subroutines

subroutine mod_fld::fld_params_read (files)
 public methods these are called in mod_rhd_phys More...
 
subroutine, public mod_fld::fld_init (He_abundance, rhd_radiation_diffusion, rhd_gamma)
 Initialising FLD-module: Read opacities Initialise Multigrid adimensionalise kappa Add extra variables to w-array, flux, kappa, eddington Tensor Lambda and R ... More...
 
subroutine, public mod_fld::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 the radiation force More...
 
subroutine, public mod_fld::fld_radforce_get_dt (w, ixIL, ixOL, dtnew, dxD, x)
 
subroutine, public mod_fld::get_fld_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 the energy exchange between gas and radiation More...
 
subroutine, public mod_fld::fld_get_opacity (w, x, ixIL, ixOL, fld_kappa)
 Sets the opacity in the w-array by calling mod_opal_opacity. More...
 
subroutine, public mod_fld::fld_get_fluxlimiter (w, x, ixIL, ixOL, fld_lambda, fld_R)
 Calculate fld flux limiter This subroutine calculates flux limiter lambda using the prescription stored in fld_fluxlimiter. It also calculates the ratio of radiation scaleheight and mean free path. More...
 
subroutine, public mod_fld::fld_get_radflux (w, x, ixIL, ixOL, rad_flux)
 Calculate Radiation Flux stores radiation flux in w-array. More...
 
subroutine mod_fld::fld_get_eddington (w, x, ixIL, ixOL, eddington_tensor)
 Calculate Eddington-tensor Stores Eddington-tensor in w-array. More...
 
subroutine, public mod_fld::fld_get_radpress (w, x, ixIL, ixOL, rad_pressure)
 Calculate Radiation Pressure Returns Radiation Pressure as tensor. More...
 
subroutine mod_fld::diffuse_e_rad_mg (dtfactor, qdt, qtC, psa, psb)
 Calling all subroutines to perform the multigrid method Communicates rad_e and diff_coeff to multigrid library. More...
 
subroutine mod_fld::evaluate_e_rad_mg (qtC, psa)
 inplace update of psa==>F_im(psa) More...
 
subroutine mod_fld::put_diffterm_onegrid (ixIL, ixOL, w)
 inplace update of psa==>F_im(psa) More...
 
subroutine mod_fld::fld_get_diffcoef_central (w, wCT, x, ixIL, ixOL)
 Calculates cell-centered diffusion coefficient to be used in multigrid. More...
 
subroutine mod_fld::update_diffcoeff (psa)
 
subroutine mod_fld::fld_smooth_diffcoef (w, ixIL, ixOL)
 Use running average on Diffusion coefficient. More...
 
subroutine mod_fld::energy_interaction (w, wCT, x, ixIL, ixOL)
 This subroutine calculates the radiation heating, radiation cooling and photon tiring using an implicit scheme. These sourceterms are applied using the root-finding of a 4th order polynomial This routine loops over every cell in the domain and computes the coefficients of the polynomials in every cell. More...
 
subroutine mod_fld::bisection_method (e_gas, E_rad, c0, c1)
 Find the root of the 4th degree polynomial using the bisection method. More...
 
subroutine mod_fld::newton_method (e_gas, E_rad, c0, c1)
 Find the root of the 4th degree polynomial using the Newton-Ralphson method. More...
 
subroutine mod_fld::halley_method (e_gas, E_rad, c0, c1)
 Find the root of the 4th degree polynomial using the Halley method. More...
 
double precision function mod_fld::polynomial_bisection (e_gas, c0, c1)
 Evaluate polynomial at argument e_gas. More...
 
double precision function mod_fld::dpolynomial_bisection (e_gas, c0, c1)
 Evaluate first derivative of polynomial at argument e_gas. More...
 
double precision function mod_fld::ddpolynomial_bisection (e_gas, c0, c1)
 Evaluate second derivative of polynomial at argument e_gas. More...
 
subroutine mod_fld::gradiento (q, x, ixIL, ixOL, idir, gradq, n)
 Calculate gradient of a scalar q within ixL in direction idir difference with gradient is gradq(ixO^S), NOT gradq(ixI^S) More...
 

Variables

logical mod_fld::fld_eint_split = .false.
 source split for energy interact and radforce: More...
 
logical mod_fld::fld_radforce_split = .false.
 
double precision, public mod_fld::fld_kappa0 = 0.34d0
 Opacity per unit of unit_density. More...
 
double precision, public mod_fld::fld_mu = 0.6d0
 mean particle mass More...
 
double precision, public mod_fld::fld_bisect_tol = 1.d-4
 Tolerance for bisection method for Energy sourceterms This is a percentage of the minimum of gas- and radiation energy. More...
 
double precision, public mod_fld::fld_diff_tol = 1.d-4
 Tolerance for adi method for radiative Energy diffusion. More...
 
double precision, public mod_fld::diff_crit
 Number for splitting the diffusion module. More...
 
character(len=8) mod_fld::fld_opacity_law = 'const'
 Use constant Opacity? More...
 
character(len=50) mod_fld::fld_opal_table = 'Y09800'
 
character(len=16) mod_fld::fld_fluxlimiter = 'Pomraning'
 Diffusion limit lambda = 0.33. More...
 
integer mod_fld::i_diff_mg
 diffusion coefficient for multigrid method More...
 
character(len=8) mod_fld::fld_diff_scheme = 'mg'
 Which method to solve diffusion part. More...
 
character(len=8) mod_fld::fld_interaction_method = 'Halley'
 Which method to find the root for the energy interaction polynomial. More...
 
logical mod_fld::diff_coef_filter = .false.
 Take running average for Diffusion coefficient. More...
 
integer mod_fld::size_d_filter = 1
 
logical mod_fld::flux_lim_filter = .false.
 Take a running average over the fluxlimiter. More...
 
integer mod_fld::size_l_filter = 1
 
logical mod_fld::lineforce_opacities = .false.
 Use or dont use lineforce opacities. More...
 
logical mod_fld::diffcrash_resume = .true.
 Resume run when multigrid returns error. More...
 
integer, dimension(:), allocatable, public mod_fld::i_opf
 Index for Flux weighted opacities. More...
 
double precision mod_fld::dt_diff = 0.d0
 running timestep for diffusion solver, initialised as zero More...