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

Helmholtz operator for equations like alpha*a(x)*phi(x) + beta*div(b(x)*grad(phi(x))) = rho. More...

#include <CD_EBHelmholtzOp.H>

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

Public Types

enum class  Smoother { NoRelax , PointJacobi , GauSaiRedBlack , GauSaiMultiColor }
 Relaxation method for the operators.
 

Public Member Functions

 EBHelmholtzOp ()=delete
 Disallowed default constructor.
 
 EBHelmholtzOp (const EBHelmholtzOp &a_other)=delete
 Disallowed copy constructor. More...
 
 EBHelmholtzOp (const EBHelmholtzOp &&a_other)=delete
 Disallowed move constructor. More...
 
 EBHelmholtzOp (const Location::Cell a_dataLocation, const EBLevelGrid &a_eblgFine, const EBLevelGrid &a_eblg, const EBLevelGrid &a_eblgCoFi, const EBLevelGrid &a_eblgCoar, const EBLevelGrid &a_eblgCoarMG, const RefCountedPtr< EBMultigridInterpolator > &a_interpolator, const RefCountedPtr< EBReflux > &a_fluxReg, const RefCountedPtr< EBCoarAve > &a_coarAve, const RefCountedPtr< EBHelmholtzDomainBC > &a_domainBC, const RefCountedPtr< EBHelmholtzEBBC > &a_ebBC, const RealVect &a_probLo, const Real &a_dx, const int &a_refToFine, const int &a_refToCoar, const bool &a_hasFine, const bool &a_hasCoar, const bool &a_hasMGObjects, const Real &a_alpha, const Real &a_beta, const RefCountedPtr< LevelData< EBCellFAB >> &a_Acoef, const RefCountedPtr< LevelData< EBFluxFAB >> &a_Bcoef, const RefCountedPtr< LevelData< BaseIVFAB< Real >>> &a_BcoIrreg, const IntVect &a_ghostCellsPhi, const IntVect &a_ghostCellsRHS, const Smoother &a_smoother)
 Full constructor. More...
 
virtual ~EBHelmholtzOp ()
 Dtor.
 
EBHelmholtzOpoperator= (const EBHelmholtzOp &a_oper)=delete
 No copy assigment allowed. More...
 
EBHelmholtzOpoperator= (const EBHelmholtzOp &&a_oper)=delete
 No move assigment allowed. More...
 
void turnOffCFInterp ()
 Turn off BCs.
 
void turnOnCFInterp ()
 Turn on BCs.
 
void turnOffExchange ()
 Turn off exchange operation. More...
 
void turnOnExchange ()
 Turn on exchange operation. More...
 
void turnOffCoarsening ()
 Turn off coarsening operation. More...
 
void turnOnCoarsening ()
 Turn on coarsening operation. More...
 
void setAcoAndBco (const RefCountedPtr< LevelData< EBCellFAB >> &a_Acoef, const RefCountedPtr< LevelData< EBFluxFAB >> &a_Bcoef, const RefCountedPtr< LevelData< BaseIVFAB< Real >>> &a_BcoefIrreg)
 Update with new A and B coefficients. More...
 
const RefCountedPtr< LevelData< EBCellFAB > > & getAcoef ()
 Get the Helmholtz A-coefficient on cell centers. More...
 
const RefCountedPtr< LevelData< EBFluxFAB > > & getBcoef ()
 Get the Helmholtz B-coefficient on faces. More...
 
const RefCountedPtr< LevelData< BaseIVFAB< Real > > > & getBcoefIrreg ()
 Get the Helmholtz B-coefficient on the EB. More...
 
void coarsenCell (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiFine)
 Coarsen data from fine to coar level. More...
 
void coarsenFlux (LevelData< EBFluxFAB > &a_flux, const LevelData< EBFluxFAB > &a_fineFlux)
 Coarsen fluxes on the fine level onto this level. More...
 
void pointJacobiKernel (EBCellFAB &a_Lcorr, EBCellFAB &a_corr, const EBCellFAB &a_resid, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit) const noexcept
 Point Jacobi kernel. More...
 
void gauSaiRedBlackKernel (EBCellFAB &a_Lcorr, EBCellFAB &a_corr, const EBCellFAB &a_resid, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit, const int &a_redBlack) const noexcept
 Red-black Gauss-Seidel kernel. More...
 
void gauSaiMultiColorKernel (EBCellFAB &a_Lcorr, EBCellFAB &a_corr, const EBCellFAB &a_resid, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit, const IntVect &a_color) const noexcept
 Multi-color Gauss-Seidel kernel. More...
 
void residual (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, const bool a_homogeneousPhysBc) override
 Compute residual on this level. More...
 
void preCond (LevelData< EBCellFAB > &a_corr, const LevelData< EBCellFAB > &a_residual) override final
 Precondition system before bottom solve. More...
 
void interpolateCF (LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > *a_phiCoar, const bool a_homogeneousCFBC)
 Apply coarse-fine boundary conditions. More...
 
void homogeneousCFInterp (LevelData< EBCellFAB > &a_phi)
 Do homogeneous coarse-fine interpolation. More...
 
void inhomogeneousCFInterp (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar)
 Inhomogeneous coarse-fine interpolation. More...
 
void gauSaiMultiColor (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_rhs, const IntVect &a_color) const
 Multi-colored Gauss-Seidel kernel. Public because MFHelmholtzOp may want to use use. More...
 
void applyOp (LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoar, const bool a_homogeneousPhysBC, const bool a_homogeneousCFBC)
 Apply operator on this level. This is a more general version which can turn on/off homogeneous and CF bcs. More...
 
void applyOp (LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousPhysBc) override final
 Apply operator. More...
 
void applyOp (EBCellFAB &a_Lphi, EBCellFAB &a_phi, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit, const bool a_homogeneousPhysBC) const noexcept
 Apply operator in a grid box. More...
 
void applyOpRegular (EBCellFAB &a_Lphi, EBCellFAB &a_phi, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const Box &a_cellBox, const DataIndex &a_dit, const bool a_homogeneousPhysBC) const noexcept
 Apply operator in regular cells. More...
 
