chombo-discharge
|
Helmholtz operator for equations like alpha*a(x)*phi(x) + beta*div(b(x)*grad(phi(x))) = rho. More...
#include <CD_EBHelmholtzOp.H>
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. | |
EBHelmholtzOp & | operator= (const EBHelmholtzOp &a_oper)=delete |
No copy assigment allowed. More... | |
EBHelmholtzOp & | operator= (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< EBHelmholtzDomainBC > | m_domainBc |
Domain bc object. | |
RefCountedPtr< EBHelmholtzEBBC > | m_ebBc |
Domain bc object. | |
RefCountedPtr< EBMultigridInterpolator > | m_interpolator |
Interpolator object. | |
RefCountedPtr< EBReflux > | m_fluxReg |
Flux register. | |
RefCountedPtr< EBCoarAve > | m_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. | |
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.
|
delete |
Disallowed copy constructor.
[in] | a_other | Other operator |
|
delete |
Disallowed move constructor.
[in] | a_other | Other operator |
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.
[in] | a_dataLocation | Data location, either cell center or cell centroid |
[in] | a_eblgFine | Fine grids |
[in] | a_eblg | Grids on this level |
[in] | a_eblgCoFi | Coarsening of fine level grids |
[in] | a_eblgCoar | Coarse grids |
[in] | a_eblgCoarMG | Multigrid-grids |
[in] | a_interpolator | Interpolator |
[in] | a_fluxReg | Flux register |
[in] | a_coarAve | Coarsener |
[in] | a_domainBC | Domain BC |
[in] | a_ebBC | Boundary conditions on EBs |
[in] | a_probLo | Lower-left corner of computational domain |
[in] | a_dx | Grid resolution |
[in] | a_refToFine | Refinement ratio to fine level |
[in] | a_refToCoar | Refinement ratio to coarse level |
[in] | a_hasFine | Has fine level or not |
[in] | a_hasCoar | Has coarse level or not |
[in] | a_hasMGObjects | Has multigrid-objects (special objects between AMR levels, or below the AMR levels) |
[in] | a_alpha | Operator alpha |
[in] | a_beta | Operator beta |
[in] | a_Acoef | Operator A-coefficient |
[in] | a_Bcoef | Operator B-coefficient |
[in] | a_BcoefIrreg | Operator B-coefficient (on EB faces) |
[in] | a_ghostCellsPhi | Ghost cells in data holders |
[in] | a_ghostCellsPhi | Ghost cells in data holders |
[in] | a_smoother | Which smoother to use |
|
finaloverride |
Apply the AMR operator, i.e. compute L(phi) in an AMR context.
[out] | a_residual | Residual on this level |
[in] | a_phiFine | Phi on fine level |
[in] | a_phi | Phi on this level |
[in] | a_phiCoar | Phi on coar level |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
[in] | a_finerOp | Finer operatator |
This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.
|
finaloverride |
Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no coarser AMR levels.
[out] | a_Lphi | L(phi) |
[in] | a_phiFine | Phi on finer level |
[in] | a_phi | Phi on this level |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.
|
finaloverride |
Apply the AMR operator, i.e. compute L(phi) in an AMR context, assuming no finer levels.
[out] | a_Lphi | L(phi) |
[in] | a_phi | Phi on this level |
[in] | a_phiCoar | Phi on coar level |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.
|
finaloverride |
Prolongation onto AMR level.
[out] | a_correction | Interpolated correction |
[in] | a_coarseCorrection | Correction on coarse level |
|
finaloverride |
Compute residual on this level. AMR version.
[out] | a_residual | Residual on this level |
[in] | a_phiFine | Phi on fine level |
[in] | a_phi | Phi on this level |
[in] | a_phiCoar | Phi on coar level |
[in] | a_rhs | Right-hand side on this level |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
[in] | a_finerOp | Finer operatator |
|
finaloverride |
Compute AMR residual on coarsest.
[out] | a_residual | Residual on this level |
[in] | a_phiFine | Phi on fine level |
[in] | a_phi | Phi on this level |
[in] | a_rhs | Right-hand side on this level |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
[in] | a_finerOp | Finer operator |
|
finaloverride |
Compute AMR residual on finest AMR level.
[out] | a_residual | Residual on this level |
[in] | a_phi | Phi on this level |
[in] | a_phiCoar | Phi on coar level |
[in] | a_rhs | Right-hand side on this level |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
|
finaloverride |
Restrict residual.
[out] | a_residualCoarse | Coarse residual |
[out] | a_residual | Residual |
[out] | a_correction | Correction on this level |
[out] | a_coarseCorrection | Coarse level correction |
[in] | a_skip_res | I have no idea what this one is supposed to do. |
|
finaloverride |
Update AMR residual.
[in] | a_residual | Residual |
[in] | a_correction | Correction |
[in] | a_coarseCorrection | Coarse-level correction |
|
noexcept |
Apply domain flux.
[in,out] | a_phi | Cell data |
[in] | a_Bcoef | Helmholtz B-coefficient |
[in] | a_cellBox | Computation box |
[in] | a_dit | Data index |
[in] | a_homogeneousPhysBC | Homogeneous BC or not |
|
noexcept |
Apply operator in a grid box.
[out] | a_Lphi | L(phi) |
[in] | a_phi | Phi |
[in] | a_cellBox | Grid box |
[in] | a_dit | Data indxe |
[in] | a_homogeneousPhysBC | Homogeneous physical BCs or not |
|
finaloverride |
Apply operator.
[out] | a_Lphi | L(phi) |
[in] | a_phi | Phi |
[in] | a_homogeneousPhysBc | Homogeneous physical BCs or not |
This computes a_Lphi = L(a_phi) using homogeneous physical BCs or not
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.
[out] | a_Lphi | L(phi) |
[out] | a_phi | Phi on this level |
[out] | a_phiCoar | Coarse-level phi. If you have a coar this |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
[in] | a_homogeneousCFBC | Use homogeneous coarse-fine bcs or not |
|
noexcept |
Apply operator in irregular cells.
[out] | a_Lphi | L(phi) |
[in] | a_phi | Cell-centered data |
[in] | a_cellBox | Cell box |
[in] | a_dit | Grid index |
[in] | a_homogeneousPhysBC | Homogeneous physical BCs or not |
|
noexcept |
Apply operator in regular cells.
[out] | a_Lphi | L(phi) |
[in] | a_phi | Phi |
[in] | a_cellBox | Grid box |
[in] | a_dit | Data indxe |
[in] | a_homogeneousPhysBC | Homogeneous physical BCs or not |
|
finaloverride |
Assign data.
This does a copy from rhs to lhs
|
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)
[out] | a_lhs | Outgoing data |
[in] | a_rhs | Incoming data |
[in] | a_copier | Copier |
|
finaloverride |
Local assignment function.
[out] | a_lhs. | Equal to a_rhs on output |
[in] | a_rhs. | Data |
|
finaloverride |
Set a_lhs = a*x + b*y.
[out] | a_lhs | Result data |
[in] | a_x | x-data |
[in] | a_y | y-data |
[in] | a_a | Scaling factor |
[in] | a_b | Scaling factor |
|
override |
Build copier.
[out] | a_copier | Copier for copying between a_lhs and a_rhs |
[in] | a_lhs | Copying from |
[in] | a_rhs | Copying to |
void EBHelmholtzOp::coarsenCell | ( | LevelData< EBCellFAB > & | a_phi, |
const LevelData< EBCellFAB > & | a_phiFine | ||
) |
Coarsen data from fine to coar level.
[in] | a_phi | Data on this level |
[in] | a_phiFine | Data on finer level |
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.
[in] | a_flux | Flux on this level. |
[in] | a_fineFlux | Flux on finer level |
|
protected |
Compute face-centered fluxes.
[out] | a_fluxCenter | Flux on center |
[in | a_phi Cell-centered data | |
[in] | a_cellBox | Cell-centered compute box. |
[in] | a_dit | Data index |
|
protected |
Compute face-centered fluxes.
[out] | a_fluxCenter | Flux on center |
[in | a_phi Cell-centered data | |
[in] | a_cellBox | Cell-centered compute box. |
[in] | a_dit | Data index |
|
protected |
Fill face centroid-centered fluxes on this level.
[in] | a_phi | Cell-centered data. Must have updated ghost cells for this to make any sense. |
|
protected |
Compute centroid fluxes in a direction.
[out] | a_fluxCentroid | Flux on face centroid |
[in | a_phi Cell-centered data | |
[in] | a_cellBox | Cell-centered compute box. |
[in] | a_dit | Data index |
|
finaloverride |
Create data which clones the layout of the other.
[out] | a_lhs | Data clone (returned data is not initialized) |
[out] | a_rhs | Data layout to be cloned |
|
finaloverride |
Create coarsening of data holder.
[out] | a_lhs | Coarsened data |
[in] | a_rhs | Fine data |
[in] | a_refRat | Coarsening factor |
|
finaloverride |
Create coarsened data.
[out] | a_coarse | Coarse data |
[in] | a_fine | Fine data |
[in] | a_ghosted | Include ghost cells or nto |
|
finaloverride |
Do diagonal scaling.
[in,out] | a_rhs | Data to be scaled |
[in] | a_kappaWeighted | Use kappa-weighted scaling or not |
|
finaloverride |
Divide by the a-coefficient.
[in,out] | a_rhs | Divided data |
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.
[in,out] | a_flux | Flux data holder. |
[in,out] | a_phi | Cell-centered data |
[in] | a_cellBox | Computation box |
[in] | a_dit | Data index |
|
finaloverride |
Not called, I think.
[in] | a_phi | Phi |
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.
[in] | a_phi | Correction |
[in] | a_phi | Lphi L(corrcetion) |
[in] | a_rhs | Residual |
[in] | a_color | "Color": |
|
noexcept |
Multi-color Gauss-Seidel kernel.
[in,out] | a_Lcorr | Storage for computing L(a_corr) |
[in,out] | a_corr | Correction |
[in] | a_resid | Residual |
[in] | a_Acoef | A-coefficient |
[in] | a_Bcoef | B-coefficient |
[in] | a_BcoefIrreg | B-coefficient on EB faces |
[in] | a_cellBox | Grid box |
[in] | a_dit | Data index |
[in] | a_color | Color |
|
noexcept |
Red-black Gauss-Seidel kernel.
[in,out] | a_Lcorr | Storage for computing L(a_corr) |
[in,out] | a_corr | Correction |
[in] | a_resid | Residual |
[in] | a_Acoef | A-coefficient |
[in] | a_Bcoef | B-coefficient |
[in] | a_BcoefIrreg | B-coefficient on EB faces |
[in] | a_redBlack | Red or black |
[in] | a_cellBox | Grid box |
[in] | a_dit | Data index |
const RefCountedPtr< LevelData< EBCellFAB > > & EBHelmholtzOp::getAcoef | ( | ) |
Get the Helmholtz A-coefficient on cell centers.
const RefCountedPtr< LevelData< EBFluxFAB > > & EBHelmholtzOp::getBcoef | ( | ) |
Get the Helmholtz B-coefficient on faces.
const RefCountedPtr< LevelData< BaseIVFAB< Real > > > & EBHelmholtzOp::getBcoefIrreg | ( | ) |
Get the Helmholtz B-coefficient on the EB.
|
protected |
Get the face-centered flux stencil.
[in] | a_face | Face |
[in] | a_dit | Data index. Need because we multiply by B-coefficient |
|
protected |
Get the face-centroid flux stencil.
[in] | a_face | Face |
[in] | a_dit | Data index. Need because we multiply by B-coefficient |
|
finaloverride |
Fill flux.
[out] | a_flux | Flux |
[in] | a_data | Data for which we will compute the flux |
[in] | a_grid | Grid |
[in] | a_dit | Corresponding data index. |
[in] | a_scale | Scaling factor |
void EBHelmholtzOp::homogeneousCFInterp | ( | LevelData< EBCellFAB > & | a_phi | ) |
Do homogeneous coarse-fine interpolation.
[in] | a_phi | Data |
|
finaloverride |
Increment operator.
[in,out] | a_lhs | Data to be incremented |
[in] | a_rhs | Incrementation data |
[in] | a_scale | Scaling factor |
void EBHelmholtzOp::inhomogeneousCFInterp | ( | LevelData< EBCellFAB > & | a_phi, |
const LevelData< EBCellFAB > & | a_phiCoar | ||
) |
Inhomogeneous coarse-fine interpolation.
[in,out] | a_phiFine | Fine data. Ghost cells will be filled. |
[in,out] | a_phiCoar | Coarse data. |
void EBHelmholtzOp::interpolateCF | ( | LevelData< EBCellFAB > & | a_phiFine, |
const LevelData< EBCellFAB > * | a_phiCoar, | ||
const bool | a_homogeneousCFBC | ||
) |
Apply coarse-fine boundary conditions.
[in,out] | a_phiFine | Fine data |
[in,out] | a_phiCoar | Coarse data |
[in] | a_homogeneousCFBC | Homogeneous CF or not (i.e. coar data is zero); |
|
finaloverride |
Compute norm of data.
[in] |
|
delete |
No move assigment allowed.
[in] | a_other | Other operator |
|
delete |
No copy assigment allowed.
[in] | a_other | Other operator |
|
noexcept |
Point Jacobi kernel.
[in,out] | a_Lcorr | Storage for computing L(a_corr) |
[in,out] | a_corr | Correction |
[in] | a_resid | Residual |
[in] | a_Acoef | A-coefficient |
[in] | a_Bcoef | B-coefficient |
[in] | a_BcoefIrreg | B-coefficient on EB faces |
[in] | a_cellBox | Grid box |
[in] | a_dit | Data index |
|
finaloverride |
Precondition system before bottom solve.
[in] | a_corr | Correction |
[in] | a_residual | Residual |
This just runs a few relaxations.
|
finaloverride |
Prolongation method.
[out] | a_phi | Correction on this level |
[out] | a_correctCoarse | Correction on coarse level |
|
protected |
Reflux algorithm.
[in,out] | a_Lphi | L(phiFine, phiCoar). This is modified along the CF interface to consistently account for fine-level fluxes. |
[in] | a_phiFine | Data on finer level |
[in] | a_phi | Data on this operator's level |
[in] | a_finerOp | Operator on fine level |
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.
[out] | a_residual | Residual on this level |
[in] | a_phiFine | Phi on fine level |
[in] | a_phi | Phi on this level |
[in] | a_phiCoar | Phi on coar level |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
[in] | a_finerOp | Finer operatator |
This involves ghost cell interpolation if there's a coarse level, and refluxing if there's a fine level.
|
finaloverride |
Relaxation method. This does smoothing for the system L(correction) = residual.
[in,out] | a_correction | Correction |
[in] | a_residual | Residual |
[in] | a_iterations | Number of iterations |
|
protected |
Multi-colored gauss-seidel relaxation.
[in,out] | a_correction | Correction |
[in] | a_residual | Residual |
[in] | a_iterations | Number of iterations |
|
protected |
Jacobi relaxation.
[in,out] | a_correction | Correction |
[in] | a_residual | Residual |
[in] | a_iterations | Number of iterations |
|
protected |
Jacobi relaxation.
[in,out] | a_correction | Correction |
[in] | a_residual | Residual |
[in] | a_iterations | Number of iterations |
|
override |
Compute residual on this level.
[out] | a_residual | Residual rhs - L(phi) |
[in] | a_phi | phi |
[in] | a_rhs | Right-hand side of system |
[in] | a_homogeneousPhysBC | Use homogeneous physical BCs or not |
|
finaloverride |
Restrict residual onto coarse level.
[in,out] | a_resCoar | Coarse residual |
[in,out] | a_phi | Phi on this level |
[in] | a_rhs | Rhs on this level |
|
finaloverride |
Scale data. Returns a_lhs = a_lhs*a_scale.
[in,out] | a_lhs | Data to be scaled |
[i] | a_scale Scaling factor |
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.
[in] | a_Acoef | Operator A-coefficient |
[in] | a_Bcoef | Operator B-coefficient |
[in] | a_BcoefIrreg | Operator B-coefficient (on EB faces) |
|
finaloverride |
Set alpha coefficient and beta coefficient (can change as diffusion solvers progress)
[in] | a_alpha | Alpha-coefficient |
[in] | a_beta | Beta-coefficient |
|
finaloverride |
Set data to zero.
[in,out] | a_lhs | Data to be set to zero. |
void EBHelmholtzOp::turnOffCoarsening | ( | ) |
Turn off coarsening operation.
void EBHelmholtzOp::turnOffExchange | ( | ) |
Turn off exchange operation.
void EBHelmholtzOp::turnOnCoarsening | ( | ) |
Turn on coarsening operation.
void EBHelmholtzOp::turnOnExchange | ( | ) |
Turn on exchange operation.
|
protected |
For making irregular stencil applications go faster.
This wraps m_relaxStencils in VCAggStencil (which computes explicit stencil offsets)
|
protected |
Turn on/off exchange operation.
|
protected |
Turn on/off exchange operation.
|
protected |
Do coarse-fine interpolation or not.
|
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.