chombo-discharge
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes | List of all members
MFHelmholtzOpFactory Class Reference

Factory class for making MFHelmholtzOp. More...

#include <CD_MFHelmholtzOpFactory.H>

Inheritance diagram for MFHelmholtzOpFactory:
Inheritance graph
[legend]
Collaboration diagram for MFHelmholtzOpFactory:
Collaboration graph
[legend]

Public Types

using MFIS = RefCountedPtr< MultiFluidIndexSpace >
 
using AmrLevelGrids = Vector< MFLevelGrid >
 
using AmrInterpolators = Vector< MFMultigridInterpolator >
 
using AmrFluxRegisters = Vector< MFReflux >
 
using AmrCoarseners = Vector< MFCoarAve >
 
using AmrResolutions = Vector< Real >
 
using AmrRefRatios = Vector< int >
 
using AmrCellData = Vector< RefCountedPtr< LevelData< MFCellFAB > >>
 
using AmrFluxData = Vector< RefCountedPtr< LevelData< MFFluxFAB > >>
 
using AmrIrreData = Vector< RefCountedPtr< LevelData< MFBaseIVFAB > >>
 
using DomainBCFactory = RefCountedPtr< MFHelmholtzDomainBCFactory >
 
using EBBCFactory = RefCountedPtr< MFHelmholtzEBBCFactory >
 
using JumpBCFactory = RefCountedPtr< MFHelmholtzJumpBCFactory >
 
using Smoother = MFHelmholtzOp::Smoother
 

Public Member Functions

 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 EBAMRIVDatagetSigma () const
 Get the BC jump factor.
 
void coarsenCoefficientsMG ()
 Go through all MG levels and coarsen the coefficients from the finer levels.
 
MFHelmholtzOpMGnewOp (const ProblemDomain &a_fineDomain, int a_depth, bool a_homogeneousOnly=true) override final
 Create multigrid operator. More...
 