void applyDomainFlux (EBCellFAB &a_phi, const EBFluxFAB &a_Bcoef, const Box &a_cellBox, const DataIndex &a_dit, const bool a_homogeneousPhysBc) const noexcept
 Apply domain flux. More...
 
void fillDomainFlux (EBFluxFAB &a_flux, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit)
 Fill domain flux. This fills the flux on the domain face using centered differencing ala applyDomainFlux. More...
 
void applyOpIrregular (EBCellFAB &a_Lphi, const EBCellFAB &a_phi, const EBCellFAB &a_Acoef, const EBFluxFAB &a_Bcoef, const BaseIVFAB< Real > &a_BcoefIrreg, const BaseIVFAB< Real > &a_alphaDiagWeight, const Box &a_cellBox, const DataIndex &a_dit, const bool a_homogeneousPhysBC) const noexcept
 Apply operator in irregular cells. More...
 
void create (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
 Create data which clones the layout of the other. More...
 
void assign (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
 Assign data. More...
 
void assignCopier (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const Copier &a_copier) override final
 Assign lhs. More...
 
void assignLocal (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
 Local assignment function. More...
 
void buildCopier (Copier &a_copier, const LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override
 Build copier. More...
 
Real dotProduct (const LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs) override final
 Compute the dot product??
 
void incr (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const Real a_scale) override final
 Increment operator. More...
 
void axby (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, const LevelData< EBCellFAB > &a_y, const Real a_a, const Real a_b) override final
 Set a_lhs = a*x + b*y. More...
 
void scale (LevelData< EBCellFAB > &a_lhs, const Real &a_scale) override final
 Scale data. Returns a_lhs = a_lhs*a_scale. More...
 
Real norm (const LevelData< EBCellFAB > &a_rhs, const int a_order) override final
 Compute norm of data. More...
 
void setToZero (LevelData< EBCellFAB > &a_lhs) override final
 Set data to zero. More...
 
void createCoarser (LevelData< EBCellFAB > &a_coarse, const LevelData< EBCellFAB > &a_fine, bool a_ghosted) override final
 Create coarsened data. More...
 
void relax (LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, int a_iterations) override final
 Relaxation method. This does smoothing for the system L(correction) = residual. More...
 
void restrictResidual (LevelData< EBCellFAB > &a_resCoar, LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs) override final
 Restrict residual onto coarse level. More...
 
void prolongIncrement (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_correctCoarse) override final
 Prolongation method. More...
 
int refToCoarser () override final
 Return coarsening factor to coarser level (1 if there is no coarser level);.
 
void AMROperator (LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB >> *a_finerOp) override final
 Apply the AMR operator, i.e. compute L(phi) in an AMR context. More...
 
void refluxFreeAMROperator (LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB >> *a_finerOp)
 Apply the AMR operator, i.e. compute L(phi) in an AMR context. More...
 
void AMRResidual (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB >> *a_finerOp) override final
 Compute residual on this level. AMR version. More...
 
void AMRResidualNF (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC) override final
 Compute AMR residual on finest AMR level. More...
 
void AMRResidualNC (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB >> *a_finerOp) override final
 Compute AMR residual on coarsest. More...
 
void AMROperatorNF (LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, bool a_homogeneousPhysBC) override final
 Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no finer levels. More...
 
void AMROperatorNC (LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoar, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< EBCellFAB >> *a_finerOp) override final
 Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no coarser AMR levels. More...
 
void AMRRestrict (LevelData< EBCellFAB > &a_residualCoarse, const LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection, bool a_skip_res) override final
 Restrict residual. More...
 
void AMRProlong (LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection) override final
 Prolongation onto AMR level. More...
 
void AMRUpdateResidual (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection) override final
 Update AMR residual. More...
 
void setAlphaAndBeta (const Real &a_alpha, const Real &a_beta) override final
 Set alpha coefficient and beta coefficient (can change as diffusion solvers progress) More...
 
void createCoarsened (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const int &a_refRat) override final
 Create coarsening of data holder. More...
 
void diagonalScale (LevelData< EBCellFAB > &a_rhs, bool a_kappaWeighted) override final
 Do diagonal scaling. More...
 
void divideByIdentityCoef (LevelData< EBCellFAB > &a_rhs) override final
 Divide by the a-coefficient. More...
 
void applyOpNoBoundary (LevelData< EBCellFAB > &a_ans, const LevelData< EBCellFAB > &a_phi) override final
 Apply operator but turn off all BCs.
 
void fillGrad (const LevelData< EBCellFAB > &a_phi) override final
 Not called, I think. More...
 
void getFlux (EBFluxFAB &a_flux, const LevelData< EBCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale) override final
 Fill flux. More...
 
void allocateFlux () const noexcept
 Allocate m_flux.
 
void deallocateFlux () const noexcept
 Deallocate m_flux.
 
LevelData< EBFluxFAB > & getFlux () const
 Returns m_flux. This is used in refluxing routines.
 

Protected Member Functions

void relaxPointJacobi (LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
 Jacobi relaxation. More...
 
void relaxGSRedBlack (LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
 Jacobi relaxation. More...
 
void relaxGSMultiColor (LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_residual, const int a_iterations)
 Multi-colored gauss-seidel relaxation. More...
 
void computeDiagWeight ()
 Calculate the weight of the diagonal term.
 
void computeRelaxationCoefficient ()
 Calculate relaxation coefficient.
 
void makeAggStencil ()
 Compute aggregated stencils.
 
void defineStencils ()
 Define stencils.
 
VoFStencil getFaceCenterFluxStencil (const FaceIndex &a_face, const DataIndex &a_dit) const
 Get the face-centered flux stencil. More...
 
VoFStencil getFaceCentroidFluxStencil (const FaceIndex &a_face, const DataIndex &a_dit) const
 Get the face-centroid flux stencil. More...
 
void computeFlux (EBFaceFAB &a_fluxCentroid, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit, const int a_dir)
 Compute centroid fluxes in a direction. More...
 
void computeFaceCenteredFlux (EBFaceFAB &a_fluxCenter, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit, const int a_dir)
 Compute face-centered fluxes. More...
 
void computeFaceCentroidFlux (EBFaceFAB &a_flux, const EBCellFAB &a_phi, const Box &a_cellBox, const DataIndex &a_dit, const int a_dir)
 Compute face-centered fluxes. More...
 
void computeFlux (const LevelData< EBCellFAB > &a_phi)
 Fill face centroid-centered fluxes on this level. More...
 
void reflux (LevelData< EBCellFAB > &a_Lphi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB >> &a_finerOp)
 Reflux algorithm. More...
 

Protected Attributes

Timer m_timer
 Timer so user can profile.
 
Interval m_interval
 Interval.
 
Smoother m_smoother
 Relaxation method.
 
Location::Cell m_dataLocation
 Data centering.
 
bool m_refluxFree
 Use reflux-free formulation or not.
 
bool m_profile
 Profile the operator.
 
bool m_hasMGObjects
 True if there is a multigrid level below this operator.
 
bool m_hasFine
 True if there's a finer level.
 
bool m_hasCoar
 True if there's a coarser level.
 
int m_refToCoar
 Refinement factor to coarse level.
 
int m_refToFine
 Refinement factor to fine level.
 
bool m_doInterpCF
 Do coarse-fine interpolation or not. More...
 
bool m_doExchange
 Turn on/off exchange operation. More...
 
bool m_doCoarsen
 Turn on/off exchange operation. More...
 
IntVect m_ghostPhi
 Ghost cells for phi.
 
IntVect m_ghostRhs
 Ghost cells for rhs (note, the operator rhs)
 
Real m_alpha
 Alpha-coefficient.
 
Real m_beta
 Beta-coefficient.
 
Real m_dx
 Grid resolution;.
 
RealVect m_vecDx
 Vector resolution.
 
RealVect m_probLo
 Lower-left corner of domain.
 
Copier m_exchangeCopier
 Pre-built exchange copier.
 
Copier m_exchangeCopierFine
 Pre-built exchange copier.
 
EBLevelGrid m_eblgFine
 Fine level grid (if the operator has a fine level)
 
EBLevelGrid m_eblg
 Grid.
 
EBLevelGrid m_eblgCoFi
 Coarsened of m_eblg.
 
EBLevelGrid m_eblgCoar
 Coarse level grid (if the operator has a coarse level)
 
EBLevelGrid m_eblgCoarMG
 Coarser grids (multigrid level)
 
EBMGRestrict m_restrictOp
 Restriction operator for AMR levels.
 
EBMGRestrict m_restrictOpMG
 Restriction operator if this is an MG level.
 
EBMGProlong m_prolongOp
 Prolongation operator for AMR levels.
 
EBMGProlong m_prolongOpMG
 Prolongation operator for MG levels.
 
RefCountedPtr< EBHelmholtzDomainBCm_domainBc
 Domain bc object.
 
RefCountedPtr< EBHelmholtzEBBCm_ebBc
 Domain bc object.
 
RefCountedPtr< EBMultigridInterpolatorm_interpolator
 Interpolator object.
 
RefCountedPtr< EBRefluxm_fluxReg
 Flux register.
 
RefCountedPtr< EBCoarAvem_coarAve
 Conservative coarsener.
 
RefCountedPtr< LevelData< EBCellFAB > > m_Acoef
 A-coefficient in Helmholtz equation.
 
RefCountedPtr< LevelData< EBFluxFAB > > m_Bcoef
 B-coefficient in Helmholtz equation.
 
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_BcoefIrreg
 B-coefficient in Helmholtz equation, but on EB faces.
 
LevelData< EBCellFAB > m_relCoef
 Relaxation coefficient.
 
LevelData< EBFluxFAB > * m_flux
 For holding fluxes.
 
LayoutData< BaseIFFAB< Real > > m_interpolant [SpaceDim]
 Interpolant for when we want centroid fluxes.
 
LayoutData< BaseIFFAB< FaceStencil > > m_interpStencil [SpaceDim]
 Face centroid interpolation stencil.
 
LayoutData< BaseIFFAB< VoFStencil > > m_centroidFluxStencil [SpaceDim]
 Face centroid flux stencil. Defined on all faces connecting one or more irregular vofs.
 
LayoutData< BaseIVFAB< VoFStencil > > m_relaxStencils
 Operator stencils in irregular cells (and ones that border irregular cells if using a centroid discretization). More...
 
LayoutData< RefCountedPtr< VCAggStencil > > m_aggRelaxStencil
 For making irregular stencil applications go faster. More...
 
LayoutData< BaseIVFAB< Real > > m_alphaDiagWeight
 Weights of diagonal alpha terms.
 
LayoutData< BaseIVFAB< Real > > m_betaDiagWeight
 Weights of diagonal beta terms.
 
LayoutData< VoFIterator > m_vofIterIrreg
 VoFIterator for irregular cells.
 
LayoutData< VoFIterator > m_vofIterMulti
 VoFIterator for "multi-cells".
 
LayoutData< VoFIterator > m_vofIterStenc
 VoFIterator which iterates over all cells that are 1) a cut-cell or 2) borders a cut-cell.
 
LayoutData< VoFIterator > m_vofIterDomLo [SpaceDim]
 VoF iterators for lo domain side.
 
LayoutData< VoFIterator > m_vofIterDomHi [SpaceDim]
 VoF iterators for hi domain side.
 
std::map< std::pair< int, Side::LoHiSide >, Box > m_sideBox
 Domain boxes on each side.
 
Vector< IntVect > m_colors
 "Colors" for the relaxation methods
 

Static Protected Attributes

static constexpr int m_nComp = 1
 Number of components that we solve for (always one..)
 
static constexpr int m_comp = 0
 Component that we solve for.
 

Detailed Description

Helmholtz operator for equations like alpha*a(x)*phi(x) + beta*div(b(x)*grad(phi(x))) = rho.

This can be used with TGA time stepping.

Constructor & Destructor Documentation

◆ EBHelmholtzOp() [1/3]

EBHelmholtzOp::EBHelmholtzOp ( const EBHelmholtzOp a_other)
delete

Disallowed copy constructor.

Parameters
[in]a_otherOther operator

◆ EBHelmholtzOp() [2/3]

EBHelmholtzOp::EBHelmholtzOp ( const EBHelmholtzOp &&  a_other)
delete

Disallowed move constructor.

Parameters
[in]a_otherOther operator

◆ EBHelmholtzOp() [3/3]

EBHelmholtzOp::EBHelmholtzOp ( const Location::Cell  a_dataLocation,
const EBLevelGrid &  a_eblgFine,
const EBLevelGrid &  a_eblg,
const EBLevelGrid &  a_eblgCoFi,
const EBLevelGrid &  a_eblgCoar,
const EBLevelGrid &  a_eblgCoarMG,
const RefCountedPtr< EBMultigridInterpolator > &  a_interpolator,
const RefCountedPtr< EBReflux > &  a_fluxReg,
const RefCountedPtr< EBCoarAve > &  a_coarAve,
const RefCountedPtr< EBHelmholtzDomainBC > &  a_domainBC,
const RefCountedPtr< EBHelmholtzEBBC > &  a_ebBC,
const RealVect &  a_probLo,
const Real &  a_dx,
const int &  a_refToFine,
const int &  a_refToCoar,
const bool &  a_hasFine,
const bool &  a_hasCoar,
const bool &  a_hasMGObjects,
const Real &  a_alpha,
const Real &  a_beta,
const RefCountedPtr< LevelData< EBCellFAB >> &  a_Acoef,
const RefCountedPtr< LevelData< EBFluxFAB >> &  a_Bcoef,
const RefCountedPtr< LevelData< BaseIVFAB< Real >>> &  a_BcoIrreg,
const IntVect &  a_ghostCellsPhi,
const IntVect &  a_ghostCellsRHS,
const Smoother a_smoother 
)

Full constructor.

Parameters
[in]a_dataLocationData location, either cell center or cell centroid
[in]a_eblgFineFine grids
[in]a_eblgGrids on this level
[in]a_eblgCoFiCoarsening of fine level grids
[in]a_eblgCoarCoarse grids
[in]a_eblgCoarMGMultigrid-grids
[in]a_interpolatorInterpolator
[in]a_fluxRegFlux register
[in]a_coarAveCoarsener
[in]a_domainBCDomain BC
[in]a_ebBCBoundary conditions on EBs
[in]a_probLoLower-left corner of computational domain
[in]a_dxGrid resolution
[in]a_refToFineRefinement ratio to fine level
[in]a_refToCoarRefinement ratio to coarse level
[in]a_hasFineHas fine level or not
[in]a_hasCoarHas coarse level or not
[in]a_hasMGObjectsHas multigrid-objects (special objects between AMR levels, or below the AMR levels)
[in]a_alphaOperator alpha
[in]a_betaOperator beta
[in]a_AcoefOperator A-coefficient
[in]a_BcoefOperator B-coefficient
[in]a_BcoefIrregOperator B-coefficient (on EB faces)
[in]a_ghostCellsPhiGhost cells in data holders
[in]a_ghostCellsPhiGhost cells in data holders
[in]a_smootherWhich smoother to use

Member Function Documentation

◆ AMROperator()

void EBHelmholtzOp::AMROperator ( LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoar,
const bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< EBCellFAB >> *  a_finerOp 
)
finaloverride

Apply the AMR operator, i.e. compute L(phi) in an AMR context.

Parameters
[out]a_residualResidual on this level
[in]a_phiFinePhi on fine level
[in]a_phiPhi on this level
[in]a_phiCoarPhi on coar level
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not
[in]a_finerOpFiner operatator

This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.

◆ AMROperatorNC()

void EBHelmholtzOp::AMROperatorNC ( LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoar,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< EBCellFAB >> *  a_finerOp 
)
finaloverride

Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no coarser AMR levels.

Parameters
[out]a_LphiL(phi)
[in]a_phiFinePhi on finer level
[in]a_phiPhi on this level
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not

This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.

◆ AMROperatorNF()

void EBHelmholtzOp::AMROperatorNF ( LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoar,
bool  a_homogeneousPhysBC 
)
finaloverride

Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no finer levels.

Parameters
[out]a_LphiL(phi)
[in]a_phiPhi on this level
[in]a_phiCoarPhi on coar level
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not

This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.

◆ AMRProlong()

void EBHelmholtzOp::AMRProlong ( LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_coarseCorrection 
)
finaloverride

Prolongation onto AMR level.

Parameters
[out]a_correctionInterpolated correction
[in]a_coarseCorrectionCorrection on coarse level

◆ AMRResidual()

void EBHelmholtzOp::AMRResidual ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoar,
const LevelData< EBCellFAB > &  a_rhs,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< EBCellFAB >> *  a_finerOp 
)
finaloverride

Compute residual on this level. AMR version.

Parameters
[out]a_residualResidual on this level
[in]a_phiFinePhi on fine level
[in]a_phiPhi on this level
[in]a_phiCoarPhi on coar level
[in]a_rhsRight-hand side on this level
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not
[in]a_finerOpFiner operatator

◆ AMRResidualNC()

void EBHelmholtzOp::AMRResidualNC ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_rhs,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< EBCellFAB >> *  a_finerOp 
)
finaloverride

Compute AMR residual on coarsest.

Parameters
[out]a_residualResidual on this level
[in]a_phiFinePhi on fine level
[in]a_phiPhi on this level
[in]a_rhsRight-hand side on this level
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not
[in]a_finerOpFiner operator

◆ AMRResidualNF()

void EBHelmholtzOp::AMRResidualNF ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoar,
const LevelData< EBCellFAB > &  a_rhs,
bool  a_homogeneousPhysBC 
)
finaloverride

Compute AMR residual on finest AMR level.

Parameters
[out]a_residualResidual on this level
[in]a_phiPhi on this level
[in]a_phiCoarPhi on coar level
[in]a_rhsRight-hand side on this level
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not

◆ AMRRestrict()

void EBHelmholtzOp::AMRRestrict ( LevelData< EBCellFAB > &  a_residualCoarse,
const LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_coarseCorrection,
bool  a_skip_res 
)
finaloverride

Restrict residual.

Parameters
[out]a_residualCoarseCoarse residual
[out]a_residualResidual
[out]a_correctionCorrection on this level
[out]a_coarseCorrectionCoarse level correction
[in]a_skip_resI have no idea what this one is supposed to do.

◆ AMRUpdateResidual()

void EBHelmholtzOp::AMRUpdateResidual ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_coarseCorrection 
)
finaloverride

Update AMR residual.

Parameters
[in]a_residualResidual
[in]a_correctionCorrection
[in]a_coarseCorrectionCoarse-level correction

◆ applyDomainFlux()

void EBHelmholtzOp::applyDomainFlux ( EBCellFAB &  a_phi,
const EBFluxFAB &  a_Bcoef,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const bool  a_homogeneousPhysBc 
) const
noexcept

Apply domain flux.

Parameters
[in,out]a_phiCell data
[in]a_BcoefHelmholtz B-coefficient
[in]a_cellBoxComputation box
[in]a_ditData index
[in]a_homogeneousPhysBCHomogeneous BC or not
Note
This sets the ghost cells of phi to phi(i-1) = phi(i) + F*dx/bco such that when we take the centered difference we get -bco*(phi(i) - phi(i-1))/dx = F.

◆ applyOp() [1/3]

void EBHelmholtzOp::applyOp ( EBCellFAB &  a_Lphi,
EBCellFAB &  a_phi,
const EBCellFAB &  a_Acoef,
const EBFluxFAB &  a_Bcoef,
const BaseIVFAB< Real > &  a_BcoefIrreg,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const bool  a_homogeneousPhysBC 
) const
noexcept

Apply operator in a grid box.

Parameters
[out]a_LphiL(phi)
[in]a_phiPhi
[in]a_cellBoxGrid box
[in]a_ditData indxe
[in]a_homogeneousPhysBCHomogeneous physical BCs or not

◆ applyOp() [2/3]

void EBHelmholtzOp::applyOp ( LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_phi,
bool  a_homogeneousPhysBc 
)
finaloverride

Apply operator.

Parameters
[out]a_LphiL(phi)
[in]a_phiPhi
[in]a_homogeneousPhysBcHomogeneous physical BCs or not

This computes a_Lphi = L(a_phi) using homogeneous physical BCs or not

◆ applyOp() [3/3]

void EBHelmholtzOp::applyOp ( LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > *const  a_phiCoar,
const bool  a_homogeneousPhysBC,
const bool  a_homogeneousCFBC 
)

Apply operator on this level. This is a more general version which can turn on/off homogeneous and CF bcs.

Parameters
[out]a_LphiL(phi)
[out]a_phiPhi on this level
[out]a_phiCoarCoarse-level phi. If you have a coar this
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not
[in]a_homogeneousCFBCUse homogeneous coarse-fine bcs or not
Note
a_phi is not really const and we do a nasty cast. Leaving it as const because we only modify the ghost cells!

◆ applyOpIrregular()

void EBHelmholtzOp::applyOpIrregular ( EBCellFAB &  a_Lphi,
const EBCellFAB &  a_phi,
const EBCellFAB &  a_Acoef,
const EBFluxFAB &  a_Bcoef,
const BaseIVFAB< Real > &  a_BcoefIrreg,
const BaseIVFAB< Real > &  a_alphaDiagWeight,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const bool  a_homogeneousPhysBC 
) const
noexcept

