#include <AMRPoissonOp.H>
Public Member Functions | |
virtual Real | dx () const |
AMRPoissonOp functions | |
AMRPoissonOp () | |
virtual | ~AMRPoissonOp () |
void | define (const DisjointBoxLayout &a_grids, const DisjointBoxLayout &a_gridsFiner, const DisjointBoxLayout &a_gridsCoarser, const Real &a_dxLevel, int a_refRatio, int a_refRatioFiner, const ProblemDomain &a_domain, BCHolder a_bc, const Copier &a_exchange, const CFRegion &a_cfregion, const int a_nComp=1) |
void | define (const DisjointBoxLayout &a_grids, const DisjointBoxLayout &a_gridsFiner, const Real &a_dxLevel, int a_refRatio, int a_refRatioFiner, const ProblemDomain &a_domain, BCHolder a_bc, const Copier &a_exchange, const CFRegion &a_cfregion, const int a_nComp=1) |
void | define (const DisjointBoxLayout &a_grids, const DisjointBoxLayout &a_baseBAPtr, const Real &a_dxLevel, int a_refRatio, const ProblemDomain &a_domain, BCHolder a_bc, const Copier &a_exchange, const CFRegion &a_cfregion, int a_numComp=1) |
void | define (const DisjointBoxLayout &a_grids, const Real &a_dx, const ProblemDomain &a_domain, BCHolder a_bc, const Copier &a_exchange, const CFRegion &a_cfregion) |
void | define (const DisjointBoxLayout &a_grids, const Real &a_dx, const ProblemDomain &a_domain, BCHolder a_bc) |
void | define (const DisjointBoxLayout &a_grids, const DisjointBoxLayout *a_baseBAPtr, Real a_dxLevel, int a_refRatio, const ProblemDomain &a_domain, BCHolder a_bc) |
Full define function that mimics the old PoissonOp. | |
virtual void | residual (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs, bool a_homogeneous=false) |
virtual void | residualNF (LevelData< FArrayBox > &a_lhs, LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > *a_phiCoarse, const LevelData< FArrayBox > &a_rhs, bool a_homogeneous=false) |
Residual which uses coarse-level boundary conditions, but ignores finer levels. | |
virtual void | residualI (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs, bool a_homogeneous=false) |
despite what you might think, the "I" here means "Ignore the coarse-fine boundary" | |
virtual void | preCond (LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_residual) |
virtual void | applyOp (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi, bool a_homogeneous=false) |
virtual void | applyOpI (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi, bool a_homogeneous=false) |
despite what you might think, the "I" here means "Ignore the coarse-fine boundary" | |
virtual void | applyOpNoBoundary (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs) |
virtual void | create (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs) |
virtual void | createCoarsened (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs, const int &a_refRat) |
virtual void | assign (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs) |
virtual void | assignLocal (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs) |
virtual void | buildCopier (Copier &a_copier, const LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs) |
virtual void | assignCopier (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs, const Copier &a_copier) |
virtual void | zeroCovered (LevelData< FArrayBox > &a_lhs, LevelData< FArrayBox > &a_rhs, const Copier &a_copier) |
virtual Real | dotProduct (const LevelData< FArrayBox > &a_1, const LevelData< FArrayBox > &a_2) |
virtual void | mDotProduct (const LevelData< FArrayBox > &a_1, const int a_sz, const LevelData< FArrayBox > a_2[], Real a_mdots[]) |
virtual void | incr (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_x, Real a_scale) |
virtual void | axby (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_x, const LevelData< FArrayBox > &a_y, Real a_a, Real a_b) |
virtual void | scale (LevelData< FArrayBox > &a_lhs, const Real &a_scale) |
virtual Real | norm (const LevelData< FArrayBox > &a_x, int a_ord) |
virtual Real | localMaxNorm (const LevelData< FArrayBox > &a_x) |
virtual void | setToZero (LevelData< FArrayBox > &a_x) |
MGLevelOp functions | |
virtual void | relax (LevelData< FArrayBox > &a_e, const LevelData< FArrayBox > &a_residual, int a_iterations) |
virtual void | relaxNF (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > *a_phiCoarse, const LevelData< FArrayBox > &a_rhs, int a_iterations) |
Do CF interpolation then relax as normal. | |
virtual void | createCoarser (LevelData< FArrayBox > &a_coarse, const LevelData< FArrayBox > &a_fine, bool a_ghosted) |
virtual void | restrictResidual (LevelData< FArrayBox > &a_resCoarse, LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_rhsFine) |
virtual void | prolongIncrement (LevelData< FArrayBox > &a_phiThisLevel, const LevelData< FArrayBox > &a_correctCoarse) |
AMRLevelOp functions | |
virtual int | refToCoarser () |
virtual void | AMRResidual (LevelData< FArrayBox > &a_residual, const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_phiCoarse, const LevelData< FArrayBox > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp) |
virtual void | AMRResidualNC (LevelData< FArrayBox > &a_residual, const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp) |
virtual void | AMRResidualNF (LevelData< FArrayBox > &a_residual, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_phiCoarse, const LevelData< FArrayBox > &a_rhs, bool a_homogeneousPhysBC) |
virtual void | AMROperator (LevelData< FArrayBox > &a_LofPhi, const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_phiCoarse, bool a_homogeneousDomBC, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp) |
virtual void | AMROperatorNC (LevelData< FArrayBox > &a_LofPhi, const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, bool a_homogeneousBC, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp) |
virtual void | AMROperatorNF (LevelData< FArrayBox > &a_LofPhi, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_phiCoarse, bool a_homogeneousBC) |
virtual void | AMRRestrict (LevelData< FArrayBox > &a_resCoarse, const LevelData< FArrayBox > &a_residual, const LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_coarseCorrection, bool a_skip_res=false) |
virtual void | AMRRestrictS (LevelData< FArrayBox > &a_resCoarse, const LevelData< FArrayBox > &a_residual, const LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_coarseCorrection, LevelData< FArrayBox > &a_scratch, bool a_skip_res=false) |
virtual void | AMRProlong (LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_coarseCorrection) |
virtual void | AMRProlongS (LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_coarseCorrection, LevelData< FArrayBox > &a_temp, const Copier &a_copier) |
virtual void | AMRProlongS_2 (LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_coarseCorrection, LevelData< FArrayBox > &a_temp, const Copier &a_copier, const Copier &a_cornerCopier, const AMRLevelOp< LevelData< FArrayBox > > *a_crsOp) |
virtual void | AMRUpdateResidual (LevelData< FArrayBox > &a_residual, const LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_coarseCorrection) |
virtual Real | AMRNorm (const LevelData< FArrayBox > &a_coarseResid, const LevelData< FArrayBox > &a_fineResid, const int &a_refRat, const int &a_ord) |
virtual void | setAlphaAndBeta (const Real &a_alpha, const Real &a_beta) |
For tga to reset stuff. | |
virtual void | setBC (const BCHolder &a_bc) |
Change boundary conditions. | |
virtual void | diagonalScale (LevelData< FArrayBox > &a_rhs, bool a_kappaWeighted) |
For tga stuff---in this case a noop. | |
virtual void | divideByIdentityCoef (LevelData< FArrayBox > &a_rhs) |
For tga stuff---in this case a noop. | |
virtual void | fillGrad (const LevelData< FArrayBox > &a_phi) |
These functions are part of the LevelTGA interface...... | |
virtual void | reflux (const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, LevelData< FArrayBox > &a_residual, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp) |
virtual void | getFlux (FluxBox &a_flux, const LevelData< FArrayBox > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale) |
virtual void | write (const LevelData< FArrayBox > *a_data, const char *a_filename) |
Public Attributes | |
Real | m_alpha |
public constants | |
Real | m_beta |
Real | m_aCoef |
Real | m_bCoef |
Real | m_dxCrse |
Vector< IntVect > | m_colors |
Static Public Attributes | |
static int | s_exchangeMode |
static int | s_relaxMode |
static int | s_maxCoarse |
static int | s_prolongType |
Protected Member Functions | |
virtual void | levelGSRB (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs) |
virtual void | levelMultiColor (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs) |
virtual void | looseGSRB (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs) |
virtual void | overlapGSRB (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs) |
virtual void | levelGSRBLazy (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs) |
virtual void | levelJacobi (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs) |
virtual void | homogeneousCFInterp (LevelData< FArrayBox > &a_phif) |
virtual void | homogeneousCFInterp (LevelData< FArrayBox > &a_phif, const DataIndex &a_datInd, int a_idir, Side::LoHiSide a_hiorlo) |
virtual void | singleBoxCFInterp (FArrayBox &a_phi) |
virtual void | interpOnIVSHomo (LevelData< FArrayBox > &a_phif, const DataIndex &a_datInd, const int a_idir, const Side::LoHiSide a_hiorlo, const IntVectSet &a_interpIVS) |
virtual void | getFlux (FArrayBox &a_flux, const FArrayBox &a_data, const Box &a_edgebox, int a_dir, int a_ref=1) const |
virtual void | getFlux (FArrayBox &a_flux, const FArrayBox &a_data, int a_dir, int a_ref=1) const |
Protected Attributes | |
Real | m_dx |
ProblemDomain | m_domain |
LevelDataOps< FArrayBox > | m_levelOps |
BCHolder | m_bc |
CFRegion | m_cfregion |
Copier | m_exchangeCopier |
QuadCFInterp | m_interpWithCoarser |
LevelFluxRegister | m_levfluxreg |
DisjointBoxLayout | m_coarsenedMGrids |
int | m_refToCoarser |
int | m_refToFiner |
AMRPoissonOp::AMRPoissonOp | ( | ) | [inline] |
virtual AMRPoissonOp::~AMRPoissonOp | ( | ) | [inline, virtual] |
void AMRPoissonOp::define | ( | const DisjointBoxLayout & | a_grids, | |
const DisjointBoxLayout & | a_gridsFiner, | |||
const DisjointBoxLayout & | a_gridsCoarser, | |||
const Real & | a_dxLevel, | |||
int | a_refRatio, | |||
int | a_refRatioFiner, | |||
const ProblemDomain & | a_domain, | |||
BCHolder | a_bc, | |||
const Copier & | a_exchange, | |||
const CFRegion & | a_cfregion, | |||
const int | a_nComp = 1 | |||
) |
full define function for AMRLevelOp with both coarser and finer levels
void AMRPoissonOp::define | ( | const DisjointBoxLayout & | a_grids, | |
const DisjointBoxLayout & | a_gridsFiner, | |||
const Real & | a_dxLevel, | |||
int | a_refRatio, | |||
int | a_refRatioFiner, | |||
const ProblemDomain & | a_domain, | |||
BCHolder | a_bc, | |||
const Copier & | a_exchange, | |||
const CFRegion & | a_cfregion, | |||
const int | a_nComp = 1 | |||
) |
full define function for AMRLevelOp with finer levels, but no coarser
void AMRPoissonOp::define | ( | const DisjointBoxLayout & | a_grids, | |
const DisjointBoxLayout & | a_baseBAPtr, | |||
const Real & | a_dxLevel, | |||
int | a_refRatio, | |||
const ProblemDomain & | a_domain, | |||
BCHolder | a_bc, | |||
const Copier & | a_exchange, | |||
const CFRegion & | a_cfregion, | |||
int | a_numComp = 1 | |||
) |
define function for AMRLevelOp which has no finer AMR level Jamie Parkinson added a_numComp parameter so that we define CF interp objects correctly
void AMRPoissonOp::define | ( | const DisjointBoxLayout & | a_grids, | |
const Real & | a_dx, | |||
const ProblemDomain & | a_domain, | |||
BCHolder | a_bc, | |||
const Copier & | a_exchange, | |||
const CFRegion & | a_cfregion | |||
) |
define function for AMRLevelOp which has no finer or coarser AMR level
void AMRPoissonOp::define | ( | const DisjointBoxLayout & | a_grids, | |
const Real & | a_dx, | |||
const ProblemDomain & | a_domain, | |||
BCHolder | a_bc | |||
) |
define function for AMRLevelOp which has no finer or coarser AMR level
void AMRPoissonOp::define | ( | const DisjointBoxLayout & | a_grids, | |
const DisjointBoxLayout * | a_baseBAPtr, | |||
Real | a_dxLevel, | |||
int | a_refRatio, | |||
const ProblemDomain & | a_domain, | |||
BCHolder | a_bc | |||
) |
Full define function that mimics the old PoissonOp.
Makes all coarse-fine information and sets internal variables
virtual void AMRPoissonOp::residual | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > & | a_rhs, | |||
bool | a_homogeneous = false | |||
) | [virtual] |
Say you are solving L(phi) = rhs. Make a_lhs = L(a_phi) - a_rhs. If a_homogeneous is true, evaluate the operator using homogeneous boundary conditions.
Implements LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::residualNF | ( | LevelData< FArrayBox > & | a_lhs, | |
LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > * | a_phiCoarse, | |||
const LevelData< FArrayBox > & | a_rhs, | |||
bool | a_homogeneous = false | |||
) | [virtual] |
Residual which uses coarse-level boundary conditions, but ignores finer levels.
useful for full-multigrid types of schemes
Reimplemented from MGLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::residualI | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > & | a_rhs, | |||
bool | a_homogeneous = false | |||
) | [virtual] |
despite what you might think, the "I" here means "Ignore the coarse-fine boundary"
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::preCond | ( | LevelData< FArrayBox > & | a_cor, | |
const LevelData< FArrayBox > & | a_residual | |||
) | [virtual] |
Given the current state of the residual the correction, apply your preconditioner to a_cor.
Implements LinearOp< LevelData< FArrayBox > >.
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::applyOp | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_phi, | |||
bool | a_homogeneous = false | |||
) | [virtual] |
In the context of solving L(phi) = rhs, set a_lhs = L(a_phi). If a_homogeneous is true, evaluate the operator using homogeneous boundary conditions.
Implements LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::applyOpI | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_phi, | |||
bool | a_homogeneous = false | |||
) | [virtual] |
despite what you might think, the "I" here means "Ignore the coarse-fine boundary"
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::applyOpNoBoundary | ( | LevelData< FArrayBox > & | a_ans, | |
const LevelData< FArrayBox > & | a_phi | |||
) | [virtual] |
Apply the differential operator without any boundary or coarse-fine boundary conditions and no finer level
a_ans | The result of the application of the operator to a_phi. | |
a_phi | The data to which the operator is applied. |
Implements TGAHelmOp< LevelData< FArrayBox > >.
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::create | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [virtual] |
Creat data holder a_lhs that mirrors a_rhs. You do not need to copy the data of a_rhs, just make a holder the same size.
Implements LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::createCoarsened | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_rhs, | |||
const int & | a_refRat | |||
) | [virtual] |
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::assign | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [virtual] |
Set a_lhs equal to a_rhs.
Implements LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::assignLocal | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [virtual] |
Reimplemented from LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::buildCopier | ( | Copier & | a_copier, | |
const LevelData< FArrayBox > & | a_lhs, | |||
const LevelData< FArrayBox > & | a_rhs | |||
) | [virtual] |
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::assignCopier | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_rhs, | |||
const Copier & | a_copier | |||
) | [virtual] |
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::zeroCovered | ( | LevelData< FArrayBox > & | a_lhs, | |
LevelData< FArrayBox > & | a_rhs, | |||
const Copier & | a_copier | |||
) | [virtual] |
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
virtual Real AMRPoissonOp::dotProduct | ( | const LevelData< FArrayBox > & | a_1, | |
const LevelData< FArrayBox > & | a_2 | |||
) | [virtual] |
Compute and return the dot product of a_1 and a_2. In most contexts, this means return the sum over all data points of a_1*a_2.
Implements LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::mDotProduct | ( | const LevelData< FArrayBox > & | a_1, | |
const int | a_sz, | |||
const LevelData< FArrayBox > | a_2[], | |||
Real | a_mdots[] | |||
) | [virtual] |
Reimplemented from LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::incr | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_x, | |||
Real | a_scale | |||
) | [virtual] |
Increment by scaled amount (a_lhs += a_scale*a_x).
Implements LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::axby | ( | LevelData< FArrayBox > & | a_lhs, | |
const LevelData< FArrayBox > & | a_x, | |||
const LevelData< FArrayBox > & | a_y, | |||
Real | a_a, | |||
Real | a_b | |||
) | [virtual] |
Set input to a scaled sum (a_lhs = a_a*a_x + a_b*a_y).
Implements LinearOp< LevelData< FArrayBox > >.
Multiply the input by a given scale (a_lhs *= a_scale).
Implements LinearOp< LevelData< FArrayBox > >.
Return the norm of a_rhs. a_ord == 0 max norm, a_ord == 1 sum(abs(a_rhs)), else, L(a_ord) norm.
Implements LinearOp< LevelData< FArrayBox > >.
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
Set a_lhs to zero.
Implements LinearOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::relax | ( | LevelData< FArrayBox > & | a_correction, | |
const LevelData< FArrayBox > & | a_residual, | |||
int | a_iterations | |||
) | [virtual] |
Use your relaxation operator to remove the high frequency wave numbers from the correction so that it may be averaged to a coarser refinement. A point relaxtion scheme, for example takes the form a_correction -= lambda*(L(a_correction) - a_residual).
Implements MGLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::relaxNF | ( | LevelData< FArrayBox > & | a_phi, | |
const LevelData< FArrayBox > * | a_phiCoarse, | |||
const LevelData< FArrayBox > & | a_rhs, | |||
int | a_iterations | |||
) | [virtual] |
Do CF interpolation then relax as normal.
useful for full-multigrid/FAS schemes
Reimplemented from MGLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::createCoarser | ( | LevelData< FArrayBox > & | a_coarse, | |
const LevelData< FArrayBox > & | a_fine, | |||
bool | ghosted | |||
) | [virtual] |
Create a coarsened (by two) version of the input data. This does not include averaging the data. So if a_fine is over a Box of (0, 0, 0) (63, 63, 63), a_fine should be over a Box (0, 0, 0) (31, 31, 31).
Implements MGLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::restrictResidual | ( | LevelData< FArrayBox > & | a_resCoarse, | |
LevelData< FArrayBox > & | a_phiFine, | |||
const LevelData< FArrayBox > & | a_rhsFine | |||
) | [virtual] |
calculate restricted residual a_resCoarse[2h] = I[h->2h] (rhsFine[h] - L[h](phiFine[h])
Implements MGLevelOp< LevelData< FArrayBox > >.
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::prolongIncrement | ( | LevelData< FArrayBox > & | a_phiThisLevel, | |
const LevelData< FArrayBox > & | a_correctCoarse | |||
) | [virtual] |
correct the fine solution based on coarse correction a_phiThisLevel += I[2h->h](a_correctCoarse)
Implements MGLevelOp< LevelData< FArrayBox > >.
virtual int AMRPoissonOp::refToCoarser | ( | ) | [inline, virtual] |
returns 1 when there are no coarser AMRLevelOp objects
Implements AMRLevelOp< LevelData< FArrayBox > >.
References m_refToCoarser.
virtual void AMRPoissonOp::AMRResidual | ( | LevelData< FArrayBox > & | a_residual, | |
const LevelData< FArrayBox > & | a_phiFine, | |||
const LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > & | a_phiCoarse, | |||
const LevelData< FArrayBox > & | a_rhs, | |||
bool | a_homogeneousPhysBC, | |||
AMRLevelOp< LevelData< FArrayBox > > * | a_finerOp | |||
) | [virtual] |
a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse)
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRResidualNC | ( | LevelData< FArrayBox > & | a_residual, | |
const LevelData< FArrayBox > & | a_phiFine, | |||
const LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > & | a_rhs, | |||
bool | a_homogeneousPhysBC, | |||
AMRLevelOp< LevelData< FArrayBox > > * | a_finerOp | |||
) | [virtual] |
residual assuming no more coarser AMR levels
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRResidualNF | ( | LevelData< FArrayBox > & | a_residual, | |
const LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > & | a_phiCoarse, | |||
const LevelData< FArrayBox > & | a_rhs, | |||
bool | a_homogeneousPhysBC | |||
) | [virtual] |
a_residual = a_rhs - L(a_phi, a_phiCoarse)
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMROperator | ( | LevelData< FArrayBox > & | a_LofPhi, | |
const LevelData< FArrayBox > & | a_phiFine, | |||
const LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > & | a_phiCoarse, | |||
bool | a_homogeneousDomBC, | |||
AMRLevelOp< LevelData< FArrayBox > > * | a_finerOp | |||
) | [virtual] |
Apply the AMR operator, including coarse-fine matching
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMROperatorNC | ( | LevelData< FArrayBox > & | a_LofPhi, | |
const LevelData< FArrayBox > & | a_phiFine, | |||
const LevelData< FArrayBox > & | a_phi, | |||
bool | a_homogeneousBC, | |||
AMRLevelOp< LevelData< FArrayBox > > * | a_finerOp | |||
) | [virtual] |
Apply the AMR operator, including coarse-fine matching assume no coarser AMR level
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMROperatorNF | ( | LevelData< FArrayBox > & | a_LofPhi, | |
const LevelData< FArrayBox > & | a_phi, | |||
const LevelData< FArrayBox > & | a_phiCoarse, | |||
bool | a_homogeneousBC | |||
) | [virtual] |
Apply the AMR operator, including coarse-fine matching. assume no finer AMR level
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRRestrict | ( | LevelData< FArrayBox > & | a_resCoarse, | |
const LevelData< FArrayBox > & | a_residual, | |||
const LevelData< FArrayBox > & | a_correction, | |||
const LevelData< FArrayBox > & | a_coarseCorrection, | |||
bool | a_skip_res = false | |||
) | [virtual] |
a_resCoarse = I[h-2h]( a_residual - L(a_correction, a_coarseCorrection)) it is assumed that a_resCoarse has already been filled in with the coarse version of AMRResidualNF and that this operation is free to overwrite in the overlap regions.
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRRestrictS | ( | LevelData< FArrayBox > & | a_resCoarse, | |
const LevelData< FArrayBox > & | a_residual, | |||
const LevelData< FArrayBox > & | a_correction, | |||
const LevelData< FArrayBox > & | a_coarseCorrection, | |||
LevelData< FArrayBox > & | a_scratch, | |||
bool | a_skip_res = false | |||
) | [virtual] |
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRProlong | ( | LevelData< FArrayBox > & | a_correction, | |
const LevelData< FArrayBox > & | a_coarseCorrection | |||
) | [virtual] |
a_correction += I[h->h](a_coarseCorrection)
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRProlongS | ( | LevelData< FArrayBox > & | a_correction, | |
const LevelData< FArrayBox > & | a_coarseCorrection, | |||
LevelData< FArrayBox > & | a_temp, | |||
const Copier & | a_copier | |||
) | [virtual] |
optimization of AMRProlong that sends in the existing temporary and copier
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRProlongS_2 | ( | LevelData< FArrayBox > & | a_correction, | |
const LevelData< FArrayBox > & | a_coarseCorrection, | |||
LevelData< FArrayBox > & | a_temp, | |||
const Copier & | a_copier, | |||
const Copier & | a_cornerCopier, | |||
const AMRLevelOp< LevelData< FArrayBox > > * | a_crsOp | |||
) | [virtual] |
optimization of AMRProlong that sends in the existing temporary and copier -- higher order
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
virtual void AMRPoissonOp::AMRUpdateResidual | ( | LevelData< FArrayBox > & | a_residual, | |
const LevelData< FArrayBox > & | a_correction, | |||
const LevelData< FArrayBox > & | a_coarseCorrection | |||
) | [virtual] |
a_residual = a_residual - L(a_correction, a_coarseCorrection)
Implements AMRLevelOp< LevelData< FArrayBox > >.
virtual Real AMRPoissonOp::AMRNorm | ( | const LevelData< FArrayBox > & | a_coarseResid, | |
const LevelData< FArrayBox > & | a_fineResid, | |||
const int & | a_refRat, | |||
const int & | a_ord | |||
) | [virtual] |
compute norm over all cells on coarse not covered by finer
Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.
For tga to reset stuff.
Implements TGAHelmOp< LevelData< FArrayBox > >.
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::setBC | ( | const BCHolder & | a_bc | ) | [virtual] |
Change boundary conditions.
virtual void AMRPoissonOp::diagonalScale | ( | LevelData< FArrayBox > & | a_rhs, | |
bool | a_kappaWeighted | |||
) | [inline, virtual] |
For tga stuff---in this case a noop.
Reimplemented from TGAHelmOp< LevelData< FArrayBox > >.
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::divideByIdentityCoef | ( | LevelData< FArrayBox > & | a_rhs | ) | [inline, virtual] |
For tga stuff---in this case a noop.
Implements TGAHelmOp< LevelData< FArrayBox > >.
Reimplemented in VCAMRPoissonOp2.
These functions are part of the LevelTGA interface......
Implements LevelTGAHelmOp< LevelData< FArrayBox >, FluxBox >.
virtual void AMRPoissonOp::reflux | ( | const LevelData< FArrayBox > & | a_phiFine, | |
const LevelData< FArrayBox > & | a_phi, | |||
LevelData< FArrayBox > & | a_residual, | |||
AMRLevelOp< LevelData< FArrayBox > > * | a_finerOp | |||
) | [virtual] |
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::getFlux | ( | FluxBox & | a_flux, | |
const LevelData< FArrayBox > & | a_data, | |||
const Box & | a_grid, | |||
const DataIndex & | a_dit, | |||
Real | a_scale | |||
) | [inline, virtual] |
Implements LevelTGAHelmOp< LevelData< FArrayBox >, FluxBox >.
Reimplemented in VCAMRPoissonOp2.
References FluxBox::define(), BoxLayoutData< T >::nComp(), and SpaceDim.
virtual void AMRPoissonOp::write | ( | const LevelData< FArrayBox > * | a, | |
const char * | filename | |||
) | [virtual] |
Debugging aid for solvers. Print out a "T" to a file named "filename" default implementation is to print out a message saying "LinearOp::write not implemented"
Reimplemented from LinearOp< LevelData< FArrayBox > >.
virtual Real AMRPoissonOp::dx | ( | ) | const [inline, virtual] |
Return dx at this level of refinement
Reimplemented from LinearOp< LevelData< FArrayBox > >.
References m_dx.
virtual void AMRPoissonOp::levelGSRB | ( | LevelData< FArrayBox > & | a_phi, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [protected, virtual] |
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::levelMultiColor | ( | LevelData< FArrayBox > & | a_phi, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [protected, virtual] |
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::looseGSRB | ( | LevelData< FArrayBox > & | a_phi, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [protected, virtual] |
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::overlapGSRB | ( | LevelData< FArrayBox > & | a_phi, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [protected, virtual] |
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::levelGSRBLazy | ( | LevelData< FArrayBox > & | a_phi, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [protected, virtual] |
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::levelJacobi | ( | LevelData< FArrayBox > & | a_phi, | |
const LevelData< FArrayBox > & | a_rhs | |||
) | [protected, virtual] |
Reimplemented in VCAMRPoissonOp2.
virtual void AMRPoissonOp::homogeneousCFInterp | ( | LevelData< FArrayBox > & | a_phif | ) | [protected, virtual] |
virtual void AMRPoissonOp::homogeneousCFInterp | ( | LevelData< FArrayBox > & | a_phif, | |
const DataIndex & | a_datInd, | |||
int | a_idir, | |||
Side::LoHiSide | a_hiorlo | |||
) | [protected, virtual] |
virtual void AMRPoissonOp::singleBoxCFInterp | ( | FArrayBox & | a_phi | ) | [protected, virtual] |
virtual void AMRPoissonOp::interpOnIVSHomo | ( | LevelData< FArrayBox > & | a_phif, | |
const DataIndex & | a_datInd, | |||
const int | a_idir, | |||
const Side::LoHiSide | a_hiorlo, | |||
const IntVectSet & | a_interpIVS | |||
) | [protected, virtual] |
virtual void AMRPoissonOp::getFlux | ( | FArrayBox & | a_flux, | |
const FArrayBox & | a_data, | |||
const Box & | a_edgebox, | |||
int | a_dir, | |||
int | a_ref = 1 | |||
) | const [protected, virtual] |
virtual void AMRPoissonOp::getFlux | ( | FArrayBox & | a_flux, | |
const FArrayBox & | a_data, | |||
int | a_dir, | |||
int | a_ref = 1 | |||
) | const [protected, virtual] |
public constants
Reimplemented in VCAMRPoissonOp2.
Reimplemented in VCAMRPoissonOp2.
int AMRPoissonOp::s_exchangeMode [static] |
int AMRPoissonOp::s_relaxMode [static] |
int AMRPoissonOp::s_maxCoarse [static] |
int AMRPoissonOp::s_prolongType [static] |
Real AMRPoissonOp::m_dx [protected] |
Referenced by dx().
ProblemDomain AMRPoissonOp::m_domain [protected] |
LevelDataOps<FArrayBox> AMRPoissonOp::m_levelOps [protected] |
BCHolder AMRPoissonOp::m_bc [protected] |
CFRegion AMRPoissonOp::m_cfregion [protected] |
Copier AMRPoissonOp::m_exchangeCopier [protected] |
QuadCFInterp AMRPoissonOp::m_interpWithCoarser [protected] |
LevelFluxRegister AMRPoissonOp::m_levfluxreg [protected] |
DisjointBoxLayout AMRPoissonOp::m_coarsenedMGrids [protected] |
int AMRPoissonOp::m_refToCoarser [protected] |
Referenced by refToCoarser().
int AMRPoissonOp::m_refToFiner [protected] |