Factory class for making MFHelmholtzOp.
More...
#include <CD_MFHelmholtzOpFactory.H>
|
| MFHelmholtzOpFactory ()=delete |
| Disallowed constructor.
|
|
| MFHelmholtzOpFactory (const MFHelmholtzOpFactory &a_otherFactory)=delete |
| Disallowed constructor.
|
|
| MFHelmholtzOpFactory (const MFIS &a_mfis, const Location::Cell a_dataLocation, const Real &a_alpha, const Real &a_beta, const RealVect &a_probLo, const AmrLevelGrids &a_amrLevelGrids, const AmrInterpolators &a_amrInterpolators, const AmrFluxRegisters &a_amrFluxRegisters, const AmrCoarseners &a_amrCoarseners, const AmrRefRatios &a_amrRefRatios, const AmrResolutions &a_amrResolutions, const AmrCellData &a_amrAcoef, const AmrFluxData &a_amrBcoef, const AmrIrreData &a_amrBcoefIrreg, const DomainBCFactory &a_domainBcFactory, const EBBCFactory &a_ebBcFactory, const JumpBCFactory &a_jumpBcFactory, const IntVect &a_ghostPhi, const IntVect &a_ghostRhs, const Smoother &a_smoother, const ProblemDomain &a_bottomDomain, const int &a_jumpOrder, const int &a_jumpWeight, const int &a_blockingFactor, const AmrLevelGrids &a_deeperLevelGrids=AmrLevelGrids()) |
| Full constructor. More...
|
|
| ~MFHelmholtzOpFactory () |
| Destructor. Does nothing.
|
|
void | setJump (const EBAMRIVData &a_sigma, const Real &a_scale) |
| Set jump condition. This is passed to the operators by reference. More...
|
|
void | setJump (const Real &a_sigma, const Real &a_scale) |
|
const EBAMRIVData & | getSigma () const |
| Get the BC jump factor.
|
|
void | coarsenCoefficientsMG () |
| Go through all MG levels and coarsen the coefficients from the finer levels.
|
|
MFHelmholtzOp * | MGnewOp (const ProblemDomain &a_fineDomain, int a_depth, bool a_homogeneousOnly=true) override final |
| Create multigrid operator. More...
|
|
MFHelmholtzOp * | AMRnewOp (const ProblemDomain &a_domain) override final |
| Create AMR operator for specified domain. More...
|
|
int | refToFiner (const ProblemDomain &a_indexspace) const override final |
| Get refinement ratio to next finest level. More...
|
|
|
void | defineMultigridLevels () |
| Function which defines the multigrid levels for this operator factory.
|
|
void | defineJump () |
| Define jump data.
|
|
bool | isCoarser (const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const |
| Check if a domain is coarser than the other. More...
|
|
bool | isFiner (const ProblemDomain &a_domainOne, const ProblemDomain &a_domainTwo) const |
| Check if a domain is finer than the other. More...
|
|
bool | getCoarserLayout (MFLevelGrid &a_coarseGrid, const MFLevelGrid &a_fineGrid, const int a_refRat, const int a_blockingFactor) const |
| Construct coarsening of a grid level. More...
|
|
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) More...
|
|
int | findAmrLevel (const ProblemDomain &a_domain) const |
| Find level corresponding to amr level. More...
|
|
|
static constexpr int | m_comp = 0 |
| Component that this operator solves for.
|
|
static constexpr int | m_nComp = 1 |
| Number of components that this operator solves for.
|
|
static constexpr int | m_mainPhase = 0 |
| Main phase.
|
|
Factory class for making MFHelmholtzOp.
- Note
- This factory is designed for coupling EBHelmholtzOps across EBIS levels
-
This factory is designed without time-dependence in BCs. Time-dependent BCs are still doable by letting the boundary condition classes carry a reference to an externally updated time. If you need new coefficients you will have to set up multigrid again.
◆ MFHelmholtzOpFactory()
MFHelmholtzOpFactory::MFHelmholtzOpFactory |
( |
const MFIS & |
a_mfis, |
|
|
const Location::Cell |
a_dataLocation, |
|
|
const Real & |
a_alpha, |
|
|
const Real & |
a_beta, |
|
|
const RealVect & |
a_probLo, |
|
|
const AmrLevelGrids & |
a_amrLevelGrids, |
|
|
const AmrInterpolators & |
a_amrInterpolators, |
|
|
const AmrFluxRegisters & |
a_amrFluxRegisters, |
|
|
const AmrCoarseners & |
a_amrCoarseners, |
|
|
const AmrRefRatios & |
a_amrRefRatios, |
|
|
const AmrResolutions & |
a_amrResolutions, |
|
|
const AmrCellData & |
a_amrAcoef, |
|
|
const AmrFluxData & |
a_amrBcoef, |
|
|
const AmrIrreData & |
a_amrBcoefIrreg, |
|
|
const DomainBCFactory & |
a_domainBcFactory, |
|
|
const EBBCFactory & |
a_ebBcFactory, |
|
|
const JumpBCFactory & |
a_jumpBcFactory, |
|
|
const IntVect & |
a_ghostPhi, |
|
|
const IntVect & |
a_ghostRhs, |
|
|
const Smoother & |
a_smoother, |
|
|
const ProblemDomain & |
a_bottomDomain, |
|
|
const int & |
a_jumpOrder, |
|
|
const int & |
a_jumpWeight, |
|
|
const int & |
a_blockingFactor, |
|
|
const AmrLevelGrids & |
a_deeperLevelGrids = AmrLevelGrids() |
|
) |
| |
Full constructor.
- Parameters
-
[in] | a_mfis | Multi-fluid index space |
[in] | a_dataLocation | Assumed data centering |
[in] | a_alpha | alpha-coefficient in Helmholtz operator. |
[in] | a_beta | beta-coefficient in Helmholtz operator. |
[in] | a_probLo | Lower-left corner of domain |
[in] | a_amrInterpolators | Interpolator objects between AMR levels. |
[in] | a_amrCoarseners | Conservative coarseners between AMR levels. |
[in] | a_amrFluxRegisters | Flux registers between AMR levels. |
[in] | a_amrResolutions | Grid resolutions for AMR levels. |
[in] | a_amrAcoef | A-coefficient in Helmholtz operator. |
[in] | a_amrBcoef | B-coefficient in Helmholtz operator. |
[in] | a_amrBcoefIrreg | B-coefficient in Helmholtz operator. This one is defined on EB faces. |
[in] | a_domainBCFactory | Factory class for making domain BC objects. |
[in] | a_ebbcFactory | Factory class for making BC objects for EB boundary conditions. |
[in] | a_jumpBCFactory | Factory class for making jump BC objects. |
[in] | a_ghostPhi | Number of ghost cells in solution vector. |
[in] | a_ghostRhs | Number of ghost cells in right-hand side. |
[in] | a_jumpOrder | Stencil order in multiphase cells |
[in] | a_jumpWeight | Equation weighting in multiphase cells |
[in] | a_relaxationMethod | Relaxation method. |
[in] | a_bottomDomain | Coarsest domain on which we run multigrid. Must be a coarsening of the AMR problem domains. |
[in] | a_deeperLevelGrids | Optional object in case you want to pre-define the deeper multigrid levels. |
◆ AMRnewOp()
MFHelmholtzOp * MFHelmholtzOpFactory::AMRnewOp |
( |
const ProblemDomain & |
a_domain | ) |
|
|
finaloverride |
Create AMR operator for specified domain.
- Parameters
-
◆ coarsenCoefficients()
void MFHelmholtzOpFactory::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 |
|
) |
| |
|
protected |
Coarsen coefficients (conservatively)
- Parameters
-
[out] | a_coarAcoef | Coarse A-coefficient |
[out] | a_coarBcoef | Coarse B-coefficient |
[out] | a_coarBcoefIrreg | Coarse B-coefficient on EB faces |
[in] | a_fineAcoef | Fine A-coefficient |
[in] | a_fineBcoef | Fine B-coefficient |
[in] | a_fineBcoefIrreg | Fine B-coefficient on EB faces |
[in] | a_eblgCoar | Coarse grids |
[in] | a_eblgFine | Fine grids |
[in] | a_refRat | Coarsening factor |
◆ findAmrLevel()
int MFHelmholtzOpFactory::findAmrLevel |
( |
const ProblemDomain & |
a_domain | ) |
const |
|
protected |
Find level corresponding to amr level.
- Parameters
-
[in] | a_domain | Problem domain. |
- Returns
- Depth in m_amrLevelGrids corresponding to a_domain.
- Note
- Run-time error if no level was found.
◆ getCoarserLayout()
bool MFHelmholtzOpFactory::getCoarserLayout |
( |
MFLevelGrid & |
a_coarseGrid, |
|
|
const MFLevelGrid & |
a_fineGrid, |
|
|
const int |
a_refRat, |
|
|
const int |
a_blockingFactor |
|
) |
| const |
|
protected |
Construct coarsening of a grid level.
- Parameters
-
[out] | a_coarseGrid | The coarse grid layout. Must be a pointer to an undefined EBLevelGrid on input |
[in] | a_fineGrid | The coarse grid layout. Must be a pointer to an undefined EBLevelGrid on input |
[in] | a_refRat | Refinement ratio |
[in] | a_blockingFactor | Blocking factor to use for grid aggregation |
- Returns
- This will return a multigrid level (i.e. one that is completely overlapping) the fine level. If we can, we use aggregation with the blocking factor. Otherwise we try to coarsen directly.
◆ isCoarser()
bool MFHelmholtzOpFactory::isCoarser |
( |
const ProblemDomain & |
a_domainOne, |
|
|
const ProblemDomain & |
a_domainTwo |
|
) |
| const |
|
protected |
Check if a domain is coarser than the other.
- Parameters
-
[in] | a_domainOne | The first domain |
[in] | a_domainTwo | The second domain |
- Returns
- Returns true of a_domainOne has fewer grid points than a_domainTwo
◆ isFiner()
bool MFHelmholtzOpFactory::isFiner |
( |
const ProblemDomain & |
a_domainOne, |
|
|
const ProblemDomain & |
a_domainTwo |
|
) |
| const |
|
protected |
Check if a domain is finer than the other.
- Parameters
-
[in] | a_domainOne | The first domain |
[in] | a_domainTwo | The second domain |
- Returns
- Returns true of a_domainOne has more grid points than a_domainTwo
◆ MGnewOp()
MFHelmholtzOp * MFHelmholtzOpFactory::MGnewOp |
( |
const ProblemDomain & |
a_fineDomain, |
|
|
int |
a_depth, |
|
|
bool |
a_homogeneousOnly = true |
|
) |
| |
|
finaloverride |
Create multigrid operator.
- Parameters
-
[in] | a_fineDomain | Domain |
[in] | a_depth | Depth. This specifies that the operator will be created at depth coarsen(a_fineDomain, 2^a_depth); |
[in] | a_homogeneousOnly | If true, only homogeneous boundary conditions will be needed. |
◆ refToFiner()
int MFHelmholtzOpFactory::refToFiner |
( |
const ProblemDomain & |
a_indexspace | ) |
const |
|
finaloverride |
Get refinement ratio to next finest level.
- Note
- Returns -1 when there are no finer levels.
◆ setJump() [1/2]
void MFHelmholtzOpFactory::setJump |
( |
const EBAMRIVData & |
a_sigma, |
|
|
const Real & |
a_scale |
|
) |
| |
Set jump condition. This is passed to the operators by reference.
- Parameters
-
[in] | a_sigma | Jump (must be defined on gas phase) |
[in] | a_scale | Scaling factor |
◆ setJump() [2/2]
void MFHelmholtzOpFactory::setJump |
( |
const Real & |
a_sigma, |
|
|
const Real & |
a_scale |
|
) |
| |
- Parameters
-
[in] | a_sigma | Jump |
[in] | a_scale | Scaling factor |
The documentation for this class was generated from the following files: