12 #ifndef CD_EddingtonSP1_H
13 #define CD_EddingtonSP1_H
19 #include <AMRMultiGrid.H>
21 #include <EBBackwardEuler.H>
22 #include <EBSimpleSolver.H>
23 #include <BiCGStabSolver.H>
24 #include <GMRESSolver.H>
29 #include <CD_EddingtonSP1DomainBc.H>
30 #include <CD_NamespaceHeader.H>
124 const Side::LoHiSide a_side,
139 preRegrid(
const int a_base,
const int a_oldFinestLevel)
override;
154 regrid(
const int a_lmin,
const int a_oldFinestLevel,
const int a_newFinestLevel)
override;
208 writeCheckpointLevel(HDF5Handle& a_handle,
const int a_level)
const override;
218 readCheckpointLevel(HDF5Handle& a_handle,
const int a_level)
override;
344 Real m_multigridExitTolerance;
349 Real m_multigridExitHang;
449 EBFluxFAB& a_helmBco,
450 BaseIVFAB<Real>& a_helmBcoIrreg,
452 const DataIndex& a_dit);
474 makeBcString(
const int a_dir,
const Side::LoHiSide a_side)
const;
538 #include <CD_NamespaceFooter.H>
Declaration of a factory class for making Poisson operators for multigrid.
Abstract parent class for various radiative transfer solvers.
Smoother
Relaxation method for the operators.
Definition: CD_EBHelmholtzOp.H:47
Class which maps boundary condition types to a side and direction.
Definition: CD_EddingtonSP1DomainBc.H:31
std::function< Real(const RealVect a_position, const Real a_time)> BcFunction
Function which maps f(R^3,t) : R. Used for setting the associated value and boundary condition type.
Definition: CD_EddingtonSP1DomainBc.H:46
BcType
Boundary condition type.
Definition: CD_EddingtonSP1DomainBc.H:37
Radiative tranfer equation solver in the SP1 (diffusion) approximation.
Definition: CD_EddingtonSP1.H:36
virtual void setupSolver()
Set up geometric multigrid.
Definition: CD_EddingtonSP1.cpp:739
virtual void computeBoundaryFlux(EBAMRIVData &a_ebFlux, const EBAMRCellData &a_phi) override
Compute the boundary flux. For Eddington, the boundary flux is = c*phi/2.
Definition: CD_EddingtonSP1.cpp:1047
EBAMRCellData m_resid
Multigrid residue.
Definition: CD_EddingtonSP1.H:394
RefCountedPtr< EBHelmholtzOpFactory > m_helmholtzOpFactory
Operator factory.
Definition: CD_EddingtonSP1.H:379
GMRESSolver< LevelData< EBCellFAB > > m_gmres
GMRES solver.
Definition: CD_EddingtonSP1.H:374
virtual void regrid(const int a_lmin, const int a_oldFinestLevel, const int a_newFinestLevel) override
Regrid function for this class.
Definition: CD_EddingtonSP1.cpp:527
MultigridType m_multigridType
GMG multigrid type.
Definition: CD_EddingtonSP1.H:269
int m_multigridBcWeight
Multigrid EBBC weight (only relevant for Dirichlet)
Definition: CD_EddingtonSP1.H:324
virtual void advanceEuler(EBAMRCellData &a_phi, const EBAMRCellData &a_source, const Real a_dt, const bool a_zeroPhi) noexcept
Advance using the Euler rule.
Definition: CD_EddingtonSP1.cpp:673
virtual void parseReflection()
Parse reflection coefficients for Robin bcs.
Definition: CD_EddingtonSP1.cpp:442
EBHelmholtzOp::Smoother m_multigridRelaxMethod
Relaxation type for gmg.
Definition: CD_EddingtonSP1.H:264
BottomSolverType m_bottomSolverType
Bottom solver type.
Definition: CD_EddingtonSP1.H:329
EBAMRFluxData m_helmBco
b-coefficient
Definition: CD_EddingtonSP1.H:404
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)
Definition: CD_EddingtonSP1.cpp:1138
virtual void setupMultigrid()
Set the GMG solver.
Definition: CD_EddingtonSP1.cpp:973
virtual void writePlotFile() override
Write plot file.
Definition: CD_EddingtonSP1.cpp:1180
int m_multigridPreSmooth
Number of smoothings before averaging.
Definition: CD_EddingtonSP1.H:294
virtual void setupHelmholtzFactory()
Set the operator factory.
Definition: CD_EddingtonSP1.cpp:895
virtual void setHelmholtzCoefficients()
Set multigrid coefficients.
Definition: CD_EddingtonSP1.cpp:754
Real m_reflectCoefTwo
Reflection coefficient.
Definition: CD_EddingtonSP1.H:359
virtual void setDefaultDomainBcFunctions()
Set default domain BC functions.
Definition: CD_EddingtonSP1.cpp:95
virtual void deallocate() override
Deallocate internal storage.
Definition: CD_EddingtonSP1.cpp:511
virtual ~EddingtonSP1()
Destructor.
Definition: CD_EddingtonSP1.cpp:57
int m_minCellsBottom
Set bottom drop depth.
Definition: CD_EddingtonSP1.H:339
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.
Definition: CD_EddingtonSP1.cpp:580
virtual void parsePlotVariables()
Parse plot variables.
Definition: CD_EddingtonSP1.cpp:327
int m_multigridMaxIterations
Maximum number of iterations.
Definition: CD_EddingtonSP1.H:309
virtual void preRegrid(const int a_base, const int a_oldFinestLevel) override
Cache state.
Definition: CD_EddingtonSP1.cpp:472
static Real s_defaultDomainBcFunction(const RealVect a_position, const Real a_time)
Default function for space-time dependence of domain boundary conditions.
Definition: CD_EddingtonSP1.cpp:36
virtual std::string makeBcString(const int a_dir, const Side::LoHiSide a_side) const
Make domain bc string.
Definition: CD_EddingtonSP1.cpp:126
EddingtonSP1(const EddingtonSP1 &&a_other)=delete
Disallowed move construction.
virtual void parseEBBC()
Parse domain BC settings.
Definition: CD_EddingtonSP1.cpp:270
std::pair< EBBCType, Real > m_ebbc
Associated boundary condition on the embedded boundary.
Definition: CD_EddingtonSP1.H:424
int m_numSmoothingsForSimpleSolver
Number of smoothing for bottom solver.
Definition: CD_EddingtonSP1.H:334
EBSimpleSolver m_simpleSolver
multi-fluid simple solver
Definition: CD_EddingtonSP1.H:384
EddingtonSP1 & operator=(const EddingtonSP1 &&a_other)=delete
Disallowed move assignment.
virtual void parseDomainBC()
Parse domain BC settings.
Definition: CD_EddingtonSP1.cpp:186
bool m_kappaScale
Use kappa scaling for source or not.
Definition: CD_EddingtonSP1.H:279
static constexpr Real m_beta
Beta-coefficient for Helmholtz operator.
Definition: CD_EddingtonSP1.H:230
EBAMRIVData m_helmBcoIrreg
b-coefficient
Definition: CD_EddingtonSP1.H:409
bool m_isSolverSetup
Needs setup.
Definition: CD_EddingtonSP1.H:274
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.
Definition: CD_EddingtonSP1.cpp:111
EBAMRCellData m_cacheSrc
For regridding the source term. This is needed when doing a stationary solve.
Definition: CD_EddingtonSP1.H:389
Real m_reflectCoefOne
Reflection coefficient.
Definition: CD_EddingtonSP1.H:354
RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > m_multigridSolver
Geometric multigrid solver.
Definition: CD_EddingtonSP1.H:364
virtual void computeDensity(EBAMRCellData &a_isotropic, const EBAMRCellData &a_phi) override
Get isotropic part.
Definition: CD_EddingtonSP1.cpp:1167
virtual void computeDomainFlux(EBAMRIFData &a_domainflux, const EBAMRCellData &a_phi) override
Compute the domain flux. For Eddington, the domain flux is = c*phi/2.
Definition: CD_EddingtonSP1.cpp:1071
EddingtonSP1DomainBc m_domainBc
Wrapper calss.
Definition: CD_EddingtonSP1.H:414
BiCGStabSolver< LevelData< EBCellFAB > > m_bicgstab
Conjugate gradient solver bottom MG level.
Definition: CD_EddingtonSP1.H:369
int m_multigridBcOrder
Multigrid EBBC order (only relevant for Dirichlet)
Definition: CD_EddingtonSP1.H:319
int m_multigridVerbosity
Verbosity for geometric multigrid.
Definition: CD_EddingtonSP1.H:289
EddingtonSP1()
Weak constructor.
Definition: CD_EddingtonSP1.cpp:41
EddingtonSP1(const EddingtonSP1 &a_other)=delete
Disallowed copy construction.
EBAMRCellData m_helmAco
a-coefficient
Definition: CD_EddingtonSP1.H:399
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.
Definition: CD_EddingtonSP1.cpp:783
int m_multigridMinIterations
Minimum number of iterations.
Definition: CD_EddingtonSP1.H:314
virtual void parseKappaScale()
Parse kappa-scaling or not.
Definition: CD_EddingtonSP1.cpp:351
static constexpr Real m_alpha
Write checkpoint data into HDF5 file.
Definition: CD_EddingtonSP1.H:225
EddingtonSP1 & operator=(const EddingtonSP1 &a_other)=delete
Disallowed copy assignment.
BottomSolverType
Enum class for supported bottom solvers in multigrid.
Definition: CD_EddingtonSP1.H:236
virtual void parseRuntimeOptions() override
Parse class options.
Definition: CD_EddingtonSP1.cpp:80
virtual EddingtonSP1DomainBc::BcType parseBcString(const std::string a_str) const
Returns BC type based on string.
Definition: CD_EddingtonSP1.cpp:160
virtual void parseOptions() override
Parse class options.
Definition: CD_EddingtonSP1.cpp:61
virtual void allocate() override
Allocate internal storage.
Definition: CD_EddingtonSP1.cpp:489
EBBCType
Enum for boundary condition types on EBs.
Definition: CD_EddingtonSP1.H:255
virtual void registerOperators() override
Register operators.
Definition: CD_EddingtonSP1.cpp:558
bool m_regridSlopes
Use slopes when regridding (or not)
Definition: CD_EddingtonSP1.H:284
std::map< EddingtonSP1DomainBc::DomainSide, EddingtonSP1DomainBc::BcFunction > m_domainBcFunctions
Actual functions on domain edges (faces).
Definition: CD_EddingtonSP1.H:419
virtual void parseStationary()
Parse whether or not this is a stationary solver.
Definition: CD_EddingtonSP1.cpp:313
int m_multigridBottomSmooth
Number of smoothing before bottom solver.
Definition: CD_EddingtonSP1.H:304
MultigridType
Enum for multigrid cycle types.
Definition: CD_EddingtonSP1.H:246
int m_multigridPostSmooth
Number of smoothings before averaging.
Definition: CD_EddingtonSP1.H:299
virtual void parseMultigridSettings()
Parse multigrid settings.
Definition: CD_EddingtonSP1.cpp:359
virtual void parseRegridSlopes()
Parse whether or not to use slopes when regridding.
Definition: CD_EddingtonSP1.cpp:459
Abstract RTE solver class for doing various kinds of radiative transfer equations....
Definition: CD_RtSolver.H:30