chombo-discharge
|
Radiative tranfer equation solver in the SP1 (diffusion) approximation. More...
#include <CD_EddingtonSP1.H>
Public Member Functions | |
EddingtonSP1 () | |
Weak constructor. | |
EddingtonSP1 (const EddingtonSP1 &a_other)=delete | |
Disallowed copy construction. | |
EddingtonSP1 (const EddingtonSP1 &&a_other)=delete | |
Disallowed move construction. | |
virtual | ~EddingtonSP1 () |
Destructor. | |
EddingtonSP1 & | operator= (const EddingtonSP1 &a_other)=delete |
Disallowed copy assignment. | |
EddingtonSP1 & | operator= (const EddingtonSP1 &&a_other)=delete |
Disallowed move assignment. | |
virtual bool | advance (const Real a_dt, EBAMRCellData &a_phi, const EBAMRCellData &a_source, const bool a_zeroPhi=false) override |
Advance RTE onto state a_phi. More... | |
virtual void | advanceEuler (EBAMRCellData &a_phi, const EBAMRCellData &a_source, const Real a_dt, const bool a_zeroPhi) noexcept |
Advance using the Euler rule. More... | |
virtual void | parseOptions () override |
Parse class options. | |
virtual void | parseRuntimeOptions () override |
Parse class options. | |
virtual void | setDomainSideBcFunction (const int a_dir, const Side::LoHiSide a_side, const EddingtonSP1DomainBc::BcFunction &a_function) |
Set the boundary condition function on a domain side. More... | |
virtual void | allocate () override |
Allocate internal storage. | |
virtual void | preRegrid (const int a_base, const int a_oldFinestLevel) override |
Cache state. More... | |
virtual void | deallocate () override |
Deallocate internal storage. | |
virtual void | regrid (const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) override |
Regrid function for this class. More... | |
virtual void | registerOperators () override |
Register operators. | |
virtual void | computeBoundaryFlux (EBAMRIVData &a_ebFlux, const EBAMRCellData &a_phi) override |
Compute the boundary flux. For Eddington, the boundary flux is = c*phi/2. More... | |
virtual void | computeDomainFlux (EBAMRIFData &a_domainflux, const EBAMRCellData &a_phi) override |
Compute the domain flux. For Eddington, the domain flux is = c*phi/2. More... | |
virtual void | computeFlux (EBAMRCellData &a_flux, const EBAMRCellData &a_phi) override |
Compute the flux. For Eddington, the flux is F = -c/(3*kappa)*grad(phi) More... | |
virtual void | computeDensity (EBAMRCellData &a_isotropic, const EBAMRCellData &a_phi) override |
Get isotropic part. More... | |
virtual void | writePlotFile () override |
Write plot file. | |
Public Member Functions inherited from RtSolver | |
RtSolver () | |
Constructor. | |
virtual | ~RtSolver () |
Constructor (does nothing) | |
virtual std::string | getName () |
Get solver name. | |
virtual const std::string | getRealm () const |
Get the realm where the solver lives. | |
virtual bool | advance (const Real a_dt, const bool a_zeroPhi=false) |
Advance equation one time step. More... | |
virtual bool | advance (const Real a_dt, EBAMRCellData &a_phi, const bool a_zeroPhi=false) |
Advance method. Advances one time step. More... | |
virtual void | setRealm (const std::string a_realm) |
Set realm where this solver lives. More... | |
virtual void | setRtSpecies (const RefCountedPtr< RtSpecies > &a_species) |
Set the radiative transfer species (RtSpecies) More... | |
virtual void | setComputationalGeometry (const RefCountedPtr< ComputationalGeometry > a_computationalGeometry) |
Set computational geometry. More... | |
virtual void | computeLoads (Vector< long long > &a_loads, const DisjointBoxLayout &a_dbl, const int a_level) const noexcept |
Get computational loads for a specific grid level. More... | |
virtual void | setAmr (const RefCountedPtr< AmrMesh > &a_amr) |
Set the amr object. More... | |
virtual void | setPhase (phase::which_phase a_phase=phase::gas) |
Set phase. More... | |
virtual void | setVerbosity (const int a_verbosity) |
Set verbosity. More... | |
virtual void | setTime (const int a_step, const Real a_time, const Real a_dt) |
Set the time for this solver. More... | |
virtual void | setStationary (const bool a_stationary) |
Set stationary solver or not. More... | |
virtual void | sanityCheck () |
Sanity check. | |
virtual bool | isStationary () |
Check if solver is stationary. | |
virtual void | initialData () |
Fill solver with initial data. By default, this sets internal data to zero. More... | |
virtual void | setSource (const EBAMRCellData &a_source) |
Set source term. More... | |
virtual void | setSource (const Real a_source) |
Set source. More... | |
virtual void | setSource (const std::function< Real(const RealVect a_pos)> a_source) |
Set source. More... | |
virtual int | getNumberOfPlotVariables () const |
Get number of output fields. More... | |
virtual Vector< std::string > | getPlotVariableNames () const |
Get output plot names. | |
virtual void | writePlotData (LevelData< EBCellFAB > &a_output, int &a_comp, const std::string a_outputRealm, const int a_level) const noexcept |
Write output data to a_output. More... | |
virtual Real | getTime () const |
Get current time. More... | |
virtual phase::which_phase | getPhase () |
Get the RTE phase. More... | |
virtual EBAMRCellData & | getPhi () |
Get solver state. More... | |
virtual EBAMRCellData & | getSource () |
Get multifluid source. More... | |
virtual EBAMRFluxData & | getKappa () |
Get the absorption length. More... | |
virtual EBAMRIVData & | getKappaEb () |
Get the absorption coefficient on irregular EB faces. More... | |
virtual RefCountedPtr< RtSpecies > & | getSpecies () |
Get species. | |
Static Public Member Functions | |
static Real | s_defaultDomainBcFunction (const RealVect a_position, const Real a_time) |
Default function for space-time dependence of domain boundary conditions. More... | |
Protected Types | |
enum class | BottomSolverType { Simple , BiCGStab , GMRES } |
Enum class for supported bottom solvers in multigrid. | |
enum class | MultigridType { VCycle , WCycle } |
Enum for multigrid cycle types. | |
enum class | EBBCType { Dirichlet , Neumann , Larsen } |
Enum for boundary condition types on EBs. | |
Protected Member Functions | |
virtual void | setupSolver () |
Set up geometric multigrid. | |
virtual void | setHelmholtzCoefficients () |
Set multigrid coefficients. | |
virtual void | setHelmholtzCoefficientsBox (EBCellFAB &a_helmAco, EBFluxFAB &a_helmBco, BaseIVFAB< Real > &a_helmBcoIrreg, const int a_lvl, const DataIndex &a_dit) |
Set EBHelmholtzOp A- and B-coefficients. More... | |
virtual void | setupHelmholtzFactory () |
Set the operator factory. | |
virtual void | setupMultigrid () |
Set the GMG solver. | |
virtual std::string | makeBcString (const int a_dir, const Side::LoHiSide a_side) const |
Make domain bc string. More... | |
virtual EddingtonSP1DomainBc::BcType | parseBcString (const std::string a_str) const |
Returns BC type based on string. More... | |
virtual void | parseDomainBC () |
Parse domain BC settings. | |
virtual void | parseEBBC () |
Parse domain BC settings. | |
virtual void | parseStationary () |
Parse whether or not this is a stationary solver. | |
virtual void | parsePlotVariables () |
Parse plot variables. | |
virtual void | parseMultigridSettings () |
Parse multigrid settings. | |
virtual void | parseKappaScale () |
Parse kappa-scaling or not. | |
virtual void | parseReflection () |
Parse reflection coefficients for Robin bcs. | |
virtual void | parseRegridSlopes () |
Parse whether or not to use slopes when regridding. | |
virtual void | setDefaultDomainBcFunctions () |
Set default domain BC functions. | |
Protected Member Functions inherited from RtSolver | |
void | setEbIndexSpace (const RefCountedPtr< EBIndexSpace > &a_ebis) |
Set ebis. | |
void | parseVerbosity () noexcept |
Parse verbosity. | |
virtual void | writeData (LevelData< EBCellFAB > &a_output, int &a_comp, const EBAMRCellData &a_data, const std::string a_outputRealm, const int a_level, const bool a_interpToCentroids, const bool a_interpGhost) const noexcept |
Write data to output. Convenience function. More... | |
Protected Attributes | |
EBHelmholtzOp::Smoother | m_multigridRelaxMethod |
Relaxation type for gmg. | |
MultigridType | m_multigridType |
GMG multigrid type. | |
bool | m_isSolverSetup |
Needs setup. | |
bool | m_kappaScale |
Use kappa scaling for source or not. | |
bool | m_regridSlopes |
Use slopes when regridding (or not) | |
int | m_multigridVerbosity |
Verbosity for geometric multigrid. | |
int | m_multigridPreSmooth |
Number of smoothings before averaging. | |
int | m_multigridPostSmooth |
Number of smoothings before averaging. | |
int | m_multigridBottomSmooth |
Number of smoothing before bottom solver. | |
int | m_multigridMaxIterations |
Maximum number of iterations. | |
int | m_multigridMinIterations |
Minimum number of iterations. | |
int | m_multigridBcOrder |
Multigrid EBBC order (only relevant for Dirichlet) | |
int | m_multigridBcWeight |
Multigrid EBBC weight (only relevant for Dirichlet) | |
BottomSolverType | m_bottomSolverType |
Bottom solver type. | |
int | m_numSmoothingsForSimpleSolver |
Number of smoothing for bottom solver. | |
int | m_minCellsBottom |
Set bottom drop depth. | |
Real | m_multigridExitTolerance |
Real | m_multigridExitHang |
Real | m_reflectCoefOne |
Reflection coefficient. | |
Real | m_reflectCoefTwo |
Reflection coefficient. | |
RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > | m_multigridSolver |
Geometric multigrid solver. | |
BiCGStabSolver< LevelData< EBCellFAB > > | m_bicgstab |
Conjugate gradient solver bottom MG level. | |
GMRESSolver< LevelData< EBCellFAB > > | m_gmres |
GMRES solver. | |
RefCountedPtr< EBHelmholtzOpFactory > | m_helmholtzOpFactory |
Operator factory. | |
EBSimpleSolver | m_simpleSolver |
multi-fluid simple solver | |
EBAMRCellData | m_cacheSrc |
For regridding the source term. This is needed when doing a stationary solve. | |
EBAMRCellData | m_resid |
Multigrid residue. | |
EBAMRCellData | m_helmAco |
a-coefficient | |
EBAMRFluxData | m_helmBco |
b-coefficient | |
EBAMRIVData | m_helmBcoIrreg |
b-coefficient | |
EddingtonSP1DomainBc | m_domainBc |
Wrapper calss. | |
std::map< EddingtonSP1DomainBc::DomainSide, EddingtonSP1DomainBc::BcFunction > | m_domainBcFunctions |
Actual functions on domain edges (faces). | |
std::pair< EBBCType, Real > | m_ebbc |
Associated boundary condition on the embedded boundary. | |
Protected Attributes inherited from RtSolver | |
Location::Cell | m_dataLocation |
Data location. | |
std::string | m_realm |
Realm where this solver lives. | |
RefCountedPtr< EBIndexSpace > | m_ebis |
EBIndexSpace for this solver. | |
RefCountedPtr< RtSpecies > | m_rtSpecies |
Radiative transfer species (contains meta-information like initial conditions) | |
RefCountedPtr< ComputationalGeometry > | m_computationalGeometry |
Computational geometry. | |
RefCountedPtr< AmrMesh > | m_amr |
AMR; needed for grid stuff. | |
phase::which_phase | m_phase |
Phase. | |
std::string | m_name |
Name for this solver. | |
std::string | m_className |
Class name – needed because inherited classes will be named different. | |
EBAMRCellData | m_cachePhi |
Cached state used for regridding. | |
EBAMRCellData | m_phi |
Internal state. More... | |
EBAMRCellData | m_source |
Source term. More... | |
EBAMRFluxData | m_kappa |
Absorption coefficient. | |
EBAMRIVData | m_kappaEB |
Absorption coefficient on EB faces. | |
Real | m_time |
Time. | |
Real | m_dt |
Time increment. | |
bool | m_stationary |
Stationary solver or not. | |
bool | m_plotPhi |
Output state. | |
bool | m_plotSource |
Output source term. | |
int | m_verbosity |
int | m_timeStep |
Time step. | |
Static Protected Attributes | |
static constexpr Real | m_alpha = 1.0 |
Write checkpoint data into HDF5 file. More... | |
static constexpr Real | m_beta = -1.0 |
Beta-coefficient for Helmholtz operator. | |
Static Protected Attributes inherited from RtSolver | |
static constexpr int | m_comp = 0 |
Default component that we solve for. | |
static constexpr int | m_nComp = 1 |
Default number of components. | |
Radiative tranfer equation solver in the SP1 (diffusion) approximation.
|
overridevirtual |
Advance RTE onto state a_phi.
[in] | a_dt | Time step |
[in,out] | a_phi | RTE solution |
[in] | a_source | Source term |
[in] | a_zeroPhi | Set phi to zero in initial guess for multigrid solve |
Implements RtSolver.
|
virtualnoexcept |
Advance using the Euler rule.
This implementation assumes that the incoming previous solution has not been weighted by the volume fraction, and that the incoming source term a_source HAS been weighted by the volume fraction.
[in,out] | a_newPhi | On exit, contains solution at time t + dt. On entry, containes solution at tim et |
[in] | a_source | Source term. Should be weighted by the volume fraction. |
[in] | a_dt | Time step. |
[in] | a_zeroPhi | Set phi to zero before multigrid solving |
|
overridevirtual |
Compute the boundary flux. For Eddington, the boundary flux is = c*phi/2.
[out] | a_ebFlux | Flux on the EBs |
[in] | a_phi | Isotropic part of the RTE solution |
Implements RtSolver.
|
overridevirtual |
Get isotropic part.
[out] | a_isotropic | Isotropic part of RTE solution |
[in] | a_phi | RTE solution |
For Eddington, only the isotropic part is solved for so this routine just copies.
Implements RtSolver.
|
overridevirtual |
Compute the domain flux. For Eddington, the domain flux is = c*phi/2.
[out] | a_domainFlux | Flux on the domain sides |
[in] | a_phi | Isotropic part of the RTE solution |
Implements RtSolver.
|
overridevirtual |
Compute the flux. For Eddington, the flux is F = -c/(3*kappa)*grad(phi)
[out] | a_flux | Cell-centered flux |
[in] | a_phi | Isotropic part of RTE solution |
Implements RtSolver.
|
protectedvirtual |
Make domain bc string.
[in] | a_dir | Coordinate direction |
[in] | a_side | Side |
Used for pasing boundary conditions from input script.
|
protectedvirtual |
Returns BC type based on string.
[in] | a_str | Boundary condition string. Must be "dirichlet <number>", "neumann <number>", "robin <number>", "dirichlet_custom", "neumann_custom", or "robin_custom". |
|
overridevirtual |
Cache state.
[in] | a_base | Coarsest level that did not change during regrid |
[in] | a_oldFinestLevel | Finest level before the regrid |
Implements RtSolver.
|
overridevirtual |
Regrid function for this class.
[in] | a_lmin | Coarsest level that changed during regrid |
[in] | a_oldFinestLevel | Finest level before the regrid |
[in] | a_newFinestLevel | Finest level after the regrid |
Implements RtSolver.
|
static |
Default function for space-time dependence of domain boundary conditions.
This is a utility function
[in] | a_position | Position in space |
[in] | a_time | Time |
|
virtual |
Set the boundary condition function on a domain side.
[in] | a_dir | Coordinate direction. |
[in] | a_side | Side (low/high) |
[in] | a_function | Boundary condition function. |
This sets a boundary condition for a particular domain side. The user must also specify how to use this BC in the input script.
|
protectedvirtual |
Set EBHelmholtzOp A- and B-coefficients.
For the B-coefficient, this also fills one of the tangential ghost faces.
[out] | a_helmAco | Helmholtz A-coefficient on on cell centers |
[out] | a_helmBco | Helmholtz B-coefficient on on face centers |
[out] | a_helmBcoIrreg | Helmholtz B-coefficient on EB faces |
[in] | a_lvl | Grid level |
[in] | a_dit | Grid index |
|
staticconstexprprotected |
Write checkpoint data into HDF5 file.
[out] | a_handle | HDF5 file |
[in] | a_level | Grid level |
Read checkpoint data from handle
[out] | a_handle | HDF5 file |
[in] | a_level | Grid level |
alpha-coefficient for Helmholtz operator