Apply operator in irregular cells.

Parameters
[out]a_LphiL(phi)
[in]a_phiCell-centered data
[in]a_cellBoxCell box
[in]a_ditGrid index
[in]a_homogeneousPhysBCHomogeneous physical BCs or not

◆ applyOpRegular()

void EBHelmholtzOp::applyOpRegular ( EBCellFAB &  a_Lphi,
EBCellFAB &  a_phi,
const EBCellFAB &  a_Acoef,
const EBFluxFAB &  a_Bcoef,
const BaseIVFAB< Real > &  a_BcoefIrreg,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const bool  a_homogeneousPhysBC 
) const
noexcept

Apply operator in regular cells.

Parameters
[out]a_LphiL(phi)
[in]a_phiPhi
[in]a_cellBoxGrid box
[in]a_ditData indxe
[in]a_homogeneousPhysBCHomogeneous physical BCs or not

◆ assign()

void EBHelmholtzOp::assign ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs 
)
finaloverride

Assign data.

This does a copy from rhs to lhs

◆ assignCopier()

void EBHelmholtzOp::assignCopier ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs,
const Copier &  a_copier 
)
finaloverride

Assign lhs.

This is the version that is called by AMRMultiGrid::VCycle. Note that the other version might be called by other operators (e.g., EBBackwardEuler)

Parameters
[out]a_lhsOutgoing data
[in]a_rhsIncoming data
[in]a_copierCopier

◆ assignLocal()

void EBHelmholtzOp::assignLocal ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs 
)
finaloverride

Local assignment function.

Parameters
[out]a_lhs.Equal to a_rhs on output
[in]a_rhs.Data

◆ axby()

void EBHelmholtzOp::axby ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_x,
const LevelData< EBCellFAB > &  a_y,
const Real  a_a,
const Real  a_b 
)
finaloverride

Set a_lhs = a*x + b*y.

Parameters
[out]a_lhsResult data
[in]a_xx-data
[in]a_yy-data
[in]a_aScaling factor
[in]a_bScaling factor

◆ buildCopier()

void EBHelmholtzOp::buildCopier ( Copier &  a_copier,
const LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs 
)
override

Build copier.

Parameters
[out]a_copierCopier for copying between a_lhs and a_rhs
[in]a_lhsCopying from
[in]a_rhsCopying to

◆ coarsenCell()

void EBHelmholtzOp::coarsenCell ( LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiFine 
)

Coarsen data from fine to coar level.

Parameters
[in]a_phiData on this level
[in]a_phiFineData on finer level
Note
This routine is used in AMROperator to ensure that EB fluxes which are described by stencils get consistent data under the fine level.

◆ coarsenFlux()

void EBHelmholtzOp::coarsenFlux ( LevelData< EBFluxFAB > &  a_flux,
const LevelData< EBFluxFAB > &  a_fineFlux 
)

Coarsen fluxes on the fine level onto this level.

User must update fluxes before calling this routine.