MFHelmholtzOpAMRnewOp (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...
 

Protected Member Functions

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...
 

Protected Attributes

MFIS m_mfis
 Index space.
 
Location::Cell m_dataLocation
 Data centering.
 
Smoother m_smoother
 Smoother.
 
int m_numAmrLevels
 Number of AMR levels.
 
IntVect m_ghostPhi
 Number of ghost cells that are used. Need because of Chombo prolongation objects.
 
IntVect m_ghostRhs
 Number of ghost cells that are used. Need because of Chombo prolongation objects.
 
Real m_alpha
 Operator alpha.
 
Real m_beta
 Operator beta.
 
RealVect m_probLo
 Lower-left corner of computational domain.
 
AmrLevelGrids m_amrLevelGrids
 AMR grids.
 
AmrInterpolators m_amrInterpolators
 Multigrid interpolators.
 
AmrFluxRegisters m_amrFluxRegisters
 Flux registers.
 
AmrCoarseners m_amrCoarseners
 Data coarseners.
 
AmrRefRatios m_amrRefRatios
 Refinement ratios.
 
AmrResolutions m_amrResolutions
 Grid resolutions.
 
AmrCellData m_amrAcoef
 Operator A-coefficient.
 
AmrFluxData m_amrBcoef
 Operator B-coefficient.
 
AmrIrreData m_amrBcoefIrreg
 Operator B-coefficient.
 
DomainBCFactory m_domainBcFactory
 Domain BC factory.
 
EBBCFactory m_ebBcFactory
 EB BC factory.
 
JumpBCFactory m_jumpBcFactory
 Jump BC factory.
 
EBAMRIVData m_amrJump
 Jump on multiphase cells.
 
ProblemDomain m_bottomDomain
 Coarsest level where we run multigrid.
 
int m_mgBlockingFactor
 Blocking factor for grid aggregration.
 
int m_jumpOrder
 Stencil order in jump cells.
 
int m_jumpWeight
 Weighting for least squares reconstruction in jump cells.
 
AmrLevelGrids m_deeperLevelGrids
 Deeper grids.
 
std::vector< bool > m_hasMgLevels
 For checking if an AMR level has multigrid levels.
 
Vector< AmrLevelGrids > m_mgLevelGrids
 Deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to the the multigrid levels below amr level 0.
 
Vector< AmrCellData > m_mgAcoef
 A-coefficient on deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to the the multigrid levels below amr level 0.
 
Vector< AmrFluxData > m_mgBcoef
 B-coefficient on deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to the the multigrid levels below amr level 0.
 
Vector< AmrIrreData > m_mgBcoefIrreg
 B-coefficient on deeper grids. Always weird to write this but e.g. m_mgLevelGrids[0] corresponds to the the multigrid levels below amr level 0.
 
Vector< EBAMRIVDatam_mgJump
 Jump coefficient on deeper grids.
 
Vector< Vector< RefCountedPtr< EBCoarAve > > > m_mgAveOp
 Coarsneing operator on deeper grids (used for coarsening jump)
 

Static Protected Attributes

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.
 

Detailed Description

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.

Constructor & Destructor Documentation

◆ 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_mfisMulti-fluid index space
[in]a_dataLocationAssumed data centering
[in]a_alphaalpha-coefficient in Helmholtz operator.
[in]a_betabeta-coefficient in Helmholtz operator.
[in]a_probLoLower-left corner of domain
[in]a_amrInterpolatorsInterpolator objects between AMR levels.
[in]a_amrCoarsenersConservative coarseners between AMR levels.
[in]a_amrFluxRegistersFlux registers between AMR levels.
[in]a_amrResolutionsGrid resolutions for AMR levels.
[in]a_amrAcoefA-coefficient in Helmholtz operator.
[in]a_amrBcoefB-coefficient in Helmholtz operator.
[in]a_amrBcoefIrregB-coefficient in Helmholtz operator. This one is defined on EB faces.
[in]a_domainBCFactoryFactory class for making domain BC objects.
[in]a_ebbcFactoryFactory class for making BC objects for EB boundary conditions.
[in]a_jumpBCFactoryFactory class for making jump BC objects.
[in]a_ghostPhiNumber of ghost cells in solution vector.
[in]a_ghostRhsNumber of ghost cells in right-hand side.
[in]a_jumpOrderStencil order in multiphase cells
[in]a_jumpWeightEquation weighting in multiphase cells
[in]a_relaxationMethodRelaxation method.
[in]a_bottomDomainCoarsest domain on which we run multigrid. Must be a coarsening of the AMR problem domains.
[in]a_deeperLevelGridsOptional object in case you want to pre-define the deeper multigrid levels.

Member Function Documentation

◆ AMRnewOp()

MFHelmholtzOp * MFHelmholtzOpFactory::AMRnewOp ( const ProblemDomain &  a_domain)
finaloverride

Create AMR operator for specified domain.

Parameters
[in]a_domainDomain

◆ 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_coarAcoefCoarse A-coefficient
[out]a_coarBcoefCoarse B-coefficient
[out]a_coarBcoefIrregCoarse B-coefficient on EB faces
[in]a_fineAcoefFine A-coefficient
[in]a_fineBcoefFine B-coefficient
[in]a_fineBcoefIrregFine B-coefficient on EB faces
[in]a_eblgCoarCoarse grids
[in]a_eblgFineFine grids
[in]a_refRatCoarsening factor

◆ findAmrLevel()

int MFHelmholtzOpFactory::findAmrLevel ( const ProblemDomain &  a_domain) const
protected

Find level corresponding to amr level.

Parameters
[in]a_domainProblem 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_coarseGridThe coarse grid layout. Must be a pointer to an undefined EBLevelGrid on input
[in]a_fineGridThe coarse grid layout. Must be a pointer to an undefined EBLevelGrid on input
[in]a_refRatRefinement ratio
[in]a_blockingFactorBlocking 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_domainOneThe first domain
[in]a_domainTwoThe 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_domainOneThe first domain
[in]a_domainTwoThe 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_fineDomainDomain
[in]a_depthDepth. This specifies that the operator will be created at depth coarsen(a_fineDomain, 2^a_depth);
[in]a_homogeneousOnlyIf 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_sigmaJump (must be defined on gas phase)
[in]a_scaleScaling factor

◆ setJump() [2/2]

void MFHelmholtzOpFactory::setJump ( const Real &  a_sigma,
const Real &  a_scale 
)
Parameters
[in]a_sigmaJump
[in]a_scaleScaling factor

The documentation for this class was generated from the following files: