chombo-discharge
Loading...
Searching...
No Matches
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
39class MFHelmholtzOpFactory : public AMRLevelOpFactory<LevelData<MFCellFAB>>
40{
41public:
42 // Various alias for cutting down on typing.
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>;
51
55
56 using DomainBCFactory = RefCountedPtr<MFHelmholtzDomainBCFactory>;
57 using EBBCFactory = RefCountedPtr<MFHelmholtzEBBCFactory>;
58 using JumpBCFactory = RefCountedPtr<MFHelmholtzJumpBCFactory>;
59
61
66
71
100 MFHelmholtzOpFactory(const MFIS& a_mfis,
102 const Real& a_alpha,
103 const Real& a_beta,
104 const RealVect& a_probLo,
105 const AmrLevelGrids& a_amrLevelGrids,
106 const AmrMask& a_validCells,
107 const AmrInterpolators& a_amrInterpolators,
108 const AmrFluxRegisters& a_amrFluxRegisters,
109 const AmrCoarseners& a_amrCoarseners,
110 const AmrRefRatios& a_amrRefRatios,
111 const AmrResolutions& a_amrResolutions,
112 const AmrCellData& a_amrAcoef,
113 const AmrFluxData& a_amrBcoef,
114 const AmrIrreData& a_amrBcoefIrreg,
115 const DomainBCFactory& a_domainBcFactory,
116 const EBBCFactory& a_ebBcFactory,
117 const JumpBCFactory& a_jumpBcFactory,
118 const IntVect& a_ghostPhi,
119 const IntVect& a_ghostRhs,
120 const Smoother& a_smoother,
121 const Real& a_relaxFactor,
123 const int& a_jumpOrder,
124 const int& a_jumpWeight,
125 const int& a_preCondSmooth,
126 const int& a_blockingFactor,
127 const AmrLevelGrids& a_deeperLevelGrids = AmrLevelGrids());
128
132 virtual ~MFHelmholtzOpFactory();
133
139 void
140 setJump(const EBAMRIVData& a_sigma, const Real& a_scale);
141
146 void
147 setJump(const Real& a_sigma, const Real& a_scale);
148
152 const EBAMRIVData&
153 getSigma() const;
154
158 void
160
169
176
181 int
183
189
194
199
203 MFIS m_mfis;
204
209
214
219
224
229
234
239
244
249
254
258 AmrLevelGrids m_amrLevelGrids;
259
264
268 AmrInterpolators m_amrInterpolators;
269
273 AmrFluxRegisters m_amrFluxRegisters;
274
278 AmrCoarseners m_amrCoarseners;
279
283 AmrRefRatios m_amrRefRatios;
284
288 AmrResolutions m_amrResolutions;
289
293 AmrCellData m_amrAcoef;
294
298 AmrFluxData m_amrBcoef;
299
303 AmrIrreData m_amrBcoefIrreg;
304
308 DomainBCFactory m_domainBcFactory;
309
313 EBBCFactory m_ebBcFactory;
314
318 JumpBCFactory m_jumpBcFactory;
319
324
329
334
339
344
348 AmrLevelGrids m_deeperLevelGrids;
349
354
358 Vector<AmrLevelGrids> m_mgLevelGrids;
359
363 Vector<AmrCellData> m_mgAcoef;
364
368 Vector<AmrFluxData> m_mgBcoef;
369
373 Vector<AmrIrreData> m_mgBcoefIrreg;
374
379
384
388 void
390
394 void
395 defineJump();
396
403 bool
405
412 bool
414
424 bool
427 const int a_refRat,
429
442 void
451 const int a_refRat);
452
459 int
461};
462
463#include <CD_NamespaceFooter.H>
464
465#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:188
ProblemDomain m_bottomDomain
Coarsest level where we run multigrid.
Definition CD_MFHelmholtzOpFactory.H:328
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:432
Vector< EBAMRIVData > m_mgJump
Jump coefficient on deeper grids.
Definition CD_MFHelmholtzOpFactory.H:378
EBBCFactory m_ebBcFactory
EB BC factory.
Definition CD_MFHelmholtzOpFactory.H:313
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:368
std::vector< bool > m_hasMgLevels
For checking if an AMR level has multigrid levels.
Definition CD_MFHelmholtzOpFactory.H:353
EBAMRIVData m_amrJump
Jump on multiphase cells.
Definition CD_MFHelmholtzOpFactory.H:323
AmrFluxData m_amrBcoef
Operator B-coefficient.
Definition CD_MFHelmholtzOpFactory.H:298
int m_mgBlockingFactor
Blocking factor for grid aggregration.
Definition CD_MFHelmholtzOpFactory.H:333
AmrFluxRegisters m_amrFluxRegisters
Flux registers.
Definition CD_MFHelmholtzOpFactory.H:273
bool isFiner(const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const
Check if a domain is finer than the other.
Definition CD_MFHelmholtzOpFactory.cpp:809
Location::Cell m_dataLocation
Data centering.
Definition CD_MFHelmholtzOpFactory.H:208
MFHelmholtzOp * AMRnewOp(const ProblemDomain &a_domain) override final
Create AMR operator for specified domain.
Definition CD_MFHelmholtzOpFactory.cpp:716
MFHelmholtzOpFactory()=delete
Disallowed constructor.
int m_numAmrLevels
Number of AMR levels.
Definition CD_MFHelmholtzOpFactory.H:218
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:358
Real m_alpha
Operator alpha.
Definition CD_MFHelmholtzOpFactory.H:238
int m_jumpWeight
Weighting for least squares reconstruction in jump cells.
Definition CD_MFHelmholtzOpFactory.H:343
int refToFiner(const ProblemDomain &a_indexspace) const override final
Get refinement ratio to next finest level.
Definition CD_MFHelmholtzOpFactory.cpp:817
AmrCoarseners m_amrCoarseners
Data coarseners.
Definition CD_MFHelmholtzOpFactory.H:278
const EBAMRIVData & getSigma() const
Get the BC jump factor.
Definition CD_MFHelmholtzOpFactory.cpp:190
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:363
AmrIrreData m_amrBcoefIrreg
Operator B-coefficient.
Definition CD_MFHelmholtzOpFactory.H:303
JumpBCFactory m_jumpBcFactory
Jump BC factory.
Definition CD_MFHelmholtzOpFactory.H:318
AmrLevelGrids m_amrLevelGrids
AMR grids.
Definition CD_MFHelmholtzOpFactory.H:258
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:496
static constexpr int m_nComp
Number of components that this operator solves for.
Definition CD_MFHelmholtzOpFactory.H:193
MFIS m_mfis
Index space.
Definition CD_MFHelmholtzOpFactory.H:203
IntVect m_ghostPhi
Number of ghost cells that are used. Need because of Chombo prolongation objects.
Definition CD_MFHelmholtzOpFactory.H:228
AmrMask m_validCells
Valid cells.
Definition CD_MFHelmholtzOpFactory.H:263
Real m_relaxFactor
Relaxation factor.
Definition CD_MFHelmholtzOpFactory.H:248
AmrResolutions m_amrResolutions
Grid resolutions.
Definition CD_MFHelmholtzOpFactory.H:288
static constexpr int m_mainPhase
Main phase.
Definition CD_MFHelmholtzOpFactory.H:198
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:118
DomainBCFactory m_domainBcFactory
Domain BC factory.
Definition CD_MFHelmholtzOpFactory.H:308
int findAmrLevel(const ProblemDomain &a_domain) const
Find level corresponding to amr level.
Definition CD_MFHelmholtzOpFactory.cpp:835
Smoother m_smoother
Smoother.
Definition CD_MFHelmholtzOpFactory.H:213
RealVect m_probLo
Lower-left corner of computational domain.
Definition CD_MFHelmholtzOpFactory.H:253
AmrRefRatios m_amrRefRatios
Refinement ratios.
Definition CD_MFHelmholtzOpFactory.H:283
AmrLevelGrids m_deeperLevelGrids
Deeper grids.
Definition CD_MFHelmholtzOpFactory.H:348
Vector< Vector< RefCountedPtr< EBCoarAve > > > m_mgAveOp
Coarsneing operator on deeper grids (used for coarsening jump)
Definition CD_MFHelmholtzOpFactory.H:383
IntVect m_ghostRhs
Number of ghost cells that are used. Need because of Chombo prolongation objects.
Definition CD_MFHelmholtzOpFactory.H:233
bool isCoarser(const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const
Check if a domain is coarser than the other.
Definition CD_MFHelmholtzOpFactory.cpp:801
MFHelmholtzOp * MGnewOp(const ProblemDomain &a_fineDomain, int a_depth, bool a_homogeneousOnly=true) override final
Create multigrid operator.
Definition CD_MFHelmholtzOpFactory.cpp:569
void coarsenCoefficientsMG()
Go through all MG levels and coarsen the coefficients from the finer levels.
Definition CD_MFHelmholtzOpFactory.cpp:389
AmrCellData m_amrAcoef
Operator A-coefficient.
Definition CD_MFHelmholtzOpFactory.H:293
AmrInterpolators m_amrInterpolators
Multigrid interpolators.
Definition CD_MFHelmholtzOpFactory.H:268
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:373
Real m_beta
Operator beta.
Definition CD_MFHelmholtzOpFactory.H:243
void defineMultigridLevels()
Function which defines the multigrid levels for this operator factory.
Definition CD_MFHelmholtzOpFactory.cpp:235
int m_numPreCondSmooth
Number of smoothings in the preconditioner. Passed directly to the operators.
Definition CD_MFHelmholtzOpFactory.H:223
void defineJump()
Define jump data.
Definition CD_MFHelmholtzOpFactory.cpp:198
virtual ~MFHelmholtzOpFactory()
Destructor. Does nothing.
Definition CD_MFHelmholtzOpFactory.cpp:112
int m_jumpOrder
Stencil order in jump cells.
Definition CD_MFHelmholtzOpFactory.H:338
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
Base class for a tracer particle solver. This solver can advance particles in a pre-defined velocity ...
Definition CD_TracerParticleSolver.H:37
Namespace for encapsulating various data centerings.
Definition CD_Location.H:24
Cell
Enum for distinguishing between cell locations.
Definition CD_Location.H:30