Parameters
[in]a_fluxFlux on this level.
[in]a_fineFluxFlux on finer level
Note
This routine is used in the reflux-free AMROperator

◆ computeFaceCenteredFlux()

void EBHelmholtzOp::computeFaceCenteredFlux ( EBFaceFAB &  a_fluxCenter,
const EBCellFAB &  a_phi,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const int  a_dir 
)
protected

Compute face-centered fluxes.

Parameters
[out]a_fluxCenterFlux on center
[ina_phi Cell-centered data
[in]a_cellBoxCell-centered compute box.
[in]a_ditData index
Note
This computes the flux for all dir-oriented faces in a_cellBox except for boundary faces.

◆ computeFaceCentroidFlux()

void EBHelmholtzOp::computeFaceCentroidFlux ( EBFaceFAB &  a_flux,
const EBCellFAB &  a_phi,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const int  a_dir 
)
protected

Compute face-centered fluxes.

Parameters
[out]a_fluxCenterFlux on center
[ina_phi Cell-centered data
[in]a_cellBoxCell-centered compute box.
[in]a_ditData index
Note
This computes the flux for all dir-oriented faces in a_cellBox except for boundary faces.

◆ computeFlux() [1/2]

void EBHelmholtzOp::computeFlux ( const LevelData< EBCellFAB > &  a_phi)
protected

Fill face centroid-centered fluxes on this level.

Parameters
[in]a_phiCell-centered data. Must have updated ghost cells for this to make any sense.
Note
This computes the flux onto m_flux everywhere except on boundary faces.

◆ computeFlux() [2/2]

void EBHelmholtzOp::computeFlux ( EBFaceFAB &  a_fluxCentroid,
const EBCellFAB &  a_phi,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const int  a_dir 
)
protected

Compute centroid fluxes in a direction.

Parameters
[out]a_fluxCentroidFlux on face centroid
[ina_phi Cell-centered data
[in]a_cellBoxCell-centered compute box.
[in]a_ditData index
Note
This computes the flux for all dir-oriented faces in a_cellBox except for boundary faces.

◆ create()

void EBHelmholtzOp::create ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs 
)
finaloverride

Create data which clones the layout of the other.

Parameters
[out]a_lhsData clone (returned data is not initialized)
[out]a_rhsData layout to be cloned

◆ createCoarsened()

void EBHelmholtzOp::createCoarsened ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs,
const int &  a_refRat 
)
finaloverride

Create coarsening of data holder.

Parameters
[out]a_lhsCoarsened data
[in]a_rhsFine data
[in]a_refRatCoarsening factor

◆ createCoarser()

void EBHelmholtzOp::createCoarser ( LevelData< EBCellFAB > &  a_coarse,
const LevelData< EBCellFAB > &  a_fine,
bool  a_ghosted 
)
finaloverride

Create coarsened data.

Parameters
[out]a_coarseCoarse data
[in]a_fineFine data
[in]a_ghostedInclude ghost cells or nto

◆ diagonalScale()

void EBHelmholtzOp::diagonalScale ( LevelData< EBCellFAB > &  a_rhs,
bool  a_kappaWeighted 
)
finaloverride

Do diagonal scaling.

Parameters
[in,out]a_rhsData to be scaled
[in]a_kappaWeightedUse kappa-weighted scaling or not

◆ divideByIdentityCoef()

void EBHelmholtzOp::divideByIdentityCoef ( LevelData< EBCellFAB > &  a_rhs)
finaloverride

Divide by the a-coefficient.

Parameters
[in,out]a_rhsDivided data

◆ fillDomainFlux()

void EBHelmholtzOp::fillDomainFlux ( EBFluxFAB &  a_flux,
const EBCellFAB &  a_phi,
const Box &  a_cellBox,
const DataIndex &  a_dit 
)

Fill domain flux. This fills the flux on the domain face using centered differencing ala applyDomainFlux.

a_flux is replaced by the user-specified flux on the domain faces.

Parameters
[in,out]a_fluxFlux data holder.
[in,out]a_phiCell-centered data
[in]a_cellBoxComputation box
[in]a_ditData index

◆ fillGrad()

void EBHelmholtzOp::fillGrad ( const LevelData< EBCellFAB > &  a_phi)
finaloverride

Not called, I think.

Parameters
[in]a_phiPhi

◆ gauSaiMultiColor()

void EBHelmholtzOp::gauSaiMultiColor ( LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_rhs,
const IntVect &  a_color 
) const

Multi-colored Gauss-Seidel kernel. Public because MFHelmholtzOp may want to use use.

Parameters
[in]a_phiCorrection
[in]a_phiLphi L(corrcetion)
[in]a_rhsResidual
[in]a_color"Color":

◆ gauSaiMultiColorKernel()

void EBHelmholtzOp::gauSaiMultiColorKernel ( EBCellFAB &  a_Lcorr,
EBCellFAB &  a_corr,
const EBCellFAB &  a_resid,
const EBCellFAB &  a_Acoef,
const EBFluxFAB &  a_Bcoef,
const BaseIVFAB< Real > &  a_BcoefIrreg,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const IntVect &  a_color 
) const
noexcept

Multi-color Gauss-Seidel kernel.

Parameters
[in,out]a_LcorrStorage for computing L(a_corr)
[in,out]a_corrCorrection
[in]a_residResidual
[in]a_AcoefA-coefficient
[in]a_BcoefB-coefficient
[in]a_BcoefIrregB-coefficient on EB faces
[in]a_cellBoxGrid box
[in]a_ditData index
[in]a_colorColor

◆ gauSaiRedBlackKernel()

void EBHelmholtzOp::gauSaiRedBlackKernel ( EBCellFAB &  a_Lcorr,
EBCellFAB &  a_corr,
const EBCellFAB &  a_resid,
const EBCellFAB &  a_Acoef,
const EBFluxFAB &  a_Bcoef,
const BaseIVFAB< Real > &  a_BcoefIrreg,
const Box &  a_cellBox,
const DataIndex &  a_dit,
const int &  a_redBlack 
) const
noexcept

Red-black Gauss-Seidel kernel.

Parameters
[in,out]a_LcorrStorage for computing L(a_corr)
[in,out]a_corrCorrection
[in]a_residResidual
[in]a_AcoefA-coefficient
[in]a_BcoefB-coefficient
[in]a_BcoefIrregB-coefficient on EB faces
[in]a_redBlackRed or black
[in]a_cellBoxGrid box
[in]a_ditData index

◆ getAcoef()

const RefCountedPtr< LevelData< EBCellFAB > > & EBHelmholtzOp::getAcoef ( )

Get the Helmholtz A-coefficient on cell centers.

Returns
m_Acoef

◆ getBcoef()

const RefCountedPtr< LevelData< EBFluxFAB > > & EBHelmholtzOp::getBcoef ( )

Get the Helmholtz B-coefficient on faces.

Returns
m_Bcoef

◆ getBcoefIrreg()

const RefCountedPtr< LevelData< BaseIVFAB< Real > > > & EBHelmholtzOp::getBcoefIrreg ( )

Get the Helmholtz B-coefficient on the EB.

Returns
m_Bcoef

◆ getFaceCenterFluxStencil()

VoFStencil EBHelmholtzOp::getFaceCenterFluxStencil ( const FaceIndex &  a_face,
const DataIndex &  a_dit 
) const
protected

Get the face-centered flux stencil.

Parameters
[in]a_faceFace
[in]a_ditData index. Need because we multiply by B-coefficient
Returns
Returns a stencil for the face-centered flux for the given face
Note
If the face is a boundary face the returned stencil is empty.

◆ getFaceCentroidFluxStencil()

VoFStencil EBHelmholtzOp::getFaceCentroidFluxStencil ( const FaceIndex &  a_face,
const DataIndex &  a_dit 
) const
protected

Get the face-centroid flux stencil.

Parameters
[in]a_faceFace
[in]a_ditData index. Need because we multiply by B-coefficient
Returns
Returns a stencil for the face-centroid flux for the given face
Note
The current version just interpolates the centered fluxes.

◆ getFlux()

void EBHelmholtzOp::getFlux ( EBFluxFAB &  a_flux,
const LevelData< EBCellFAB > &  a_data,
const Box &  a_grid,
const DataIndex &  a_dit,
Real  a_scale 
)
finaloverride

Fill flux.

Parameters
[out]a_fluxFlux
[in]a_dataData for which we will compute the flux
[in]a_gridGrid
[in]a_ditCorresponding data index.
[in]a_scaleScaling factor

◆ homogeneousCFInterp()

void EBHelmholtzOp::homogeneousCFInterp ( LevelData< EBCellFAB > &  a_phi)

Do homogeneous coarse-fine interpolation.

Parameters
[in]a_phiData

◆ incr()

void EBHelmholtzOp::incr ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs,
const Real  a_scale 
)
finaloverride

