chombo-discharge
CD_CdrMultigrid.H
Go to the documentation of this file.
1 /* chombo-discharge
2  * Copyright © 2021 SINTEF Energy Research.
3  * Please refer to Copyright.txt and LICENSE in the chombo-discharge root directory.
4  */
5 
12 #ifndef CD_CdrMultigrid_H
13 #define CD_CdrMultigrid_H
14 
15 // Std includes
16 #include <random>
17 #include <time.h>
18 #include <chrono>
19 
20 // Chombo includes
21 #include <AMRMultiGrid.H>
22 #include <BiCGStabSolver.H>
23 #include <EBSimpleSolver.H>
24 #include <GMRESSolver.H>
25 
26 // Our includes
28 #include <CD_CdrSolver.H>
29 #include <CD_NamespaceHeader.H>
30 
34 class CdrMultigrid : public CdrSolver
35 {
36 protected:
40  enum class BottomSolverType
41  {
42  Simple,
43  BiCGStab,
44  GMRES
45  };
46 
50  enum class MultigridType
51  {
52  VCycle,
53  WCycle,
54  };
55 
56 public:
60  CdrMultigrid();
61 
65  virtual ~CdrMultigrid();
66 
70  virtual void
71  parseOptions() override = 0;
72 
76  virtual void
77  parseRuntimeOptions() override = 0;
78 
82  virtual void
83  registerOperators() override;
84 
88  virtual void
89  allocate() override;
90 
97  virtual void
98  preRegrid(const int a_lbase, const int a_oldFinestLevel) override;
99 
109  virtual void
110  computeDivJ(EBAMRCellData& a_divJ,
111  EBAMRCellData& a_phi,
112  const Real a_extrapDt,
113  const bool a_conservativeOnly,
114  const bool a_ebFlux,
115  const bool a_domainFlux) override final;
116 
127  virtual void
128  computeDivF(EBAMRCellData& a_divF,
129  EBAMRCellData& a_phi,
130  const Real a_extrapDt,
131  const bool a_conservativeOnly,
132  const bool a_ebFlux,
133  const bool a_domainFlux) override final;
134 
144  virtual void
145  computeDivD(EBAMRCellData& a_divD,
146  EBAMRCellData& a_phi,
147  const bool a_conservativeOnly,
148  const bool a_ebFlux,
149  const bool a_domainFlux) override final;
150 
160  virtual void
161  advanceEuler(EBAMRCellData& a_newPhi,
162  const EBAMRCellData& a_oldPhi,
163  const EBAMRCellData& a_source,
164  const Real a_dt) override final;
165 
175  virtual void
177  const EBAMRCellData& a_oldPhi,
178  const EBAMRCellData& a_source,
179  const Real a_dt) override final;
180 
181 protected:
186 
191 
195  RefCountedPtr<AMRMultiGrid<LevelData<EBCellFAB>>> m_multigridSolver;
196 
200  RefCountedPtr<EBHelmholtzOpFactory> m_helmholtzOpFactory;
201 
205  BiCGStabSolver<LevelData<EBCellFAB>> m_bicgstab;
206 
210  EBSimpleSolver m_simpleSolver;
211 
215  GMRESSolver<LevelData<EBCellFAB>> m_gmres;
216 
221 
226 
231 
236 
241 
246 
251 
256 
261 
266 
271 
276 
281 
286 
290  virtual void
291  advectToFaces(EBAMRFluxData& a_facePhi, const EBAMRCellData& a_phi, const Real a_extrapDt) override = 0;
292 
296  virtual void
298 
302  virtual void
304 
308  virtual void
309  setupMultigrid();
310 
315  virtual void
317 
323  virtual void
324  resetAlphaAndBeta(const Real a_alpha, const Real a_beta);
325 
332  virtual void
333  computeKappaLphi(EBAMRCellData& a_kappaLphi, const EBAMRCellData& a_phi);
334 
338  virtual void
340 };
341 
342 #include <CD_NamespaceFooter.H>
343 
344 #endif
Declaration of an abstract class for solving scalar convection-diffusion-reaction problems.
Declaration of a factory class for making Poisson operators for multigrid.
Extension class of CdrSolver that uses multigrid for diffusion part. Can also solve for stochastic di...
Definition: CD_CdrMultigrid.H:35
virtual void registerOperators() override
Register operator.
Definition: CD_CdrMultigrid.cpp:38
MultigridType
Enum for multigrid cycle types.
Definition: CD_CdrMultigrid.H:51
virtual void setMultigridSolverCoefficients()
Set the multigrid solver coefficients.
Definition: CD_CdrMultigrid.cpp:96
virtual void advanceEuler(EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const EBAMRCellData &a_source, const Real a_dt) override final
Implicit diffusion Euler advance with source term.
Definition: CD_CdrMultigrid.cpp:164
CdrMultigrid()
Constructor.
Definition: CD_CdrMultigrid.cpp:24
EBSimpleSolver m_simpleSolver
Simple solver.
Definition: CD_CdrMultigrid.H:210
MultigridType m_multigridType
GMG multigrid type.
Definition: CD_CdrMultigrid.H:190
RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > m_multigridSolver
Geometric multigrid solver.
Definition: CD_CdrMultigrid.H:195
virtual void computeDivD(EBAMRCellData &a_divD, EBAMRCellData &a_phi, const bool a_conservativeOnly, const bool a_ebFlux, const bool a_domainFlux) override final
Compute div(D*grad(phi)) explicitly.
Definition: CD_CdrMultigrid.cpp:636
bool m_hasMultigridSolver
Is multigrid solver defined or not?
Definition: CD_CdrMultigrid.H:220
BottomSolverType
Enum class for supported bottom solvers in multigrid.
Definition: CD_CdrMultigrid.H:41
int m_multigridVerbosity
Verbosity for geometric multigrid.
Definition: CD_CdrMultigrid.H:235
int m_numSmoothingsForSimpleSolver
Number of smoothing for bottom solver.
Definition: CD_CdrMultigrid.H:270
virtual void computeKappaLphi(EBAMRCellData &a_kappaLphi, const EBAMRCellData &a_phi)
Compute kappa * L(phi) using the multigrid operator. This is different from computeDivD.
Definition: CD_CdrMultigrid.cpp:141
virtual void parseMultigridSettings()
Parse solver settings for geometric multigrid.
Definition: CD_CdrMultigrid.cpp:677
EBAMRCellData m_residual
Storage for multigrid residual.
Definition: CD_CdrMultigrid.H:230
int m_multigridBottomSmooth
Number of smoothing before bottom solver.
Definition: CD_CdrMultigrid.H:250
int m_minCellsBottom
Set bottom drop depth.
Definition: CD_CdrMultigrid.H:275
EBAMRCellData m_helmAcoef
Storage for Helmholtz a-coefficient. Always 1.
Definition: CD_CdrMultigrid.H:225
virtual void preRegrid(const int a_lbase, const int a_oldFinestLevel) override
Perform pre-regrid operations.
Definition: CD_CdrMultigrid.cpp:64
Real m_multigridExitHang
Multigrid hand.
Definition: CD_CdrMultigrid.H:285
BiCGStabSolver< LevelData< EBCellFAB > > m_bicgstab
Conjugate gradient solver bottom MG level.
Definition: CD_CdrMultigrid.H:205
virtual void parseRuntimeOptions() override=0
Parse class options.
virtual void resetAlphaAndBeta(const Real a_alpha, const Real a_beta)
Reset alpha and beta-coefficients in the multigrid solvers.
Definition: CD_CdrMultigrid.cpp:77
virtual void setupMultigrid()
Setup multigrid.
Definition: CD_CdrMultigrid.cpp:433
BottomSolverType m_bottomSolverType
Bottom solver type.
Definition: CD_CdrMultigrid.H:265
virtual void allocate() override
Allocate internal storage.
Definition: CD_CdrMultigrid.cpp:53
int m_multigridPreSmooth
Number of smoothings before averaging.
Definition: CD_CdrMultigrid.H:240
GMRESSolver< LevelData< EBCellFAB > > m_gmres
GMRES solver.
Definition: CD_CdrMultigrid.H:215
virtual void setupDiffusionSolver()
Set up diffusion solver.
Definition: CD_CdrMultigrid.cpp:343
virtual void setupHelmholtzFactory()
Setup the operator factory.
Definition: CD_CdrMultigrid.cpp:368
virtual void computeDivF(EBAMRCellData &a_divF, EBAMRCellData &a_phi, const Real a_extrapDt, const bool a_conservativeOnly, const bool a_ebFlux, const bool a_domainFlux) override final
Compute div(v*phi) explicitly.
Definition: CD_CdrMultigrid.cpp:590
int m_multigridPostSmooth
Number of smoothings before averaging.
Definition: CD_CdrMultigrid.H:245
EBHelmholtzOp::Smoother m_smoother
Relaxation type for gmg.
Definition: CD_CdrMultigrid.H:185
virtual void computeDivJ(EBAMRCellData &a_divJ, EBAMRCellData &a_phi, const Real a_extrapDt, const bool a_conservativeOnly, const bool a_ebFlux, const bool a_domainFlux) override final
Compute div(J) explicitly, where J = nV - D*grad(n)
Definition: CD_CdrMultigrid.cpp:525
virtual void advectToFaces(EBAMRFluxData &a_facePhi, const EBAMRCellData &a_phi, const Real a_extrapDt) override=0
Advection-only extrapolation to faces.
RefCountedPtr< EBHelmholtzOpFactory > m_helmholtzOpFactory
Operator factory.
Definition: CD_CdrMultigrid.H:200
virtual void parseOptions() override=0
Parse class options.
Real m_multigridExitTolerance
Multigrid tolerance.
Definition: CD_CdrMultigrid.H:280
int m_multigridMinIterations
Minimum number of iterations.
Definition: CD_CdrMultigrid.H:260
int m_multigridMaxIterations
Maximum number of iterations.
Definition: CD_CdrMultigrid.H:255
virtual ~CdrMultigrid()
Destructor.
Definition: CD_CdrMultigrid.cpp:34
virtual void advanceCrankNicholson(EBAMRCellData &a_newPhi, const EBAMRCellData &a_oldPhi, const EBAMRCellData &a_source, const Real a_dt) override final
Implicit diffusion Crank-Nicholson advance with source term.
Definition: CD_CdrMultigrid.cpp:244
Base class for solving convection-diffusion-reaction equations.
Definition: CD_CdrSolver.H:34
Smoother
Relaxation method for the operators.
Definition: CD_EBHelmholtzOp.H:47