chombo-discharge
Loading...
Searching...
No Matches
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_LoadBalancing.H>
23#include <CD_MFLevelGrid.H>
24#include <CD_EBCoarAve.H>
25#include <CD_EBReflux.H>
26#include <CD_DomainFluxIFFAB.H>
30#include <CD_EBGradient.H>
39#include <CD_NamespaceHeader.H>
40
41// These are operators that can be defined.
42static const std::string s_outer_cf_region = "outer_cf_region"; // Outside coarse-fine region
43static const std::string s_inner_cf_region = "inner_cf_region"; // Inside coarse-fine region
44static const std::string s_outer_particle_halo = "outer_particle_halo"; // Outside particle halo region
45static const std::string s_inner_particle_halo = "inner_particle_halo"; // Inside particle halo region
46static const std::string s_cfivs = "cfivs"; // Coarse-fine interface
47static const std::string s_eb_gradient = "eb_gradient"; // For computing gradients
48static const std::string s_eb_irreg_interp = "eb_irreg_interp"; // For data recentering
49static const std::string s_eb_coar_ave = "eb_coar_ave"; // For coarsening
50static const std::string s_eb_fill_patch = "eb_fill_patch"; // For regridding data
51static const std::string s_eb_fine_interp = "eb_fine_interp"; // For linearly filling ghost cells
52static const std::string s_eb_flux_reg = "eb_flux_reg"; // For flux registeration
53static const std::string s_eb_redist = "eb_redist"; // For redistribution
54static const std::string s_noncons_div = "eb_non_cons_div"; // For computing non-conservative divergences
55static const std::string s_eb_multigrid = "eb_multigrid"; // For multigrid interpolation
56static const std::string s_levelset = "levelset"; // For putting level-set on mesh
57static const std::string s_particle_mesh = "particle_mesh"; // For doing particle-mesh operations.
58
70{
71public:
75 PhaseRealm();
76
81 PhaseRealm(const PhaseRealm& a_other) = delete;
82
87 PhaseRealm(const PhaseRealm&& a_other) = delete;
88
94 operator=(const PhaseRealm& a_other) = delete;
95
101 operator=(const PhaseRealm&& a_other) = delete;
102
106 virtual ~PhaseRealm();
107
128 void
131 const Vector<int>& a_refRat,
132 const Vector<Real>& a_dx,
133 const RealVect a_probLo,
134 const int a_finestLevel,
135 const int a_ebGhost,
136 const int a_numGhost,
137 const int a_lsfGhost,
138 const int a_redistRad,
139 const int a_mgInterpOrder,
140 const int a_mgInterpRadius,
141 const int a_mgInterpWeight,
146
153 void
155
159 void
160 preRegrid();
161
167 void
168 regridBase(const int a_lmin);
169
175 void
176 regridOperators(const int a_lmin);
177
183 void
185
191 bool
193
198 getEBIndexSpace() const;
199
204 const Vector<int>&
205 getRefinementRatios() const;
206
211 const Vector<Real>&
212 getDx() const;
213
219 getGrids() const;
220
226 getDomains() const;
227
232 const Vector<EBISLayout>&
233 getEBISLayout() const;
234
239 getEBLevelGrid() const;
240
245 getEBLevelGridCoFi() const;
246
251 getEBLevelGridFiCo() const;
252
258 getVofIterator() const;
259
265
271
277
282 getGradientOp() const;
283
288 getParticleMesh() const;
289
294 getSurfaceDeposition() const;
295
300 getCoarseAverage() const;
301
307
313
319
324 getFineInterp() const;
325
330 getFluxRegister() const;
331
337 getRedistributionOp() const;
338
342 const EBAMRFAB&
343 getLevelset() const;
344
345protected:
350
355
360
365
370
375
380
385
390
395
400
405
410
415
420
425
430
435
440
445
450
455
460
465
471
477
482
487
492
497
502
507
512
517
522
527
532
537
542
547
552 void
553 defineEBLevelGrid(const int a_lmin);
554
559 void
560 defineVofIterator(const int a_lmin);
561
566 void
567 defineEBCoarAve(const int a_lmin);
568
573 void
574 defineEBMultigrid(const int a_lmin);
575
580 void
581 defineFillPatch(const int a_lmin);
582
587 void
589
594 void
596
601 void
602 defineFluxReg(const int a_lmin, const int a_regsize);
603
608 void
609 defineRedistOper(const int a_lmin, const int a_regsize);
610
615 void
617
622 void
624
629 void
631
635 void
637
642 void
643 defineGradSten(const int a_lmin);
644
649 void
651
656 void
658
663 void
664 defineLevelSet(const int a_lmin, const int a_numGhost);
665};
666
667#include <CD_NamespaceFooter.H>
668
669#endif
Declaration of class for interpolating cell-centered data 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 class for interpolating cell-centered data to cell centroids.
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 for computing non-conservative divergences.
Declaration of a class which can reflux over the coarse-fine interface.
Declaration of a static class for various load balancing operations.
Declaration of a wrapper for wrapping multifluid EBLevelGrids.
Multi-fluid index space.
Declaration of a class for holding particles on an AMR hierarchy.
Type
Supported interpolation types.
Definition CD_CellCentroidInterpolation.H:43
Class for handling particle-mesh operations with AMR.
Definition CD_EBAMRParticleMesh.H:52
class for handling surface deposition of particles with EB and AMR.
Definition CD_EBAMRSurfaceDeposition.H:29
Type
Supported interpolation types.
Definition CD_EBCentroidInterpolation.H:43
Class that holds important things for doing AMR over a specific phase and processor distribution.
Definition CD_PhaseRealm.H:70
PhaseRealm & operator=(const PhaseRealm &&a_other)=delete
Disallowed move assignment operator.
Vector< RefCountedPtr< EBLevelGrid > > m_eblgFiCo
Coarsened fine-level EB grids.
Definition CD_PhaseRealm.H:476
Vector< RefCountedPtr< EBCoarseToFineInterp > > & getFineInterp() const
Get piecewise linear ghost cell interpolation utility.
Definition CD_PhaseRealm.cpp:970
void defineNonConservativeDivergence(const int a_lmin)
Define non-conservative divergence stencils.
Definition CD_PhaseRealm.cpp:820
virtual ~PhaseRealm()
Destructor.
Definition CD_PhaseRealm.cpp:46
int m_finestLevel
Finest grid level.
Definition CD_PhaseRealm.H:369
EBAMRSurfaceDeposition m_surfaceDeposition
For doing particle deposition onto surfaces.
Definition CD_PhaseRealm.H:531
Vector< RefCountedPtr< EBMultigridInterpolator > > m_multigridInterpolator
Multigrid interpolation utility.
Definition CD_PhaseRealm.H:496
const Vector< EBISLayout > & getEBISLayout() const
Get EBIS layout.
Definition CD_PhaseRealm.cpp:870
const Vector< RefCountedPtr< EBNonConservativeDivergence > > & getNonConservativeDivergence() const
Get objects for computing non-conservative divergences.
Definition CD_PhaseRealm.cpp:930
std::map< std::string, bool > m_operatorMap
Operator map for checking which ones are registered.
Definition CD_PhaseRealm.H:444
Vector< RefCountedPtr< EBCoarseFineParticleMesh > > & getEBCoarseFineParticleMesh() const
Get deposition operator.
Vector< RefCountedPtr< EBGhostCellInterpolator > > m_ghostCellInterpolator
Ghost cell interpolation utility.
Definition CD_PhaseRealm.H:501
Vector< RefCountedPtr< EBCoarAve > > m_coarAve
Coarsening operator.
Definition CD_PhaseRealm.H:491
Vector< RefCountedPtr< EBCoarseToFineInterp > > m_ebFineInterp
Regridding utility (for filling new grid patches)
Definition CD_PhaseRealm.H:506
Vector< RefCountedPtr< EBLevelGrid > > m_eblgCoFi
Coarsened fine-level EB grids.
Definition CD_PhaseRealm.H:470
Vector< DisjointBoxLayout > m_grids
AMR grids.
Definition CD_PhaseRealm.H:449
void defineIrregSten()
Define data recentering stencils.
Definition CD_PhaseRealm.cpp:796
PhaseRealm()
Default constructor. Must subsequently call define.
Definition CD_PhaseRealm.cpp:28
const Vector< DisjointBoxLayout > & getGrids() const
Get AMR grids.
Definition CD_PhaseRealm.cpp:858
Vector< RefCountedPtr< EBFluxRedistribution > > & getRedistributionOp() const
Get the redistribution operators.
Definition CD_PhaseRealm.cpp:990
void defineEBCoarAve(const int a_lmin)
Define coarsening utility.
Definition CD_PhaseRealm.cpp:514
EBCentroidInterpolation::Type m_ebCentroidInterpolationType
Cell-center to cell-centroid interpolation method.
Definition CD_PhaseRealm.H:419
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:333
EBAMRFAB m_levelset
Level-set function.
Definition CD_PhaseRealm.H:481
Vector< int > m_refinementRatios
Refinement ratios between levels.
Definition CD_PhaseRealm.H:429
const Vector< RefCountedPtr< CellCentroidInterpolation > > & getCellCentroidInterpolation() const
Get objects for interpolation from cell centers to centroids.
Definition CD_PhaseRealm.cpp:910
void defineRedistOper(const int a_lmin, const int a_regsize)
Define operator redistribution utility.
Definition CD_PhaseRealm.cpp:672
Vector< RefCountedPtr< EBCoarAve > > & getCoarseAverage() const
Get coarsening operators.
Definition CD_PhaseRealm.cpp:940
const Vector< RefCountedPtr< EBLevelGrid > > & getEBLevelGridCoFi() const
Get coarsened fine grids.
Definition CD_PhaseRealm.cpp:882
void defineEBMultigrid(const int a_lmin)
Define multigrid interpolation utility.
Definition CD_PhaseRealm.cpp:539
PhaseRealm(const PhaseRealm &&a_other)=delete
Disallowed move ctor.
Vector< RefCountedPtr< EBReflux > > & getFluxRegister() const
Get flux register utility.
Definition CD_PhaseRealm.cpp:980
Vector< RefCountedPtr< EBFluxRedistribution > > m_redistributionOp
Redistribution utilities.
Definition CD_PhaseRealm.H:516
Vector< RefCountedPtr< EBLevelGrid > > m_eblg
EB level grids.
Definition CD_PhaseRealm.H:464
EBAMRParticleMesh & getParticleMesh() const
Get particle-mesh operator.
Definition CD_PhaseRealm.cpp:1000
Vector< RefCountedPtr< CellCentroidInterpolation > > m_cellCentroidInterpolation
For interpolating data from cell centers to cell centroids.
Definition CD_PhaseRealm.H:536
Vector< RefCountedPtr< LayoutData< VoFIterator > > > & getVofIterator() const
Return vof iterator for iterating over irregular cells in each grid patch.
Definition CD_PhaseRealm.cpp:894
Vector< EBISLayout > m_ebisl
EBIS layouts.
Definition CD_PhaseRealm.H:459
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:379
int m_multigridInterpolationWeight
Multigrid interpolator weight (for least squares)
Definition CD_PhaseRealm.H:404
Vector< Real > m_dx
Grid resolutions.
Definition CD_PhaseRealm.H:424
void defineLevelSet(const int a_lmin, const int a_numGhost)
Put the level-set on the mesh.
Definition CD_PhaseRealm.cpp:461
RefCountedPtr< BaseIF > m_baseif
Implicit/SD function.
Definition CD_PhaseRealm.H:439
void setGrids(const Vector< DisjointBoxLayout > &a_grids, const int a_finestLevel)
Set grid method.
Definition CD_PhaseRealm.cpp:101
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 CellCentroidInterpolation::Type a_centroidStencil, const EBCentroidInterpolation::Type a_ebStencil, const RefCountedPtr< BaseIF > &a_baseif, const RefCountedPtr< EBIndexSpace > &a_ebis)
Full define function.
Definition CD_PhaseRealm.cpp:50
int m_multigridInterpolationOrder
Multigrid interpolator order.
Definition CD_PhaseRealm.H:394
void defineFluxReg(const int a_lmin, const int a_regsize)
Define flux register utility.
Definition CD_PhaseRealm.cpp:639
void defineMultigridInjection(const int a_lmin)
Define multigrid injection utility.
EBAMRSurfaceDeposition & getSurfaceDeposition() const
Get the surface deposition operator.
Definition CD_PhaseRealm.cpp:1010
bool m_verbose
Verbose or not.
Definition CD_PhaseRealm.H:364
Vector< RefCountedPtr< LayoutData< VoFIterator > > > m_vofIter
Cut-cell iterator.
Definition CD_PhaseRealm.H:486
Vector< RefCountedPtr< EBMultigridInterpolator > > & getMultigridInterpolator() const
Get multigrid interpolator.
Definition CD_PhaseRealm.cpp:950
EBAMRParticleMesh m_particleMesh
For doing particle-mesh operations.
Definition CD_PhaseRealm.H:526
PhaseRealm(const PhaseRealm &a_other)=delete
Disallowed copy ctor.
const Vector< Real > & getDx() const
Returns resolutions.
Definition CD_PhaseRealm.cpp:852
CellCentroidInterpolation::Type m_cellCentroidInterpolationType
Cell-center to cell-centroid interpolation method.
Definition CD_PhaseRealm.H:414
RealVect m_probLo
Lower-left corner of computational domain.
Definition CD_PhaseRealm.H:409
Vector< RefCountedPtr< EBCentroidInterpolation > > m_ebCentroidInterpolation
For interpolating data from cell centers to EB centroids.
Definition CD_PhaseRealm.H:541
const EBAMRFAB & getLevelset() const
Get level-set function.
Definition CD_PhaseRealm.cpp:1020
int m_numGhostCells
Number of ghost cells in data holders.
Definition CD_PhaseRealm.H:374
Vector< RefCountedPtr< EBNonConservativeDivergence > > m_nonConservativeDivergence
For computing non-conservative divergences.
Definition CD_PhaseRealm.H:546
void regridBase(const int a_lmin)
Regrid method for EBAMR base.
Definition CD_PhaseRealm.cpp:141
Vector< RefCountedPtr< EBGhostCellInterpolator > > & getGhostCellInterpolator() const
Get ghost cell interpolation utility.
Definition CD_PhaseRealm.cpp:960
void defineEBCoarseToFineInterp(const int a_lmin)
Define ghost cell interpolation utility.
Definition CD_PhaseRealm.cpp:613
const Vector< RefCountedPtr< EBCentroidInterpolation > > & getEBCentroidInterpolation() const
Get objects for interpolation from cell centers to EB centroids.
Definition CD_PhaseRealm.cpp:920
void defineGradSten(const int a_lmin)
Define gradient stencils.
Definition CD_PhaseRealm.cpp:743
PhaseRealm & operator=(const PhaseRealm &a_other)=delete
Disallowed copy assignment operator.
Vector< RefCountedPtr< EBReflux > > m_ebReflux
Flux register utility.
Definition CD_PhaseRealm.H:511
const Vector< RefCountedPtr< EBLevelGrid > > & getEBLevelGrid() const
Get EBLevelGrids.
Definition CD_PhaseRealm.cpp:876
const Vector< RefCountedPtr< EBGradient > > & getGradientOp() const
Get gradient object.
Definition CD_PhaseRealm.cpp:900
const Vector< RefCountedPtr< EBLevelGrid > > & getEBLevelGridFiCo() const
Get refined grids.
Definition CD_PhaseRealm.cpp:888
const Vector< int > & getRefinementRatios() const
Get refinment ratios.
Definition CD_PhaseRealm.cpp:846
void defineFillPatch(const int a_lmin)
Define regrid utility.
Definition CD_PhaseRealm.cpp:578
bool m_hasEbCf
Relic of an ancient past. Always true.
Definition CD_PhaseRealm.H:354
void defineParticleMesh()
Define particle mesh operator.
Definition CD_PhaseRealm.cpp:731
bool queryOperator(const std::string a_operator) const
Query if an AMR operator has been registered.
Definition CD_PhaseRealm.cpp:358
int m_multigridInterpolationRadius
Multigrid interpolator radius.
Definition CD_PhaseRealm.H:399
const RefCountedPtr< EBIndexSpace > & getEBIndexSpace() const
Return the EB index space.
Definition CD_PhaseRealm.cpp:840
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:349
Vector< RefCountedPtr< EBGradient > > m_gradientOp
Gradient operator.
Definition CD_PhaseRealm.H:521
void preRegrid()
Perform pre-regrid operations.
Definition CD_PhaseRealm.cpp:115
RefCountedPtr< EBIndexSpace > m_ebis
EB index space.
Definition CD_PhaseRealm.H:434
int m_numLsfGhostCells
Number of ghost cells in level-set function (on the mesh)
Definition CD_PhaseRealm.H:384
bool m_profile
Profile or not.
Definition CD_PhaseRealm.H:359
void defineVofIterator(const int a_lmin)
Define vof iterators.
Definition CD_PhaseRealm.cpp:427
int m_redistributionRadius
Redistribution radius.
Definition CD_PhaseRealm.H:389
Vector< ProblemDomain > m_domains
Problem domains.
Definition CD_PhaseRealm.H:454
const Vector< ProblemDomain > & getDomains() const
Get problem domains.
Definition CD_PhaseRealm.cpp:864
int m_numEbGhostsCells
Number of ghost cells in EBLevelGrid/EBISLayout/EBISBox.
Definition CD_PhaseRealm.H:379
void regridOperators(const int a_lmin)
Regrid method for EBAMR operators.
Definition CD_PhaseRealm.cpp:183
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
TracerParticleSolver()
Default constructor.
Definition CD_TracerParticleSolverImplem.H:25