Increment operator.

Parameters
[in,out]a_lhsData to be incremented
[in]a_rhsIncrementation data
[in]a_scaleScaling factor

◆ inhomogeneousCFInterp()

void EBHelmholtzOp::inhomogeneousCFInterp ( LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoar 
)

Inhomogeneous coarse-fine interpolation.

Parameters
[in,out]a_phiFineFine data. Ghost cells will be filled.
[in,out]a_phiCoarCoarse data.

◆ interpolateCF()

void EBHelmholtzOp::interpolateCF ( LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > *  a_phiCoar,
const bool  a_homogeneousCFBC 
)

Apply coarse-fine boundary conditions.

Parameters
[in,out]a_phiFineFine data
[in,out]a_phiCoarCoarse data
[in]a_homogeneousCFBCHomogeneous CF or not (i.e. coar data is zero);

◆ norm()

Real EBHelmholtzOp::norm ( const LevelData< EBCellFAB > &  a_rhs,
const int  a_order 
)
finaloverride

Compute norm of data.

Parameters
[in]

◆ operator=() [1/2]

EBHelmholtzOp& EBHelmholtzOp::operator= ( const EBHelmholtzOp &&  a_oper)
delete

No move assigment allowed.

Parameters
[in]a_otherOther operator

◆ operator=() [2/2]

EBHelmholtzOp& EBHelmholtzOp::operator= ( const EBHelmholtzOp a_oper)
delete

No copy assigment allowed.

Parameters
[in]a_otherOther operator

◆ pointJacobiKernel()

void EBHelmholtzOp::pointJacobiKernel ( EBCellFAB &  a_Lcorr,
EBCellFAB &  a_corr,
const EBCellFAB &  a_resid,
const EBCellFAB &  a_Acoef,
const EBFluxFAB &  a_Bcoef,
const BaseIVFAB< Real > &  a_BcoefIrreg,
const Box &  a_cellBox,
const DataIndex &  a_dit 
) const
noexcept

Point Jacobi kernel.

Parameters
[in,out]a_LcorrStorage for computing L(a_corr)
[in,out]a_corrCorrection
[in]a_residResidual
[in]a_AcoefA-coefficient
[in]a_BcoefB-coefficient
[in]a_BcoefIrregB-coefficient on EB faces
[in]a_cellBoxGrid box
[in]a_ditData index

◆ preCond()

void EBHelmholtzOp::preCond ( LevelData< EBCellFAB > &  a_corr,
const LevelData< EBCellFAB > &  a_residual 
)
finaloverride

Precondition system before bottom solve.

Parameters
[in]a_corrCorrection
[in]a_residualResidual

This just runs a few relaxations.

◆ prolongIncrement()

void EBHelmholtzOp::prolongIncrement ( LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_correctCoarse 
)
finaloverride

Prolongation method.

Parameters
[out]a_phiCorrection on this level
[out]a_correctCoarseCorrection on coarse level

◆ reflux()

void EBHelmholtzOp::reflux ( LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
AMRLevelOp< LevelData< EBCellFAB >> &  a_finerOp 
)
protected

Reflux algorithm.

Parameters
[in,out]a_LphiL(phiFine, phiCoar). This is modified along the CF interface to consistently account for fine-level fluxes.
[in]a_phiFineData on finer level
[in]a_phiData on this operator's level
[in]a_finerOpOperator on fine level

◆ refluxFreeAMROperator()

void EBHelmholtzOp::refluxFreeAMROperator ( LevelData< EBCellFAB > &  a_Lphi,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoar,
const bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< EBCellFAB >> *  a_finerOp 
)

