chombo-discharge
|
Operator for solving multifluid Helmholtz on a grid level. More...
#include <CD_MFHelmholtzOp.H>
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. | |
MFHelmholtzOp & | operator= (const MFHelmholtzOp &a_oper)=delete |
No copy assigment allowed. | |
MFHelmholtzOp & | operator= (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< MFHelmholtzJumpBC > | m_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. | |
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.
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.
[in] | a_dataLocation | Data location, either cell center or cell centroid |
[in] | a_mflgFine | Fine grids |
[in] | a_mflg | Grids on this level |
[in] | a_mflgCoFi | Coarsening of fine level grids |
[in] | a_mflgCoar | Coarse grids |
[in] | a_mflgCoarMG | 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_jumpBcFactory | Factory class for making jump BCs |
[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_isMGOperator | Is MG operator or not |
[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_jumpOrder | Stencil order to use on jump cells |
[in] | a_jumpWeights | Weights for least squares stencils on jump cells |
[in] | a_smoother | Which smoother to use |
|
overrideprotected |
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 |
|
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
|
protected |
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 |
|
finaloverride |
Assignment fucntion.
[out] | a_lhs. | Equal to a_rhs on output |
[in] | a_rhs. | Data |
|
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 |
|
noexcept |
Time the applyOp routine. Template parameter is std::chrono duration. E.g. std::chrono::microseconds.
[in] | a_phi | Dummy cell-centered data. must have the correct number of ghost cells |
[in] | a_numApply | Number of times we apply the operator. |
|
finaloverride |
Create method.
[out] | a_lhs | Clone |
[out] | a_rhs | Original data |
|
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 |
Divide by the a-coefficient.
[in,out] | a_rhs | Divided data |
|
protected |
Perform an exchange operation, event if the data is const.
[in] | a_phi | Data |
|
finaloverride |
Not called, I think.
[in] | a_phi | Phi |
const RefCountedPtr< LevelData< MFCellFAB > > & MFHelmholtzOp::getAcoef | ( | ) |
Get the Helmholtz A-coefficient on cell centers.
const RefCountedPtr< LevelData< MFFluxFAB > > & MFHelmholtzOp::getBcoef | ( | ) |
Get the Helmholtz B-coefficient on faces.
const RefCountedPtr< LevelData< MFBaseIVFAB > > & MFHelmholtzOp::getBcoefIrreg | ( | ) |
Get the Helmholtz B-coefficient on the EB.
|
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 |
|
finaloverride |
Increment function.
[in,out] | a_lhs | On output, a_lhs = a_lhs + a_rhs*a_scale |
[in] | a_rhs | Data |
[in] | a_scale | Scaling factor |
|
protected |
Do coarse-fine interpolation.
[in,out] | a_phi | Fine-level data |
[in] | a_phiCoar | Coarse-level data |
[in] | a_homogeneousCF | Homogeneous interpolation or not. |
|
finaloverride |
Compute solution norm.
[in] | a_lhs | Data |
[in] | a_order | Norm order. Not used. |
|
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 |
|
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 |
void MFHelmholtzOp::relaxGSMultiColor | ( | LevelData< MFCellFAB > & | a_correction, |
const LevelData< MFCellFAB > & | a_residual, | ||
const int | a_iterations | ||
) |
Multi-colored gauss-seidel relaxation.
[in,out] | a_correction | Correction |
[in] | a_residual | Residual |
[in] | a_iterations | Number of iterations |
void MFHelmholtzOp::relaxGSRedBlack | ( | LevelData< MFCellFAB > & | a_correction, |
const LevelData< MFCellFAB > & | a_residual, | ||
const int | a_iterations | ||
) |
Jacobi relaxation.
[in,out] | a_correction | Correction |
[in] | a_residual | Residual |
[in] | a_iterations | Number of iterations |
void MFHelmholtzOp::relaxPointJacobi | ( | LevelData< MFCellFAB > & | a_correction, |
const LevelData< MFCellFAB > & | a_residual, | ||
const int | a_iterations | ||
) |
Jacobi relaxation.
[in,out] | a_correction | Correction |
[in] | a_residual | Residual |
[in] | a_iterations | Number of iterations |
|
finaloverride |
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 function.
[in,out] | a_lhs | On output, a_lhs = a_lhs*a_scale |
[in] | a_scale | Scaling factor |
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.
[in] | a_Acoef | Operator A-coefficient |
[in] | a_Bcoef | Operator B-coefficient |
[in] | a_BcoefIrreg | Operator B-coefficient (on EB faces) |
|
finaloverride |
Set alpha and beta.
[in] | a_alpha | Alpha |
[in] | a_beta | Beta |
|
finaloverride |
Set to zero.
[out] | a_lhs | Data |
|
protected |
Update the jump condition.
[in] | a_phi | Data |
[in] | a_homogeneousPhysBC | homogeneous physical BC or not |