chombo-discharge
CD_MFHelmholtzOpFactory.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_MFHelmholtzOpFactory_H
13 #define CD_MFHelmholtzOpFactory_H
14 
15 // Chombo includes
16 #include <MFCellFAB.H>
17 #include <MFFluxFAB.H>
18 
19 // Our includes
21 #include <CD_Location.H>
22 #include <CD_EBAMRData.H>
23 #include <CD_MFHelmholtzOp.H>
24 #include <CD_MFCoarAve.H>
25 #include <CD_MFReflux.H>
27 #include <CD_MFLevelGrid.H>
28 #include <CD_MFBaseIVFAB.H>
31 #include <CD_NamespaceHeader.H>
32 
39 class MFHelmholtzOpFactory : public AMRLevelOpFactory<LevelData<MFCellFAB>>
40 {
41 public:
42  // Various alias for cutting down on typing.
43  using MFIS = RefCountedPtr<MultiFluidIndexSpace>;
44  using AmrLevelGrids = Vector<MFLevelGrid>;
45  using AmrInterpolators = Vector<MFMultigridInterpolator>;
46  using AmrFluxRegisters = Vector<MFReflux>;
47  using AmrCoarseners = Vector<MFCoarAve>;
48  using AmrResolutions = Vector<Real>;
49  using AmrRefRatios = Vector<int>;
50 
51  using AmrCellData = Vector<RefCountedPtr<LevelData<MFCellFAB>>>;
52  using AmrFluxData = Vector<RefCountedPtr<LevelData<MFFluxFAB>>>;
53  using AmrIrreData = Vector<RefCountedPtr<LevelData<MFBaseIVFAB>>>;
54 
55  using DomainBCFactory = RefCountedPtr<MFHelmholtzDomainBCFactory>;
56  using EBBCFactory = RefCountedPtr<MFHelmholtzEBBCFactory>;
57  using JumpBCFactory = RefCountedPtr<MFHelmholtzJumpBCFactory>;
58 
60 
65 
69  MFHelmholtzOpFactory(const MFHelmholtzOpFactory& a_otherFactory) = delete;
70 
96  MFHelmholtzOpFactory(const MFIS& a_mfis,
97  const Location::Cell a_dataLocation,
98  const Real& a_alpha,
99  const Real& a_beta,
100  const RealVect& a_probLo,
101  const AmrLevelGrids& a_amrLevelGrids,
102  const AmrInterpolators& a_amrInterpolators,
103  const AmrFluxRegisters& a_amrFluxRegisters,
104  const AmrCoarseners& a_amrCoarseners,
105  const AmrRefRatios& a_amrRefRatios,
106  const AmrResolutions& a_amrResolutions,
107  const AmrCellData& a_amrAcoef,
108  const AmrFluxData& a_amrBcoef,
109  const AmrIrreData& a_amrBcoefIrreg,
110  const DomainBCFactory& a_domainBcFactory,
111  const EBBCFactory& a_ebBcFactory,
112  const JumpBCFactory& a_jumpBcFactory,
113  const IntVect& a_ghostPhi,
114  const IntVect& a_ghostRhs,
115  const Smoother& a_smoother,
116  const ProblemDomain& a_bottomDomain,
117  const int& a_jumpOrder,
118  const int& a_jumpWeight,
119  const int& a_blockingFactor,
120  const AmrLevelGrids& a_deeperLevelGrids = AmrLevelGrids());
121 
126 
132  void
133  setJump(const EBAMRIVData& a_sigma, const Real& a_scale);
134 
139  void
140  setJump(const Real& a_sigma, const Real& a_scale);
141 
145  const EBAMRIVData&
146  getSigma() const;
147 
151  void
153 
161  MGnewOp(const ProblemDomain& a_fineDomain, int a_depth, bool a_homogeneousOnly = true) override final;
162 
168  AMRnewOp(const ProblemDomain& a_domain) override final;
169 
174  int
175  refToFiner(const ProblemDomain& a_indexspace) const override final;
176 
177 protected:
181  static constexpr int m_comp = 0;
182 
186  static constexpr int m_nComp = 1;
187 
191  static constexpr int m_mainPhase = 0;
192 
196  MFIS m_mfis;
197 
202 
207 
212 
216  IntVect m_ghostPhi;
217 
221  IntVect m_ghostRhs;
222 
226  Real m_alpha;
227 
231  Real m_beta;
232 
236  RealVect m_probLo;
237 
241  AmrLevelGrids m_amrLevelGrids;
242 
246  AmrInterpolators m_amrInterpolators;
247 
251  AmrFluxRegisters m_amrFluxRegisters;
252 
256  AmrCoarseners m_amrCoarseners;
257 
261  AmrRefRatios m_amrRefRatios;
262 
266  AmrResolutions m_amrResolutions;
267 
271  AmrCellData m_amrAcoef;
272 
276  AmrFluxData m_amrBcoef;
277 
281  AmrIrreData m_amrBcoefIrreg;
282 
286  DomainBCFactory m_domainBcFactory;
287 
291  EBBCFactory m_ebBcFactory;
292 
296  JumpBCFactory m_jumpBcFactory;
297 
302 
306  ProblemDomain m_bottomDomain;
307 
312 
317 
322 
326  AmrLevelGrids m_deeperLevelGrids;
327 
331  std::vector<bool> m_hasMgLevels;
332 
336  Vector<AmrLevelGrids> m_mgLevelGrids;
337 
341  Vector<AmrCellData> m_mgAcoef;
342 
346  Vector<AmrFluxData> m_mgBcoef;
347 
351  Vector<AmrIrreData> m_mgBcoefIrreg;
352 
357 
361  Vector<Vector<RefCountedPtr<EBCoarAve>>> m_mgAveOp;
362 
366  void
368 
372  void
373  defineJump();
374 
381  bool
382  isCoarser(const ProblemDomain& a_domainOne, const ProblemDomain& a_domainTwo) const;
383 
390  bool
391  isFiner(const ProblemDomain& a_domainOne, const ProblemDomain& a_domainTwo) const;
392 
402  bool
403  getCoarserLayout(MFLevelGrid& a_coarseGrid,
404  const MFLevelGrid& a_fineGrid,
405  const int a_refRat,
406  const int a_blockingFactor) const;
407 
420  void
421  coarsenCoefficients(LevelData<MFCellFAB>& a_coarAcoef,
422  LevelData<MFFluxFAB>& a_coarBcoef,
423  LevelData<MFBaseIVFAB>& a_coarBcoefIrreg,
424  const LevelData<MFCellFAB>& a_fineAcoef,
425  const LevelData<MFFluxFAB>& a_fineBcoef,
426  const LevelData<MFBaseIVFAB>& a_fineBcoefIrreg,
427  const MFLevelGrid& a_eblgCoar,
428  const MFLevelGrid& a_eblgFine,
429  const int a_refRat);
430 
437  int
438  findAmrLevel(const ProblemDomain& a_domain) const;
439 };
440 
441 #include <CD_NamespaceFooter.H>
442 
443 #endif
Class for holding data across EBAMR hierarchies.
Declaration of cell positions.
Declaration of a multiphase BaseIVFAB<Real>
Wrapper class for holding multifluid EBCoarAves.
Declaration of an EB boundary condition factory class for MFHelmholtzOp.
Declaration of a factory class for making MFHelmholtzJumpBC objects for use in a multifluid MFHelmhol...
Declaration of a class for solving multiphase Helmholtz equations.
Declaration of a wrapper for wrapping multifluid EBLevelGrids.
Declaration of a wrapper class for holding multifluid EBMultigridInterpolators.
Declaration of a class for refluxing in a multiphase context.
Multi-fluid index space.
Class which replaces data at coarse level of refinement with average at fine level of refinement.
Definition: CD_EBCoarAve.H:31
Multiphase BaseIVFAB<Real>.
Definition: CD_MFBaseIVFAB.H:29
Factory class for making MFHelmholtzOp.
Definition: CD_MFHelmholtzOpFactory.H:40
static constexpr int m_comp
Component that this operator solves for.
Definition: CD_MFHelmholtzOpFactory.H:181
ProblemDomain m_bottomDomain
Coarsest level where we run multigrid.
Definition: CD_MFHelmholtzOpFactory.H:306
void coarsenCoefficients(LevelData< MFCellFAB > &a_coarAcoef, LevelData< MFFluxFAB > &a_coarBcoef, LevelData< MFBaseIVFAB > &a_coarBcoefIrreg, const LevelData< MFCellFAB > &a_fineAcoef, const LevelData< MFFluxFAB > &a_fineBcoef, const LevelData< MFBaseIVFAB > &a_fineBcoefIrreg, const MFLevelGrid &a_eblgCoar, const MFLevelGrid &a_eblgFine, const int a_refRat)
Coarsen coefficients (conservatively)
Definition: CD_MFHelmholtzOpFactory.cpp:426
Vector< EBAMRIVData > m_mgJump
Jump coefficient on deeper grids.
Definition: CD_MFHelmholtzOpFactory.H:356
EBBCFactory m_ebBcFactory
EB BC factory.
Definition: CD_MFHelmholtzOpFactory.H:291
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_MFHelmholtzOpFactory.H:346
std::vector< bool > m_hasMgLevels
For checking if an AMR level has multigrid levels.
Definition: CD_MFHelmholtzOpFactory.H:331
EBAMRIVData m_amrJump
Jump on multiphase cells.
Definition: CD_MFHelmholtzOpFactory.H:301
AmrFluxData m_amrBcoef
Operator B-coefficient.
Definition: CD_MFHelmholtzOpFactory.H:276
int m_mgBlockingFactor
Blocking factor for grid aggregration.
Definition: CD_MFHelmholtzOpFactory.H:311
AmrFluxRegisters m_amrFluxRegisters
Flux registers.
Definition: CD_MFHelmholtzOpFactory.H:251
bool isFiner(const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const
Check if a domain is finer than the other.
Definition: CD_MFHelmholtzOpFactory.cpp:795
Location::Cell m_dataLocation
Data centering.
Definition: CD_MFHelmholtzOpFactory.H:201
MFHelmholtzOp * AMRnewOp(const ProblemDomain &a_domain) override final
Create AMR operator for specified domain.
Definition: CD_MFHelmholtzOpFactory.cpp:705
MFHelmholtzOpFactory()=delete
Disallowed constructor.
int m_numAmrLevels
Number of AMR levels.
Definition: CD_MFHelmholtzOpFactory.H:211
Vector< AmrLevelGrids > m_mgLevelGrids
Deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to the the multigrid ...
Definition: CD_MFHelmholtzOpFactory.H:336
Real m_alpha
Operator alpha.
Definition: CD_MFHelmholtzOpFactory.H:226
int m_jumpWeight
Weighting for least squares reconstruction in jump cells.
Definition: CD_MFHelmholtzOpFactory.H:321
int refToFiner(const ProblemDomain &a_indexspace) const override final
Get refinement ratio to next finest level.
Definition: CD_MFHelmholtzOpFactory.cpp:803
AmrCoarseners m_amrCoarseners
Data coarseners.
Definition: CD_MFHelmholtzOpFactory.H:256
const EBAMRIVData & getSigma() const
Get the BC jump factor.
Definition: CD_MFHelmholtzOpFactory.cpp:184
MFHelmholtzOpFactory(const MFHelmholtzOpFactory &a_otherFactory)=delete
Disallowed 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_MFHelmholtzOpFactory.H:341
AmrIrreData m_amrBcoefIrreg
Operator B-coefficient.
Definition: CD_MFHelmholtzOpFactory.H:281
JumpBCFactory m_jumpBcFactory
Jump BC factory.
Definition: CD_MFHelmholtzOpFactory.H:296
AmrLevelGrids m_amrLevelGrids
AMR grids.
Definition: CD_MFHelmholtzOpFactory.H:241
bool getCoarserLayout(MFLevelGrid &a_coarseGrid, const MFLevelGrid &a_fineGrid, const int a_refRat, const int a_blockingFactor) const
Construct coarsening of a grid level.
Definition: CD_MFHelmholtzOpFactory.cpp:490
static constexpr int m_nComp
Number of components that this operator solves for.
Definition: CD_MFHelmholtzOpFactory.H:186
MFIS m_mfis
Index space.
Definition: CD_MFHelmholtzOpFactory.H:196
IntVect m_ghostPhi
Number of ghost cells that are used. Need because of Chombo prolongation objects.
Definition: CD_MFHelmholtzOpFactory.H:216
AmrResolutions m_amrResolutions
Grid resolutions.
Definition: CD_MFHelmholtzOpFactory.H:266
static constexpr int m_mainPhase
Main phase.
Definition: CD_MFHelmholtzOpFactory.H:191
void setJump(const EBAMRIVData &a_sigma, const Real &a_scale)
Set jump condition. This is passed to the operators by reference.
Definition: CD_MFHelmholtzOpFactory.cpp:112
DomainBCFactory m_domainBcFactory
Domain BC factory.
Definition: CD_MFHelmholtzOpFactory.H:286
int findAmrLevel(const ProblemDomain &a_domain) const
Find level corresponding to amr level.
Definition: CD_MFHelmholtzOpFactory.cpp:821
Smoother m_smoother
Smoother.
Definition: CD_MFHelmholtzOpFactory.H:206
RealVect m_probLo
Lower-left corner of computational domain.
Definition: CD_MFHelmholtzOpFactory.H:236
AmrRefRatios m_amrRefRatios
Refinement ratios.
Definition: CD_MFHelmholtzOpFactory.H:261
AmrLevelGrids m_deeperLevelGrids
Deeper grids.
Definition: CD_MFHelmholtzOpFactory.H:326
Vector< Vector< RefCountedPtr< EBCoarAve > > > m_mgAveOp
Coarsneing operator on deeper grids (used for coarsening jump)
Definition: CD_MFHelmholtzOpFactory.H:361
IntVect m_ghostRhs
Number of ghost cells that are used. Need because of Chombo prolongation objects.
Definition: CD_MFHelmholtzOpFactory.H:221
bool isCoarser(const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const
Check if a domain is coarser than the other.
Definition: CD_MFHelmholtzOpFactory.cpp:787
MFHelmholtzOp * MGnewOp(const ProblemDomain &a_fineDomain, int a_depth, bool a_homogeneousOnly=true) override final
Create multigrid operator.
Definition: CD_MFHelmholtzOpFactory.cpp:563
void coarsenCoefficientsMG()
Go through all MG levels and coarsen the coefficients from the finer levels.
Definition: CD_MFHelmholtzOpFactory.cpp:383
AmrCellData m_amrAcoef
Operator A-coefficient.
Definition: CD_MFHelmholtzOpFactory.H:271
AmrInterpolators m_amrInterpolators
Multigrid interpolators.
Definition: CD_MFHelmholtzOpFactory.H:246
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_MFHelmholtzOpFactory.H:351
Real m_beta
Operator beta.
Definition: CD_MFHelmholtzOpFactory.H:231
void defineMultigridLevels()
Function which defines the multigrid levels for this operator factory.
Definition: CD_MFHelmholtzOpFactory.cpp:229
void defineJump()
Define jump data.
Definition: CD_MFHelmholtzOpFactory.cpp:192
~MFHelmholtzOpFactory()
Destructor. Does nothing.
Definition: CD_MFHelmholtzOpFactory.cpp:106
int m_jumpOrder
Stencil order in jump cells.
Definition: CD_MFHelmholtzOpFactory.H:316
Operator for solving multifluid Helmholtz on a grid level.
Definition: CD_MFHelmholtzOp.H:43
Smoother
Relaxation methods for this operator.
Definition: CD_MFHelmholtzOp.H:49
Wrapper class for holding multifluid EBLevelGrids.
Definition: CD_MFLevelGrid.H:29
Namespace for encapsulating various data centerings.
Definition: CD_Location.H:24
Cell
Enum for distinguishing between cell locations.
Definition: CD_Location.H:30