Apply the AMR operator, i.e. compute L(phi) in an AMR context.

This is the reflux-free version of AMROperator.

Parameters
[out]a_residualResidual on this level
[in]a_phiFinePhi on fine level
[in]a_phiPhi on this level
[in]a_phiCoarPhi on coar level
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not
[in]a_finerOpFiner operatator

This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.

◆ relax()

void EBHelmholtzOp::relax ( LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_residual,
int  a_iterations 
)
finaloverride

Relaxation method. This does smoothing for the system L(correction) = residual.

Parameters
[in,out]a_correctionCorrection
[in]a_residualResidual
[in]a_iterationsNumber of iterations

◆ relaxGSMultiColor()

void EBHelmholtzOp::relaxGSMultiColor ( LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_residual,
const int  a_iterations 
)
protected

Multi-colored gauss-seidel relaxation.

Parameters
[in,out]a_correctionCorrection
[in]a_residualResidual
[in]a_iterationsNumber of iterations

◆ relaxGSRedBlack()

void EBHelmholtzOp::relaxGSRedBlack ( LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_residual,
const int  a_iterations 
)
protected

Jacobi relaxation.

Parameters
[in,out]a_correctionCorrection
[in]a_residualResidual
[in]a_iterationsNumber of iterations

◆ relaxPointJacobi()

void EBHelmholtzOp::relaxPointJacobi ( LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_residual,
const int  a_iterations 
)
protected

Jacobi relaxation.

Parameters
[in,out]a_correctionCorrection
[in]a_residualResidual
[in]a_iterationsNumber of iterations

◆ residual()

void EBHelmholtzOp::residual ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_rhs,
const bool  a_homogeneousPhysBc 
)
override

Compute residual on this level.

Parameters
[out]a_residualResidual rhs - L(phi)
[in]a_phiphi
[in]a_rhsRight-hand side of system
[in]a_homogeneousPhysBCUse homogeneous physical BCs or not

◆ restrictResidual()

void EBHelmholtzOp::restrictResidual ( LevelData< EBCellFAB > &  a_resCoar,
LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_rhs 
)
finaloverride

Restrict residual onto coarse level.

Parameters
[in,out]a_resCoarCoarse residual
[in,out]a_phiPhi on this level
[in]a_rhsRhs on this level

◆ scale()

void EBHelmholtzOp::scale ( LevelData< EBCellFAB > &  a_lhs,
const Real &  a_scale 
)
finaloverride

Scale data. Returns a_lhs = a_lhs*a_scale.

Parameters
[in,out]a_lhsData to be scaled
[i]a_scale Scaling factor

◆ setAcoAndBco()

void EBHelmholtzOp::setAcoAndBco ( const RefCountedPtr< LevelData< EBCellFAB >> &  a_Acoef,
const RefCountedPtr< LevelData< EBFluxFAB >> &  a_Bcoef,
const RefCountedPtr< LevelData< BaseIVFAB< Real >>> &  a_BcoefIrreg 
)

Update with new A and B coefficients.

Note
This will call defineStencils but will not recompute the stencils from the BC objects.
Parameters
[in]a_AcoefOperator A-coefficient
[in]a_BcoefOperator B-coefficient
[in]a_BcoefIrregOperator B-coefficient (on EB faces)

◆ setAlphaAndBeta()

void EBHelmholtzOp::setAlphaAndBeta ( const Real &  a_alpha,
const Real &  a_beta 
)
finaloverride

Set alpha coefficient and beta coefficient (can change as diffusion solvers progress)

Parameters
[in]a_alphaAlpha-coefficient
[in]a_betaBeta-coefficient

◆ setToZero()

void EBHelmholtzOp::setToZero ( LevelData< EBCellFAB > &  a_lhs)
finaloverride

Set data to zero.

Parameters
[in,out]a_lhsData to be set to zero.

◆ turnOffCoarsening()

void EBHelmholtzOp::turnOffCoarsening ( )

Turn off coarsening operation.

Note
Done when MFHelmholtzOp does the coarsening and we don't want EBHelmholtzOp to redo it.

◆ turnOffExchange()

void EBHelmholtzOp::turnOffExchange ( )

Turn off exchange operation.

Note
Done when MFHelmholtzOp does the exchange and we don't want EBHelmholtzOp to redo it.

◆ turnOnCoarsening()

void EBHelmholtzOp::turnOnCoarsening ( )

Turn on coarsening operation.

Note
Done when MFHelmholtzOp does the coarsening and we don't want EBHelmholtzOp to redo it.

◆ turnOnExchange()

void EBHelmholtzOp::turnOnExchange ( )

Turn on exchange operation.

Note
Done when MFHelmholtzOp does the exchange and we don't want EBHelmholtzOp to redo it.

Member Data Documentation

◆ m_aggRelaxStencil

LayoutData<RefCountedPtr<VCAggStencil> > EBHelmholtzOp::m_aggRelaxStencil
protected

For making irregular stencil applications go faster.

This wraps m_relaxStencils in VCAggStencil (which computes explicit stencil offsets)

◆ m_doCoarsen

bool EBHelmholtzOp::m_doCoarsen
protected

Turn on/off exchange operation.

Note
Used by MFHelmholtzOp for optimizing how calls in EBHelmholtzOp are made

◆ m_doExchange

bool EBHelmholtzOp::m_doExchange
protected

Turn on/off exchange operation.

Note
Used by MFHelmholtzOp for optimizing how calls in EBHelmholtzOp are made

◆ m_doInterpCF

bool EBHelmholtzOp::m_doInterpCF
protected

Do coarse-fine interpolation or not.

Note
Used by MFHelmholtzOp for optimizing how calls in EBHelmholtzOp are made

◆ m_relaxStencils

LayoutData<BaseIVFAB<VoFStencil> > EBHelmholtzOp::m_relaxStencils
protected

Operator stencils in irregular cells (and ones that border irregular cells if using a centroid discretization).

This stencil is => sum(fluxes)/dx, not including boundary faces or EB faces. I.e. this is the same as kappa*div(F) with the exclusion of faces where we have boundary conditions.


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