chombo-discharge
CD_PhaseRealm.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_PhaseRealm_H
13 #define CD_PhaseRealm_H
14 
15 // Chombo includes
16 #include <DisjointBoxLayout.H>
17 #include <ProblemDomain.H>
18 
19 // Our includes
22 #include <CD_IrregAmrStencil.H>
23 #include <CD_LoadBalancing.H>
24 #include <CD_MFLevelGrid.H>
25 #include <CD_EBCoarAve.H>
26 #include <CD_EBReflux.H>
27 #include <CD_DomainFluxIFFAB.H>
29 #include <CD_ParticleContainer.H>
31 #include <CD_EBGradient.H>
33 #include <CD_EBAMRParticleMesh.H>
40 #include <CD_NamespaceHeader.H>
41 
42 // These are operator that can be defined.
43 static const std::string s_eb_gradient = "eb_gradient"; // For computing gradients
44 static const std::string s_eb_irreg_interp = "eb_irreg_interp"; // For data recentering
45 static const std::string s_eb_coar_ave = "eb_coar_ave"; // For coarsening
46 static const std::string s_eb_fill_patch = "eb_fill_patch"; // For regridding data
47 static const std::string s_eb_fine_interp = "eb_fine_interp"; // For linearly filling ghost cells
48 static const std::string s_eb_flux_reg = "eb_flux_reg"; // For flux registeration
49 static const std::string s_eb_redist = "eb_redist"; // For redistribution
50 static const std::string s_noncons_div = "eb_non_cons_div"; // For computing non-conservative divergences
51 static const std::string s_eb_multigrid = "eb_multigrid"; // For multigrid interpolation
52 static const std::string s_levelset = "levelset"; // For putting level-set on mesh
53 static const std::string s_particle_mesh = "particle_mesh"; // For doing particle-mesh operations.
54 
66 {
67 public:
71  PhaseRealm();
72 
77  PhaseRealm(const PhaseRealm& a_other) = delete;
78 
83  PhaseRealm(const PhaseRealm&& a_other) = delete;
84 
89  PhaseRealm&
90  operator=(const PhaseRealm& a_other) = delete;
91 
96  PhaseRealm&
97  operator=(const PhaseRealm&& a_other) = delete;
98 
102  virtual ~PhaseRealm();
103 
124  void
125  define(const Vector<DisjointBoxLayout>& a_grids,
126  const Vector<ProblemDomain>& a_domains,
127  const Vector<int>& a_refRat,
128  const Vector<Real>& a_dx,
129  const RealVect a_probLo,
130  const int a_finestLevel,
131  const int a_ebGhost,
132  const int a_numGhost,
133  const int a_lsfGhost,
134  const int a_redistRad,
135  const int a_mgInterpOrder,
136  const int a_mgInterpRadius,
137  const int a_mgInterpWeight,
138  const IrregStencil::StencilType a_centroidStencil,
139  const IrregStencil::StencilType a_ebStencil,
140  const RefCountedPtr<BaseIF>& a_baseif,
141  const RefCountedPtr<EBIndexSpace>& a_ebis);
142 
149  void
150  setGrids(const Vector<DisjointBoxLayout>& a_grids, const int a_finestLevel);
151 
155  void
156  preRegrid();
157 
163  void
164  regridBase(const int a_lmin);
165 
171  void
172  regridOperators(const int a_lmin);
173 
179  void
180  registerOperator(const std::string a_operator);
181 
187  bool
188  queryOperator(const std::string a_operator) const;
189 
193  const RefCountedPtr<EBIndexSpace>&
194  getEBIndexSpace() const;
195 
200  const Vector<int>&
201  getRefinementRatios() const;
202 
207  const Vector<Real>&
208  getDx() const;
209 
214  const Vector<DisjointBoxLayout>&
215  getGrids() const;
216 
221  const Vector<ProblemDomain>&
222  getDomains() const;
223 
228  const Vector<EBISLayout>&
229  getEBISLayout() const;
230 
235  const Vector<RefCountedPtr<EBLevelGrid>>&
236  getEBLevelGrid() const;
237 
242  const Vector<RefCountedPtr<EBLevelGrid>>&
243  getEBLevelGridCoFi() const;
244 
249  Vector<RefCountedPtr<LayoutData<VoFIterator>>>&
250  getVofIterator() const;
251 
257 
263 
269 
273  const Vector<RefCountedPtr<EBGradient>>&
274  getGradientOp() const;
275 
280  getParticleMesh() const;
281 
286  getSurfaceDeposition() const;
287 
291  Vector<RefCountedPtr<EBCoarAve>>&
292  getCoarseAverage() const;
293 
297  Vector<RefCountedPtr<EBCoarseFineParticleMesh>>&
299 
303  Vector<RefCountedPtr<EBMultigridInterpolator>>&
304  getMultigridInterpolator() const;
305 
309  Vector<RefCountedPtr<EBGhostCellInterpolator>>&
310  getGhostCellInterpolator() const;
311 
315  Vector<RefCountedPtr<EBCoarseToFineInterp>>&
316  getFineInterp() const;
317 
321  Vector<RefCountedPtr<EBReflux>>&
322  getFluxRegister() const;
323 
328  Vector<RefCountedPtr<EBFluxRedistribution>>&
329  getRedistributionOp() const;
330 
334  const EBAMRFAB&
335  getLevelset() const;
336 
337 protected:
342 
346  bool m_hasEbCf;
347 
351  bool m_profile;
352 
356  bool m_verbose;
357 
362 
367 
372 
377 
382 
387 
392 
397 
401  RealVect m_probLo;
402 
407 
412 
416  Vector<Real> m_dx;
417 
421  Vector<int> m_refinementRatios;
422 
426  RefCountedPtr<EBIndexSpace> m_ebis;
427 
431  RefCountedPtr<BaseIF> m_baseif;
432 
436  std::map<std::string, bool> m_operatorMap;
437 
441  Vector<DisjointBoxLayout> m_grids;
442 
446  Vector<ProblemDomain> m_domains;
447 
451  Vector<EBISLayout> m_ebisl;
452 
456  Vector<RefCountedPtr<EBLevelGrid>> m_eblg;
457 
462  Vector<RefCountedPtr<EBLevelGrid>> m_eblgCoFi;
463 
468  Vector<RefCountedPtr<EBLevelGrid>> m_eblgFiCo;
469 
474 
478  mutable Vector<RefCountedPtr<LayoutData<VoFIterator>>> m_vofIter;
479 
483  mutable Vector<RefCountedPtr<EBCoarAve>> m_coarAve;
484 
488  mutable Vector<RefCountedPtr<EBMultigridInterpolator>> m_multigridInterpolator;
489 
493  mutable Vector<RefCountedPtr<EBGhostCellInterpolator>> m_ghostCellInterpolator;
494 
498  mutable Vector<RefCountedPtr<EBCoarseToFineInterp>> m_ebFineInterp;
499 
503  mutable Vector<RefCountedPtr<EBReflux>> m_ebReflux;
504 
508  mutable Vector<RefCountedPtr<EBFluxRedistribution>> m_redistributionOp;
509 
513  mutable Vector<RefCountedPtr<EBGradient>> m_gradientOp;
514 
519 
524 
528  RefCountedPtr<IrregAmrStencil<CentroidInterpolationStencil>> m_centroidInterpolationStencil;
529 
533  RefCountedPtr<IrregAmrStencil<EbCentroidInterpolationStencil>> m_ebCentroidInterpolationStencil;
534 
538  RefCountedPtr<IrregAmrStencil<NonConservativeDivergenceStencil>> m_NonConservativeDivergenceStencil;
539 
544  void
545  defineEBLevelGrid(const int a_lmin);
546 
551  void
552  defineVofIterator(const int a_lmin);
553 
558  void
559  defineEBCoarAve(const int a_lmin);
560 
565  void
566  defineEBMultigrid(const int a_lmin);
567 
572  void
573  defineFillPatch(const int a_lmin);
574 
579  void
580  defineEBCoarseToFineInterp(const int a_lmin);
581 
586  void
587  defineMultigridInjection(const int a_lmin);
588 
593  void
594  defineFluxReg(const int a_lmin, const int a_regsize);
595 
600  void
601  defineRedistOper(const int a_lmin, const int a_regsize);
602 
607  void
608  defineFineToCoarRedistOper(const int a_lmin, const int a_regsize);
609 
614  void
615  defineCoarToFineRedistOper(const int a_lmin, const int a_regsize);
616 
621  void
622  defineCoarToCoarRedistOper(const int a_lmin, const int a_regsize);
623 
627  void
629 
634  void
635  defineGradSten(const int a_lmin);
636 
641  void
642  defineIrregSten();
643 
648  void
650 
655  void
656  defineLevelSet(const int a_lmin, const int a_numGhost);
657 };
658 
659 #include <CD_NamespaceFooter.H>
660 
661 #endif
Implementation of IrregStencil that interpolates from cell centers to cell centroids.
Declaration of base class for defining geometries.
Declaration of factory class for DomainFluxIFFAB.
Declaration of a BaseIFFAB wrapper that holds domain fluxes.
Declaration of a class for handling particle-mesh operations with AMR.
Declaration of a class for handling surface deposition of particles with EB and AMR.
Declaration of conservative coarsening utility.
Declaration of an aggregated class for regrid operations.
Declaration of a class which can perform redistribution in an EB-AMR context.
Declaration of a class which interpolates ghost cells across the coarse-fine interface.
Declaration of a class which computes gradients in an EBAMR context.
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 IrregStencil-derived class that can from cell centers to embedded boundary centroids.
Class for holding stencils on irregular over an entire AMR hierarchy.
Declaration of a static class for various load balancing operations.
Declaration of a wrapper for wrapping multifluid EBLevelGrids.
Multi-fluid index space.
Implementation of IrregStencil that can perform the nonconservative divergence averaging.
Declaration of a class for holding particles on an AMR hierarchy.
Class for handling particle-mesh operations with AMR.
Definition: CD_EBAMRParticleMesh.H:47
class for handling surface deposition of particles with EB and AMR.
Definition: CD_EBAMRSurfaceDeposition.H:29
Class for holding VoFStencils on irregular cells over an entire AMR hierarchy. The template parameter...
Definition: CD_IrregAmrStencil.H:29
StencilType
Enum for identifying stencil – only meant for enhancing code visibility.
Definition: CD_IrregStencil.H:44
Class that holds important things for doing AMR over a specific phase and processor distribution.
Definition: CD_PhaseRealm.H:66
const IrregAmrStencil< NonConservativeDivergenceStencil > & getNonConservativeDivergenceStencils() const
Get stencils for computing "non-conservative divergences".
Definition: CD_PhaseRealm.cpp:939
Vector< RefCountedPtr< EBLevelGrid > > m_eblgFiCo
Coarsened fine-level EB grids.
Definition: CD_PhaseRealm.H:468
Vector< RefCountedPtr< EBCoarseToFineInterp > > & getFineInterp() const
Get piecewise linear ghost cell interpolation utility.
Definition: CD_PhaseRealm.cpp:979
RefCountedPtr< IrregAmrStencil< NonConservativeDivergenceStencil > > m_NonConservativeDivergenceStencil
Stencils for non-conservative divergences.
Definition: CD_PhaseRealm.H:538
virtual ~PhaseRealm()
Destructor.
Definition: CD_PhaseRealm.cpp:46
int m_finestLevel
Finest grid level.
Definition: CD_PhaseRealm.H:361
EBAMRSurfaceDeposition m_surfaceDeposition
For doing particle deposition onto surfaces.
Definition: CD_PhaseRealm.H:523
Vector< RefCountedPtr< EBMultigridInterpolator > > m_multigridInterpolator
Multigrid interpolation utility.
Definition: CD_PhaseRealm.H:488
const Vector< EBISLayout > & getEBISLayout() const
Get EBIS layout.
Definition: CD_PhaseRealm.cpp:893
std::map< std::string, bool > m_operatorMap
Operator map for checking which ones are registered.
Definition: CD_PhaseRealm.H:436
Vector< RefCountedPtr< EBGhostCellInterpolator > > m_ghostCellInterpolator
Ghost cell interpolation utility.
Definition: CD_PhaseRealm.H:493
Vector< RefCountedPtr< EBCoarAve > > m_coarAve
Coarsening operator.
Definition: CD_PhaseRealm.H:483
Vector< RefCountedPtr< EBCoarseToFineInterp > > m_ebFineInterp
Regridding utility (for filling new grid patches)
Definition: CD_PhaseRealm.H:498
Vector< RefCountedPtr< EBLevelGrid > > m_eblgCoFi
Coarsened fine-level EB grids.
Definition: CD_PhaseRealm.H:462
Vector< DisjointBoxLayout > m_grids
AMR grids.
Definition: CD_PhaseRealm.H:441
void defineIrregSten()
Define data recentering stencils.
Definition: CD_PhaseRealm.cpp:797
PhaseRealm()
Default constructor. Must subsequently call define.
Definition: CD_PhaseRealm.cpp:28
PhaseRealm & operator=(const PhaseRealm &a_other)=delete
Disallowed copy assignment operator.
const Vector< DisjointBoxLayout > & getGrids() const
Get AMR grids.
Definition: CD_PhaseRealm.cpp:881
Vector< RefCountedPtr< EBFluxRedistribution > > & getRedistributionOp() const
Get the redistribution operators.
Definition: CD_PhaseRealm.cpp:999
void defineEBCoarAve(const int a_lmin)
Define coarsening utility.
Definition: CD_PhaseRealm.cpp:515
void defineCoarToFineRedistOper(const int a_lmin, const int a_regsize)
Define operator redistribution utility.
void registerOperator(const std::string a_operator)
Register an AMR operator.
Definition: CD_PhaseRealm.cpp:334
EBAMRFAB m_levelset
Level-set function.
Definition: CD_PhaseRealm.H:473
Vector< int > m_refinementRatios
Refinement ratios between levels.
Definition: CD_PhaseRealm.H:421
void defineRedistOper(const int a_lmin, const int a_regsize)
Define operator redistribution utility.
Definition: CD_PhaseRealm.cpp:673
Vector< RefCountedPtr< EBCoarAve > > & getCoarseAverage() const
Get coarsening operators.
Definition: CD_PhaseRealm.cpp:949
const Vector< RefCountedPtr< EBLevelGrid > > & getEBLevelGridCoFi() const
Get coarsened fine grids.
Definition: CD_PhaseRealm.cpp:905
void defineEBMultigrid(const int a_lmin)
Define multigrid interpolation utility.
Definition: CD_PhaseRealm.cpp:540
PhaseRealm(const PhaseRealm &&a_other)=delete
Disallowed move ctor.
Vector< RefCountedPtr< EBReflux > > & getFluxRegister() const
Get flux register utility.
Definition: CD_PhaseRealm.cpp:989
Vector< RefCountedPtr< EBFluxRedistribution > > m_redistributionOp
Redistribution utilities.
Definition: CD_PhaseRealm.H:508
Vector< RefCountedPtr< EBCoarseFineParticleMesh > > & getEBCoarseFineParticleMesh() const
Get deposition operator.
Vector< RefCountedPtr< EBLevelGrid > > m_eblg
EB level grids.
Definition: CD_PhaseRealm.H:456
EBAMRParticleMesh & getParticleMesh() const
Get particle-mesh operator.
Definition: CD_PhaseRealm.cpp:1009
Vector< RefCountedPtr< LayoutData< VoFIterator > > > & getVofIterator() const
Return vof iterator for iterating over irregular cells in each grid patch.
Definition: CD_PhaseRealm.cpp:911
Vector< EBISLayout > m_ebisl
EBIS layouts.
Definition: CD_PhaseRealm.H:451
void defineCoarToCoarRedistOper(const int a_lmin, const int a_regsize)
Define operator redistribution utility.
void defineEBLevelGrid(const int a_lmin)
Define EBLevelGrids.
Definition: CD_PhaseRealm.cpp:380
int m_multigridInterpolationWeight
Multigrid interpolator weight (for least squares)
Definition: CD_PhaseRealm.H:396
Vector< Real > m_dx
Grid resolutions.
Definition: CD_PhaseRealm.H:416
void defineLevelSet(const int a_lmin, const int a_numGhost)
Put the level-set on the mesh.
Definition: CD_PhaseRealm.cpp:462
RefCountedPtr< BaseIF > m_baseif
Implicit/SD function.
Definition: CD_PhaseRealm.H:431
void defineNonConsDivSten()
Define non-conservative divergence stencils.
Definition: CD_PhaseRealm.cpp:837
void setGrids(const Vector< DisjointBoxLayout > &a_grids, const int a_finestLevel)
Set grid method.
Definition: CD_PhaseRealm.cpp:101
int m_multigridInterpolationOrder
Multigrid interpolator order.
Definition: CD_PhaseRealm.H:386
void defineFluxReg(const int a_lmin, const int a_regsize)
Define flux register utility.
Definition: CD_PhaseRealm.cpp:640
void defineMultigridInjection(const int a_lmin)
Define multigrid injection utility.
EBAMRSurfaceDeposition & getSurfaceDeposition() const
Get the surface deposition operator.
Definition: CD_PhaseRealm.cpp:1019
bool m_verbose
Verbose or not.
Definition: CD_PhaseRealm.H:356
Vector< RefCountedPtr< LayoutData< VoFIterator > > > m_vofIter
Cut-cell iterator.
Definition: CD_PhaseRealm.H:478
Vector< RefCountedPtr< EBMultigridInterpolator > > & getMultigridInterpolator() const
Get multigrid interpolator.
Definition: CD_PhaseRealm.cpp:959
EBAMRParticleMesh m_particleMesh
For doing particle-mesh operations.
Definition: CD_PhaseRealm.H:518
PhaseRealm(const PhaseRealm &a_other)=delete
Disallowed copy ctor.
const Vector< Real > & getDx() const
Returns resolutions.
Definition: CD_PhaseRealm.cpp:875
RealVect m_probLo
Lower-left corner of computational domain.
Definition: CD_PhaseRealm.H:401
const EBAMRFAB & getLevelset() const
Get level-set function.
Definition: CD_PhaseRealm.cpp:1029
int m_numGhostCells
Number of ghost cells in data holders.
Definition: CD_PhaseRealm.H:366
void regridBase(const int a_lmin)
Regrid method for EBAMR base.
Definition: CD_PhaseRealm.cpp:142
Vector< RefCountedPtr< EBGhostCellInterpolator > > & getGhostCellInterpolator() const
Get ghost cell interpolation utility.
Definition: CD_PhaseRealm.cpp:969
void defineEBCoarseToFineInterp(const int a_lmin)
Define ghost cell interpolation utility.
Definition: CD_PhaseRealm.cpp:614
void defineGradSten(const int a_lmin)
Define gradient stencils.
Definition: CD_PhaseRealm.cpp:744
PhaseRealm & operator=(const PhaseRealm &&a_other)=delete
Disallowed move assignment operator.
Vector< RefCountedPtr< EBReflux > > m_ebReflux
Flux register utility.
Definition: CD_PhaseRealm.H:503
const Vector< RefCountedPtr< EBLevelGrid > > & getEBLevelGrid() const
Get EBLevelGrids.
Definition: CD_PhaseRealm.cpp:899
void define(const Vector< DisjointBoxLayout > &a_grids, const Vector< ProblemDomain > &a_domains, const Vector< int > &a_refRat, const Vector< Real > &a_dx, const RealVect a_probLo, const int a_finestLevel, const int a_ebGhost, const int a_numGhost, const int a_lsfGhost, const int a_redistRad, const int a_mgInterpOrder, const int a_mgInterpRadius, const int a_mgInterpWeight, const IrregStencil::StencilType a_centroidStencil, const IrregStencil::StencilType a_ebStencil, const RefCountedPtr< BaseIF > &a_baseif, const RefCountedPtr< EBIndexSpace > &a_ebis)
Full define function.
Definition: CD_PhaseRealm.cpp:50
IrregStencil::StencilType m_ebCentroidStencilType
Stencil recentering type.
Definition: CD_PhaseRealm.H:411
const Vector< RefCountedPtr< EBGradient > > & getGradientOp() const
Get gradient object.
Definition: CD_PhaseRealm.cpp:929
const Vector< int > & getRefinementRatios() const
Get refinment ratios.
Definition: CD_PhaseRealm.cpp:869
void defineFillPatch(const int a_lmin)
Define regrid utility.
Definition: CD_PhaseRealm.cpp:579
bool m_hasEbCf
Relic of an ancient past. Always true.
Definition: CD_PhaseRealm.H:346
void defineParticleMesh()
Define particle mesh operator.
Definition: CD_PhaseRealm.cpp:732
bool queryOperator(const std::string a_operator) const
Query if an AMR operator has been registered.
Definition: CD_PhaseRealm.cpp:359
IrregStencil::StencilType m_centroidStencilType
Stencil recentering type.
Definition: CD_PhaseRealm.H:406
int m_multigridInterpolationRadius
Multigrid interpolator radius.
Definition: CD_PhaseRealm.H:391
const RefCountedPtr< EBIndexSpace > & getEBIndexSpace() const
Return the EB index space.
Definition: CD_PhaseRealm.cpp:863
void defineFineToCoarRedistOper(const int a_lmin, const int a_regsize)
Define operator redistribution utility.
bool m_isDefined
True if things on this phase can be defined. False otherwise. Only used internally.
Definition: CD_PhaseRealm.H:341
Vector< RefCountedPtr< EBGradient > > m_gradientOp
Gradient operator.
Definition: CD_PhaseRealm.H:513
void preRegrid()
Perform pre-regrid operations.
Definition: CD_PhaseRealm.cpp:115
RefCountedPtr< IrregAmrStencil< EbCentroidInterpolationStencil > > m_ebCentroidInterpolationStencil
Data recentering utility (center to eb centroid)
Definition: CD_PhaseRealm.H:533
RefCountedPtr< IrregAmrStencil< CentroidInterpolationStencil > > m_centroidInterpolationStencil
Data recentering utility (center to centroid)
Definition: CD_PhaseRealm.H:528
RefCountedPtr< EBIndexSpace > m_ebis
EB index space.
Definition: CD_PhaseRealm.H:426
const IrregAmrStencil< EbCentroidInterpolationStencil > & getEbCentroidInterpolationStencilStencils() const
Get stencils for interpolating from cell center to eb centroids.
Definition: CD_PhaseRealm.cpp:923
int m_numLsfGhostCells
Number of ghost cells in level-set function (on the mesh)
Definition: CD_PhaseRealm.H:376
bool m_profile
Profile or not.
Definition: CD_PhaseRealm.H:351
void defineVofIterator(const int a_lmin)
Define vof iterators.
Definition: CD_PhaseRealm.cpp:428
int m_redistributionRadius
Redistribution radius.
Definition: CD_PhaseRealm.H:381
Vector< ProblemDomain > m_domains
Problem domains.
Definition: CD_PhaseRealm.H:446
const Vector< ProblemDomain > & getDomains() const
Get problem domains.
Definition: CD_PhaseRealm.cpp:887
int m_numEbGhostsCells
Number of ghost cells in EBLevelGrid/EBISLayout/EBISBox.
Definition: CD_PhaseRealm.H:371
const IrregAmrStencil< CentroidInterpolationStencil > & getCentroidInterpolationStencils() const
Get stencils for interpolating from cell center to centroids.
Definition: CD_PhaseRealm.cpp:917
void regridOperators(const int a_lmin)
Regrid method for EBAMR operators.
Definition: CD_PhaseRealm.cpp:184