chombo-discharge
CD_EBHelmholtzOpFactory.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_EBHelmholtzOpFactory_H
13 #define CD_EBHelmholtzOpFactory_H
14 
15 // Chombo includes
16 #include <RealVect.H>
17 #include <EBLevelGrid.H>
18 #include <EBCellFAB.H>
19 #include <EBFluxFAB.H>
20 #include <BaseIVFAB.H>
21 #include <ConductivityBaseDomainBC.H>
22 #include <BaseEBBC.H>
23 
24 // Our includes
25 #include <CD_Location.H>
26 #include <CD_EBCoarAve.H>
27 #include <CD_EBReflux.H>
29 #include <CD_EBHelmholtzOp.H>
30 #include <CD_EBHelmholtzEBBCFactory.H>
31 #include <CD_EBHelmholtzDomainBCFactory.H>
32 #include <CD_NamespaceHeader.H>
33 
40 class EBHelmholtzOpFactory : public AMRLevelOpFactory<LevelData<EBCellFAB>>
41 {
42 public:
43  // Various alias for cutting down on typing.
45 
46  using AmrLevelGrids = Vector<RefCountedPtr<EBLevelGrid>>;
47  using AmrInterpolators = Vector<RefCountedPtr<EBMultigridInterpolator>>;
48  using AmrFluxRegisters = Vector<RefCountedPtr<EBReflux>>;
49  using AmrCoarseners = Vector<RefCountedPtr<EBCoarAve>>;
50 
51  using AmrCellData = Vector<RefCountedPtr<LevelData<EBCellFAB>>>;
52  using AmrFluxData = Vector<RefCountedPtr<LevelData<EBFluxFAB>>>;
53  using AmrIrreData = Vector<RefCountedPtr<LevelData<BaseIVFAB<Real>>>>;
54 
55  using AmrResolutions = Vector<Real>;
56  using AmrRefRatios = Vector<int>;
57 
58  using DomainBCFactory = RefCountedPtr<EBHelmholtzDomainBCFactory>; // Will be replaced with better BC class...
59  using EBBCFactory = RefCountedPtr<EBHelmholtzEBBCFactory>;
60 
65 
69  EBHelmholtzOpFactory(const EBHelmholtzOpFactory& a_otherFactory) = delete;
70 
74  EBHelmholtzOpFactory(const EBHelmholtzOpFactory&& a_otherFactory) = delete;
75 
101  EBHelmholtzOpFactory(const Location::Cell a_dataLocation,
102  const Real& a_alpha,
103  const Real& a_beta,
104  const RealVect& a_probLo,
105  const AmrLevelGrids& a_amrLevelGrids,
106  const AmrInterpolators& a_amrInterpolators,
107  const AmrFluxRegisters& a_amrFluxRegisters,
108  const AmrCoarseners& a_amrCoarseners,
109  const AmrRefRatios& a_amrRefRatios,
110  const AmrResolutions& a_amrResolutions,
111  const AmrCellData& a_amrAcoef,
112  const AmrFluxData& a_amrBcoef,
113  const AmrIrreData& a_amrBcoefIrreg,
114  const DomainBCFactory& a_domainBcFactory,
115  const EBBCFactory& a_ebbcFactory,
116  const IntVect& a_ghostPhi,
117  const IntVect& a_ghostRHS,
118  const Smoother& a_relaxationMethod,
119  const ProblemDomain& a_bottomDomain,
120  const int& a_mgBlockingFactor,
121  const AmrLevelGrids& a_deeperLevelGrids = AmrLevelGrids());
122 
127 
131  void
132  operator=(const EBHelmholtzOpFactory& a_opin) = delete;
133 
137  void
138  operator=(const EBHelmholtzOpFactory&& a_opin) = delete;
139 
143  void
145 
153  MGnewOp(const ProblemDomain& a_fineDomain, int a_depth, bool a_homogeneousOnly = true) override final;
154 
160  AMRnewOp(const ProblemDomain& a_domain) override final;
161 
166  int
167  refToFiner(const ProblemDomain& a_indexspace) const override final;
168 
169 protected:
173  static constexpr int m_comp = 0;
174 
178  static constexpr int m_nComp = 1;
179 
184 
189 
194 
198  IntVect m_ghostPhi;
199 
203  IntVect m_ghostRhs;
204 
208  Real m_alpha;
209 
213  Real m_beta;
214 
218  RealVect m_probLo;
219 
220  // Things that pertain to AMR levels. The first entry corresponds to the coarsest AMR level.
224  AmrLevelGrids m_amrLevelGrids;
225 
229  AmrInterpolators m_amrInterpolators;
230 
234  AmrFluxRegisters m_amrFluxRegisters;
235 
239  AmrCoarseners m_amrCoarseners;
240 
244  AmrRefRatios m_amrRefRatios;
245 
249  AmrResolutions m_amrResolutions;
250 
254  AmrCellData m_amrAcoef;
255 
259  AmrFluxData m_amrBcoef;
260 
264  AmrIrreData m_amrBcoefIrreg;
265 
269  DomainBCFactory m_domainBcFactory;
270 
274  EBBCFactory m_ebBcFactory;
275 
279  ProblemDomain m_bottomDomain;
280 
285 
289  AmrLevelGrids m_deeperLevelGrids;
290 
294  std::vector<bool> m_hasMgLevels;
295 
299  Vector<AmrLevelGrids> m_mgLevelGrids;
300 
304  Vector<AmrCellData> m_mgAcoef;
305 
309  Vector<AmrFluxData> m_mgBcoef;
310 
314  Vector<AmrIrreData> m_mgBcoefIrreg;
315 
319  void
321 
328  bool
329  isCoarser(const ProblemDomain& a_domainOne, const ProblemDomain& a_domainTwo) const;
330 
337  bool
338  isFiner(const ProblemDomain& a_domainOne, const ProblemDomain& a_domainTwo) const;
339 
348  bool
349  getCoarserLayout(EBLevelGrid& a_coarseGrid,
350  const EBLevelGrid& a_fineGrid,
351  const int a_refRat,
352  const int a_blockingFactor) const;
353 
366  void
367  coarsenCoefficients(LevelData<EBCellFAB>& a_coarAcoef,
368  LevelData<EBFluxFAB>& a_coarBcoef,
369  LevelData<BaseIVFAB<Real>>& a_coarBcoefIrreg,
370  const LevelData<EBCellFAB>& a_fineAcoef,
371  const LevelData<EBFluxFAB>& a_fineBcoef,
372  const LevelData<BaseIVFAB<Real>>& a_fineBcoefIrreg,
373  const EBLevelGrid& a_eblgCoar,
374  const EBLevelGrid& a_eblgFine,
375  const int a_refRat);
376 
383  int
384  findAmrLevel(const ProblemDomain& a_domain) const;
385 };
386 
387 #include <CD_NamespaceFooter.H>
388 
389 #endif
Declaration of conservative coarsening utility.
Declaration of Helmholtz multigrid operators.
Declaration of a class that can interpolate more ghost cells near the coarse-fine boundary near the E...
Declaration of a class which can reflux over the coarse-fine interface.
Declaration of cell positions.
Factory class for making variable-coefficient Helmholtz operators.
Definition: CD_EBHelmholtzOpFactory.H:41
Vector< AmrIrreData > m_mgBcoefIrreg
B-coefficient on deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to t...
Definition: CD_EBHelmholtzOpFactory.H:314
int m_mgBlockingFactor
Blocking factor for when we create intermediate and deep multigrid levels.
Definition: CD_EBHelmholtzOpFactory.H:284
AmrIrreData m_amrBcoefIrreg
Helmholtz B-coefficient (on EB faces.
Definition: CD_EBHelmholtzOpFactory.H:264
AmrCoarseners m_amrCoarseners
Data coarseners.
Definition: CD_EBHelmholtzOpFactory.H:239
AmrResolutions m_amrResolutions
Resolutions one each level.
Definition: CD_EBHelmholtzOpFactory.H:249
bool isFiner(const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const
Check if a domain is finer than the other.
Definition: CD_EBHelmholtzOpFactory.cpp:341
std::vector< bool > m_hasMgLevels
For checking if an AMR level has multigrid levels.
Definition: CD_EBHelmholtzOpFactory.H:294
EBHelmholtzOpFactory()=delete
Disallowed constructor. Use the full constructor.
Vector< AmrCellData > m_mgAcoef
A-coefficient on deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to t...
Definition: CD_EBHelmholtzOpFactory.H:304
int findAmrLevel(const ProblemDomain &a_domain) const
Find level corresponding to amr level.
Definition: CD_EBHelmholtzOpFactory.cpp:654
bool getCoarserLayout(EBLevelGrid &a_coarseGrid, const EBLevelGrid &a_fineGrid, const int a_refRat, const int a_blockingFactor) const
Construct coarsening of a grid level.
Definition: CD_EBHelmholtzOpFactory.cpp:349
IntVect m_ghostRhs
Number of ghost cells that are used. Need because of Chombo prolongation objects.
Definition: CD_EBHelmholtzOpFactory.H:203
EBBCFactory m_ebBcFactory
EB BC factory.
Definition: CD_EBHelmholtzOpFactory.H:274
EBHelmholtzOpFactory(const EBHelmholtzOpFactory &&a_otherFactory)=delete
Disallowed move constructor. Use the full constructor.
AmrFluxRegisters m_amrFluxRegisters
Flux registers.
Definition: CD_EBHelmholtzOpFactory.H:234
void operator=(const EBHelmholtzOpFactory &&a_opin)=delete
Disallowed move assignment.
static constexpr int m_comp
Component number that is solved for.
Definition: CD_EBHelmholtzOpFactory.H:173
int refToFiner(const ProblemDomain &a_indexspace) const override final
Get refinement ratio to next finest level.
Definition: CD_EBHelmholtzOpFactory.cpp:631
Vector< AmrLevelGrids > m_mgLevelGrids
Deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to the the multigrid ...
Definition: CD_EBHelmholtzOpFactory.H:299
~EBHelmholtzOpFactory()
Destructor. Does nothing.
Definition: CD_EBHelmholtzOpFactory.cpp:89
AmrLevelGrids m_deeperLevelGrids
This is for using pre-defined grids for the deeper multigrid levels, i.e. for the levels that are coa...
Definition: CD_EBHelmholtzOpFactory.H:289
void operator=(const EBHelmholtzOpFactory &a_opin)=delete
Disallowed assignment operator.
void coarsenCoefficients(LevelData< EBCellFAB > &a_coarAcoef, LevelData< EBFluxFAB > &a_coarBcoef, LevelData< BaseIVFAB< Real >> &a_coarBcoefIrreg, const LevelData< EBCellFAB > &a_fineAcoef, const LevelData< EBFluxFAB > &a_fineBcoef, const LevelData< BaseIVFAB< Real >> &a_fineBcoefIrreg, const EBLevelGrid &a_eblgCoar, const EBLevelGrid &a_eblgFine, const int a_refRat)
Coarsen coefficients (conservatively)
Definition: CD_EBHelmholtzOpFactory.cpp:292
ProblemDomain m_bottomDomain
Bottom domain, i.e. the coarsest domain which will be used in multigrid.
Definition: CD_EBHelmholtzOpFactory.H:279
Smoother m_smoother
Smoother.
Definition: CD_EBHelmholtzOpFactory.H:188
int m_numAmrLevels
Number of AMR levels.
Definition: CD_EBHelmholtzOpFactory.H:193
EBHelmholtzOpFactory(const EBHelmholtzOpFactory &a_otherFactory)=delete
Disallowed constructor. Use the full constructor.
Location::Cell m_dataLocation
Data location.
Definition: CD_EBHelmholtzOpFactory.H:183
Vector< AmrFluxData > m_mgBcoef
B-coefficient on deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to t...
Definition: CD_EBHelmholtzOpFactory.H:309
IntVect m_ghostPhi
Number of ghost cells that are used. Need because of Chombo prolongation objects.
Definition: CD_EBHelmholtzOpFactory.H:198
EBHelmholtzOp * AMRnewOp(const ProblemDomain &a_domain) override final
Create AMR operator for specified domain.
Definition: CD_EBHelmholtzOpFactory.cpp:551
static constexpr int m_nComp
Number of components that we solve for.
Definition: CD_EBHelmholtzOpFactory.H:178
AmrFluxData m_amrBcoef
Helmholtz B-coefficient.
Definition: CD_EBHelmholtzOpFactory.H:259
AmrLevelGrids m_amrLevelGrids
AMR grids.
Definition: CD_EBHelmholtzOpFactory.H:224
bool isCoarser(const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const
Check if a domain is coarser than the other.
Definition: CD_EBHelmholtzOpFactory.cpp:333
AmrCellData m_amrAcoef
Helmholtz A-coefficient.
Definition: CD_EBHelmholtzOpFactory.H:254
void coarsenCoefficientsMG()
Go through all MG levels and coarsen the coefficients from the finer levels.
Definition: CD_EBHelmholtzOpFactory.cpp:249
EBHelmholtzOp * MGnewOp(const ProblemDomain &a_fineDomain, int a_depth, bool a_homogeneousOnly=true) override final
Create multigrid operator.
Definition: CD_EBHelmholtzOpFactory.cpp:418
void defineMultigridLevels()
Function which defines the multigrid levels for this operator factory.
Definition: CD_EBHelmholtzOpFactory.cpp:95
Real m_beta
Operator beta.
Definition: CD_EBHelmholtzOpFactory.H:213
DomainBCFactory m_domainBcFactory
Domain BC factory.
Definition: CD_EBHelmholtzOpFactory.H:269
AmrRefRatios m_amrRefRatios
Refinement ratios.
Definition: CD_EBHelmholtzOpFactory.H:244
RealVect m_probLo
Lower-left corner of computational domain.
Definition: CD_EBHelmholtzOpFactory.H:218
Real m_alpha
Operator alpha.
Definition: CD_EBHelmholtzOpFactory.H:208
AmrInterpolators m_amrInterpolators
Ghost cell interpolations.
Definition: CD_EBHelmholtzOpFactory.H:229
Helmholtz operator for equations like alpha*a(x)*phi(x) + beta*div(b(x)*grad(phi(x))) = rho.
Definition: CD_EBHelmholtzOp.H:41
Smoother
Relaxation method for the operators.
Definition: CD_EBHelmholtzOp.H:47
Namespace for encapsulating various data centerings.
Definition: CD_Location.H:24
Cell
Enum for distinguishing between cell locations.
Definition: CD_Location.H:30