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

Operator for solving multifluid Helmholtz on a grid level. More...

#include <CD_MFHelmholtzOp.H>

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

Public Types

enum class  Smoother { PointJacobi , GauSaiRedBlack , GauSaiMultiColor }
 Relaxation methods for this operator.
 

Public Member Functions

 MFHelmholtzOp ()=delete
 Constructor. Must subsequently call define(...)
 
 MFHelmholtzOp (const MFHelmholtzOp &a_op)=delete
 No copy construction allowed.
 
 MFHelmholtzOp (const MFHelmholtzOp &&a_op)=delete
 No move construction allowed.
 
 MFHelmholtzOp (const Location::Cell a_dataLocation, const MFLevelGrid &a_mflgFine, const MFLevelGrid &a_mflg, const MFLevelGrid &a_mflgCoFi, const MFLevelGrid &a_mflgCoar, const MFLevelGrid &a_mflgCoarMG, const MFMultigridInterpolator &a_interpolator, const MFReflux &a_fluxReg, const MFCoarAve &a_coarAve, const RefCountedPtr< MFHelmholtzDomainBCFactory > &a_domainBC, const RefCountedPtr< MFHelmholtzEBBCFactory > &a_ebBC, const RefCountedPtr< MFHelmholtzJumpBCFactory > &a_jumpBcFactory, 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 bool &a_isMGOperator, const Real &a_alpha, const Real &a_beta, const RefCountedPtr< LevelData< MFCellFAB >> &a_Acoef, const RefCountedPtr< LevelData< MFFluxFAB >> &a_Bcoef, const RefCountedPtr< LevelData< MFBaseIVFAB >> &a_BcoefIrreg, const IntVect &a_ghostPhi, const IntVect &a_ghostRhs, const int &a_jumpOrder, const int &a_jumpWeight, const Smoother &a_smoother)
 Full constructor. More...
 
 ~MFHelmholtzOp ()
 Destructor.
 
MFHelmholtzOpoperator= (const MFHelmholtzOp &a_oper)=delete
 No copy assigment allowed.
 
MFHelmholtzOpoperator= (const MFHelmholtzOp &&a_oper)=delete
 No move assigment allowed.
 
void setAcoAndBco (const RefCountedPtr< LevelData< MFCellFAB >> &a_Acoef, const RefCountedPtr< LevelData< MFFluxFAB >> &a_Bcoef, const RefCountedPtr< LevelData< MFBaseIVFAB >> &a_BcoefIrreg)
 Update operators with new coefficients. More...
 
const RefCountedPtr< LevelData< MFCellFAB > > & getAcoef ()
 Get the Helmholtz A-coefficient on cell centers. More...
 
const RefCountedPtr< LevelData< MFFluxFAB > > & getBcoef ()
 Get the Helmholtz B-coefficient on faces. More...
 
const RefCountedPtr< LevelData< MFBaseIVFAB > > & getBcoefIrreg ()
 Get the Helmholtz B-coefficient on the EB. More...
 
void setJump (RefCountedPtr< LevelData< BaseIVFAB< Real >>> &a_jump)
 Set the jump boundary condition.
 
RefCountedPtr< EBHelmholtzOp > & getHelmOp (const int a_phase)
 Get Helmholtz operator.
 
int refToCoarser () override final
 Return coarsening factor to coarser level (1 if there is no coarser level);.
 
void setAlphaAndBeta (const Real &a_alpha, const Real &a_beta) override final
 Set alpha and beta. More...
 
void divideByIdentityCoef (LevelData< MFCellFAB > &a_rhs) override final
 Divide by the a-coefficient. More...
 
void applyOpNoBoundary (LevelData< MFCellFAB > &a_ans, const LevelData< MFCellFAB > &a_phi) override final
 Apply operator but turn off all BCs.
 
void fillGrad (const LevelData< MFCellFAB > &a_phi) override final
 Not called, I think. More...
 
void getFlux (MFFluxFAB &a_flux, const LevelData< MFCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale) override final
 Fill flux. More...
 
void residual (LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_rhs, const bool a_homogeneousPhysBc) override final
 Compute residual on this level. More...
 
void preCond (LevelData< MFCellFAB > &a_corr, const LevelData< MFCellFAB > &a_residual) override final
 Precondition system before bottom solve. More...
 
void applyOp (LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phi, bool a_homogeneousPhysBc) override final
 Apply operator. More...
 
template<typename Duration = std::chrono::microseconds>
Vector< long long > computeOperatorLoads (LevelData< MFCellFAB > &a_phi, const int a_numApply) noexcept
 Time the applyOp routine. Template parameter is std::chrono duration. E.g. std::chrono::microseconds. More...
 
void relax (LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, int a_iterations) override final
 Relaxation method. This does smoothing for the system L(correction) = residual. More...
 
void relaxPointJacobi (LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, const int a_iterations)
 Jacobi relaxation. More...
 
void relaxGSRedBlack (LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, const int a_iterations)
 Jacobi relaxation. More...
 
void relaxGSMultiColor (LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_residual, const int a_iterations)
 Multi-colored gauss-seidel relaxation. More...
 
void create (LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override final
 Create method. More...
 
void createCoarser (LevelData< MFCellFAB > &a_coarse, const LevelData< MFCellFAB > &a_fine, bool a_ghosted) override final
 Create coarsened data. More...
 
void createCoarsened (LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs, const int &a_refRat) override final
 Create coarsening of data holder. More...
 
void incr (LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs, Real a_scale) override final
 Increment function. More...
 
Real dotProduct (const LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_2) override final
 Dot product.
 
void scale (LevelData< MFCellFAB > &a_lhs, const Real &a_scale) override final
 Scale function. More...
 
void setToZero (LevelData< MFCellFAB > &a_lhs) override final
 Set to zero. More...
 
void assign (LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override final
 Assignment fucntion. More...
 
void assignCopier (LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs, const Copier &a_copier) override final
 Assign lhs. More...
 
void assignLocal (LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override final
 Local assignment function. More...
 
void buildCopier (Copier &a_copier, const LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_rhs) override
 Build copier. More...
 
Real norm (const LevelData< MFCellFAB > &a_lhs, int a_order) override final
 Compute solution norm. More...
 
void axby (LevelData< MFCellFAB > &a_lhs, const LevelData< MFCellFAB > &a_x, const LevelData< MFCellFAB > &a_y, const Real a_a, const Real a_b) override final
 Set a_lhs = a*x + b*y. More...
 
void restrictResidual (LevelData< MFCellFAB > &a_resCoar, LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_rhs) override final
 Restrict residual onto coarse level. More...
 
void prolongIncrement (LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_correctCoarse) override final
 Prolongation method. More...
 
void AMRUpdateResidual (LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_coarseCorrection) override final
 Update AMR residual. More...
 
void AMRRestrict (LevelData< MFCellFAB > &a_residualCoarse, const LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_coarseCorrection, bool a_skip_res) override final
 Restrict residual. More...
 
void AMRProlong (LevelData< MFCellFAB > &a_correction, const LevelData< MFCellFAB > &a_coarseCorrection) override final
 Prolongation onto AMR level. More...
 
void AMRResidual (LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phiFine, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, const LevelData< MFCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override final
 Compute residual on this level. AMR version. More...
 
void AMRResidualNF (LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, const LevelData< MFCellFAB > &a_rhs, bool a_homogeneousPhysBC) override final
 Compute AMR residual on finest AMR level. More...
 
void AMRResidualNC (LevelData< MFCellFAB > &a_residual, const LevelData< MFCellFAB > &a_phiFine, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override final
 Compute AMR residual on coarsest. More...
 
void AMROperatorNF (LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &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< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override final
 Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no coarser AMR levels. More...
 

Protected Member Functions

void updateJumpBC (const LevelData< MFCellFAB > &a_phi, const bool a_homogeneousPhysBC)
 Update the jump condition. More...
 
void exchangeGhost (const LevelData< MFCellFAB > &a_phi) const
 Perform an exchange operation, event if the data is const. More...
 
void interpolateCF (const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > *a_phiCoar, const bool a_homogeneousCF)
 Do coarse-fine interpolation. More...
 
void applyOp (LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > *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 AMROperator (LevelData< MFCellFAB > &a_Lphi, const LevelData< MFCellFAB > &a_phiFine, const LevelData< MFCellFAB > &a_phi, const LevelData< MFCellFAB > &a_phiCoar, const bool a_homogeneousPhysBC, AMRLevelOp< LevelData< MFCellFAB >> *a_finerOp) override
 Apply the AMR operator, i.e. compute L(phi) in an AMR context. More...
 

Protected Attributes

Location::Cell m_dataLocation
 Interpretation of data. Either on cell center or on cell centroid.
 
Smoother m_smoother
 Relaxation method.
 
Vector< IntVect > m_colors
 "Colors" for the multi-coloered relaxation method
 
std::map< int, RefCountedPtr< EBHelmholtzOp > > m_helmOps
 Helmholtz operators on each phase. Note that I'm using int rather than Phase as identifier because that is the standard terminology for MFCellFAB.
 
RefCountedPtr< MFHelmholtzJumpBCm_jumpBC
 BC jump object. This is the one that has the stencils and can compute derivatives.
 
RefCountedPtr< LevelData< BaseIVFAB< Real > > > m_jump
 Actual BC jump in data-based format. This is the right-hand side of dphi/dn1 + dphi/dn2 = jump.
 
std::map< int, RefCountedPtr< LevelData< BaseIVFAB< Real > > > > m_dirichletBcValues
 Dirichlet BC values for each phase.
 
MFLevelGrid m_mflg
 Level grid.
 
MFLevelGrid m_mflgFine
 Fine level grid.
 
MFLevelGrid m_mflgCoar
 Coarse grid.
 
MFLevelGrid m_mflgCoFi
 Coarsened version of this grid.
 
MFLevelGrid m_mflgCoarMG
 Coarse multigrid-grid.
 
MFMultigridInterpolator m_interpolator
 Multi grid interpolator.
 
MFCoarAve m_coarAve
 Coarsener.
 
RefCountedPtr< LevelData< MFCellFAB > > m_Acoef
 Helmholtz A-coefficeint.
 
RefCountedPtr< LevelData< MFFluxFAB > > m_Bcoef
 Helmholtz B-coefficeint.
 
RefCountedPtr< LevelData< MFBaseIVFAB > > m_BcoefIrreg
 Helmholtz B-coefficeint on EB.
 
Copier m_exchangeCopier
 Copier for exchange operation.
 
bool m_multifluid
 Multifluid operator or not.
 
bool m_hasMGObjects
 Has MG objects or not.
 
int m_numPhases
 Number of phases.
 
int m_refToCoar
 Refinement factor to coarser AMR level.
 
IntVect m_ghostPhi
 Number of ghost cells.
 
IntVect m_ghostRhs
 Number of ghost cells.
 
bool m_hasCoar
 True if there is a coarser AMR level.
 
bool m_hasFine
 True if there is a finer AMR level.
 
bool m_hasMGObjcts
 True if there are multigrid levels.
 

Static Protected Attributes

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

Detailed Description

Operator for solving multifluid Helmholtz on a grid level.

This operator should typically not be constructed directly by the user but used together with Chombo's AMRMultiGrid. In this case the user will want to constructor these operators through the factor.

Constructor & Destructor Documentation

◆ MFHelmholtzOp()

MFHelmholtzOp::MFHelmholtzOp ( const Location::Cell  a_dataLocation,
const MFLevelGrid a_mflgFine,
const MFLevelGrid a_mflg,
const MFLevelGrid a_mflgCoFi,
const MFLevelGrid a_mflgCoar,
const MFLevelGrid a_mflgCoarMG,
const MFMultigridInterpolator a_interpolator,
const MFReflux a_fluxReg,
const MFCoarAve a_coarAve,
const RefCountedPtr< MFHelmholtzDomainBCFactory > &  a_domainBC,
const RefCountedPtr< MFHelmholtzEBBCFactory > &  a_ebBC,
const RefCountedPtr< MFHelmholtzJumpBCFactory > &  a_jumpBcFactory,
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 bool &  a_isMGOperator,
const Real &  a_alpha,
const Real &  a_beta,
const RefCountedPtr< LevelData< MFCellFAB >> &  a_Acoef,
const RefCountedPtr< LevelData< MFFluxFAB >> &  a_Bcoef,
const RefCountedPtr< LevelData< MFBaseIVFAB >> &  a_BcoefIrreg,
const IntVect &  a_ghostPhi,
const IntVect &  a_ghostRhs,
const int &  a_jumpOrder,
const int &  a_jumpWeight,
const Smoother a_smoother 
)

Full constructor.

Parameters
[in]a_dataLocationData location, either cell center or cell centroid
[in]a_mflgFineFine grids
[in]a_mflgGrids on this level
[in]a_mflgCoFiCoarsening of fine level grids
[in]a_mflgCoarCoarse grids
[in]a_mflgCoarMGMultigrid-grids
[in]a_interpolatorInterpolator
[in]a_fluxRegFlux register
[in]a_coarAveCoarsener
[in]a_domainBCDomain BC
[in]a_ebBCBoundary conditions on EBs
[in]a_jumpBcFactoryFactory class for making jump BCs
[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_isMGOperatorIs MG operator or not
[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_jumpOrderStencil order to use on jump cells
[in]a_jumpWeightsWeights for least squares stencils on jump cells
[in]a_smootherWhich smoother to use

Member Function Documentation

◆ AMROperator()

void MFHelmholtzOp::AMROperator ( LevelData< MFCellFAB > &  a_Lphi,
const LevelData< MFCellFAB > &  a_phiFine,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  a_phiCoar,
const bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< MFCellFAB >> *  a_finerOp 
)
overrideprotected

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 MFHelmholtzOp::AMROperatorNC ( LevelData< MFCellFAB > &  a_Lphi,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  a_phiCoar,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< MFCellFAB >> *  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 MFHelmholtzOp::AMROperatorNF ( LevelData< MFCellFAB > &  a_Lphi,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::AMRProlong ( LevelData< MFCellFAB > &  a_correction,
const LevelData< MFCellFAB > &  a_coarseCorrection 
)
finaloverride

Prolongation onto AMR level.

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

◆ AMRResidual()

void MFHelmholtzOp::AMRResidual ( LevelData< MFCellFAB > &  a_residual,
const LevelData< MFCellFAB > &  a_phiFine,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  a_phiCoar,
const LevelData< MFCellFAB > &  a_rhs,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< MFCellFAB >> *  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 MFHelmholtzOp::AMRResidualNC ( LevelData< MFCellFAB > &  a_residual,
const LevelData< MFCellFAB > &  a_phiFine,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  a_rhs,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< MFCellFAB >> *  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 MFHelmholtzOp::AMRResidualNF ( LevelData< MFCellFAB > &  a_residual,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  a_phiCoar,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::AMRRestrict ( LevelData< MFCellFAB > &  a_residualCoarse,
const LevelData< MFCellFAB > &  a_residual,
const LevelData< MFCellFAB > &  a_correction,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::AMRUpdateResidual ( LevelData< MFCellFAB > &  a_residual,
const LevelData< MFCellFAB > &  a_correction,
const LevelData< MFCellFAB > &  a_coarseCorrection 
)
finaloverride

Update AMR residual.

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

◆ applyOp() [1/2]

void MFHelmholtzOp::applyOp ( LevelData< MFCellFAB > &  a_Lphi,
const LevelData< MFCellFAB > &  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() [2/2]

void MFHelmholtzOp::applyOp ( LevelData< MFCellFAB > &  a_Lphi,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > *const  a_phiCoar,
const bool  a_homogeneousPhysBC,
const bool  a_homogeneousCFBC 
)
protected

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!

◆ assign()

void MFHelmholtzOp::assign ( LevelData< MFCellFAB > &  a_lhs,
const LevelData< MFCellFAB > &  a_rhs 
)
finaloverride

Assignment fucntion.

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

◆ assignCopier()

void MFHelmholtzOp::assignCopier ( LevelData< MFCellFAB > &  a_lhs,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::assignLocal ( LevelData< MFCellFAB > &  a_lhs,
const LevelData< MFCellFAB > &  a_rhs 
)
finaloverride

Local assignment function.

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

◆ axby()

void MFHelmholtzOp::axby ( LevelData< MFCellFAB > &  a_lhs,
const LevelData< MFCellFAB > &  a_x,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::buildCopier ( Copier &  a_copier,
const LevelData< MFCellFAB > &  a_lhs,
const LevelData< MFCellFAB > &  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

◆ computeOperatorLoads()

template<typename Duration >
Vector< long long > MFHelmholtzOp::computeOperatorLoads ( LevelData< MFCellFAB > &  a_phi,
const int  a_numApply 
)
noexcept

Time the applyOp routine. Template parameter is std::chrono duration. E.g. std::chrono::microseconds.

Parameters
[in]a_phiDummy cell-centered data. must have the correct number of ghost cells
[in]a_numApplyNumber of times we apply the operator.
Returns
List of execution times in microseconds, ordered along the input DBLs boxArray()

◆ create()

void MFHelmholtzOp::create ( LevelData< MFCellFAB > &  a_lhs,
const LevelData< MFCellFAB > &  a_rhs 
)
finaloverride

Create method.

Parameters
[out]a_lhsClone
[out]a_rhsOriginal data

◆ createCoarsened()

void MFHelmholtzOp::createCoarsened ( LevelData< MFCellFAB > &  a_lhs,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::createCoarser ( LevelData< MFCellFAB > &  a_coarse,
const LevelData< MFCellFAB > &  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

◆ divideByIdentityCoef()

void MFHelmholtzOp::divideByIdentityCoef ( LevelData< MFCellFAB > &  a_rhs)
finaloverride

Divide by the a-coefficient.

Parameters
[in,out]a_rhsDivided data

◆ exchangeGhost()

void MFHelmholtzOp::exchangeGhost ( const LevelData< MFCellFAB > &  a_phi) const
protected

Perform an exchange operation, event if the data is const.

Parameters
[in]a_phiData

◆ fillGrad()

void MFHelmholtzOp::fillGrad ( const LevelData< MFCellFAB > &  a_phi)
finaloverride

Not called, I think.

Parameters
[in]a_phiPhi

◆ getAcoef()

const RefCountedPtr< LevelData< MFCellFAB > > & MFHelmholtzOp::getAcoef ( )

Get the Helmholtz A-coefficient on cell centers.

Returns
m_Acoef

◆ getBcoef()

const RefCountedPtr< LevelData< MFFluxFAB > > & MFHelmholtzOp::getBcoef ( )

Get the Helmholtz B-coefficient on faces.

Returns
m_Bcoef

◆ getBcoefIrreg()

const RefCountedPtr< LevelData< MFBaseIVFAB > > & MFHelmholtzOp::getBcoefIrreg ( )

Get the Helmholtz B-coefficient on the EB.

Returns
m_Bcoef

◆ getFlux()

void MFHelmholtzOp::getFlux ( MFFluxFAB &  a_flux,
const LevelData< MFCellFAB > &  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

◆ incr()

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

Increment function.

Parameters
[in,out]a_lhsOn output, a_lhs = a_lhs + a_rhs*a_scale
[in]a_rhsData
[in]a_scaleScaling factor

◆ interpolateCF()

void MFHelmholtzOp::interpolateCF ( const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > *  a_phiCoar,
const bool  a_homogeneousCF 
)
protected

Do coarse-fine interpolation.

Parameters
[in,out]a_phiFine-level data
[in]a_phiCoarCoarse-level data
[in]a_homogeneousCFHomogeneous interpolation or not.
Note
Only ghost cells in a_phi are touched.

◆ norm()

Real MFHelmholtzOp::norm ( const LevelData< MFCellFAB > &  a_lhs,
int  a_order 
)
finaloverride

Compute solution norm.

Parameters
[in]a_lhsData
[in]a_orderNorm order. Not used.

◆ preCond()

void MFHelmholtzOp::preCond ( LevelData< MFCellFAB > &  a_corr,
const LevelData< MFCellFAB > &  a_residual 
)
finaloverride

Precondition system before bottom solve.

Parameters
[in]a_corrCorrection
[in]a_residualResidual

This just runs a few relaxations.

◆ prolongIncrement()

void MFHelmholtzOp::prolongIncrement ( LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  a_correctCoarse 
)
finaloverride

Prolongation method.

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

◆ relax()

void MFHelmholtzOp::relax ( LevelData< MFCellFAB > &  a_correction,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::relaxGSMultiColor ( LevelData< MFCellFAB > &  a_correction,
const LevelData< MFCellFAB > &  a_residual,
const int  a_iterations 
)

Multi-colored gauss-seidel relaxation.

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

◆ relaxGSRedBlack()

void MFHelmholtzOp::relaxGSRedBlack ( LevelData< MFCellFAB > &  a_correction,
const LevelData< MFCellFAB > &  a_residual,
const int  a_iterations 
)

Jacobi relaxation.

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

◆ relaxPointJacobi()

void MFHelmholtzOp::relaxPointJacobi ( LevelData< MFCellFAB > &  a_correction,
const LevelData< MFCellFAB > &  a_residual,
const int  a_iterations 
)

Jacobi relaxation.

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

◆ residual()

void MFHelmholtzOp::residual ( LevelData< MFCellFAB > &  a_residual,
const LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  a_rhs,
const bool  a_homogeneousPhysBc 
)
finaloverride

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 MFHelmholtzOp::restrictResidual ( LevelData< MFCellFAB > &  a_resCoar,
LevelData< MFCellFAB > &  a_phi,
const LevelData< MFCellFAB > &  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 MFHelmholtzOp::scale ( LevelData< MFCellFAB > &  a_lhs,
const Real &  a_scale 
)
finaloverride

Scale function.

Parameters
[in,out]a_lhsOn output, a_lhs = a_lhs*a_scale
[in]a_scaleScaling factor

◆ setAcoAndBco()

void MFHelmholtzOp::setAcoAndBco ( const RefCountedPtr< LevelData< MFCellFAB >> &  a_Acoef,
const RefCountedPtr< LevelData< MFFluxFAB >> &  a_Bcoef,
const RefCountedPtr< LevelData< MFBaseIVFAB >> &  a_BcoefIrreg 
)

Update operators with new coefficients.

Parameters
[in]a_AcoefOperator A-coefficient
[in]a_BcoefOperator B-coefficient
[in]a_BcoefIrregOperator B-coefficient (on EB faces)

◆ setAlphaAndBeta()

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

Set alpha and beta.

Parameters
[in]a_alphaAlpha
[in]a_betaBeta

◆ setToZero()

void MFHelmholtzOp::setToZero ( LevelData< MFCellFAB > &  a_lhs)
finaloverride

Set to zero.

Parameters
[out]a_lhsData

◆ updateJumpBC()

void MFHelmholtzOp::updateJumpBC ( const LevelData< MFCellFAB > &  a_phi,
const bool  a_homogeneousPhysBC 
)
protected

Update the jump condition.

Parameters
[in]a_phiData
[in]a_homogeneousPhysBChomogeneous physical BC or not

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