EBConductivityOp Class Reference

#include <EBConductivityOp.H>

Inheritance diagram for EBConductivityOp:

Inheritance graph
[legend]
Collaboration diagram for EBConductivityOp:

Collaboration graph
[legend]

List of all members.


Detailed Description

Operator to solve (alpha a + beta div (b grad) )phi = rhs. This follows the AMRLevelOp interface.

Public Member Functions

virtual void setAlphaAndBeta (const Real &a_alpha, const Real &a_beta)
 set the constants in the equation
virtual void diagonalScale (LevelData< EBCellFAB > &a_rhs)
virtual void fillGrad (const LevelData< EBCellFAB > &a_phi)
 a leveltgaism
virtual void getFlux (EBFluxFAB &a_flux, const LevelData< EBCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale)
 a leveltgaism
void getFlux (EBFaceFAB &a_fluxCentroid, const EBCellFAB &a_phi, const Box &a_ghostedBox, const Box &a_fabBox, const ProblemDomain &a_domain, const EBISBox &a_ebisBox, const Real &a_dx, const DataIndex &a_datInd, const int &a_idir)
virtual ~EBConductivityOp ()
 EBConductivityOp ()
void AMRResidualNC (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
void AMROperatorNC (LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
 EBConductivityOp (const EBLevelGrid &a_eblgFine, const EBLevelGrid &a_eblg, const EBLevelGrid &a_eblgCoar, const EBLevelGrid &a_eblgCoarMG, const RefCountedPtr< EBQuadCFInterp > &a_quadCFI, const RefCountedPtr< ConductivityBaseDomainBC > &a_domainBC, const RefCountedPtr< ConductivityBaseEBBC > &a_ebBC, const Real &a_dx, const Real &a_dxCoar, const int &a_refToFine, const int &a_refToCoar, const bool &a_hasFine, const bool &a_hasCoar, const bool &a_hasMGObjects, const bool &a_layoutChanged, 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)
virtual void residual (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousPhysBC=false)
virtual void preCond (LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
virtual void applyOp (LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoarse, const bool &a_homogeneousPhysBC, const bool &a_homogeneousCFBC)
virtual void applyOpNoBoundary (LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi)
 virtual function called by LevelTGA
virtual void applyOp (LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousPhysBC)
virtual void create (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
virtual void createCoarsened (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs, const int &a_refRat)
Real AMRNorm (const LevelData< EBCellFAB > &a_coarResid, const LevelData< EBCellFAB > &a_fineResid, const int &a_refRat, const int &a_ord)
virtual void assign (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_rhs)
virtual Real dotProduct (const LevelData< EBCellFAB > &a_1, const LevelData< EBCellFAB > &a_2)
virtual void incr (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, Real a_scale)
virtual void axby (LevelData< EBCellFAB > &a_lhs, const LevelData< EBCellFAB > &a_x, const LevelData< EBCellFAB > &a_y, Real a_a, Real a_b)
virtual void scale (LevelData< EBCellFAB > &a_lhs, const Real &a_scale)
virtual Real norm (const LevelData< EBCellFAB > &a_rhs, int a_ord)
virtual void setToZero (LevelData< EBCellFAB > &a_lhs)
virtual void setVal (LevelData< EBCellFAB > &a_lhs, const Real &a_value)
virtual void createCoarser (LevelData< EBCellFAB > &a_coarse, const LevelData< EBCellFAB > &a_fine, bool a_ghosted)
virtual void relax (LevelData< EBCellFAB > &a_e, const LevelData< EBCellFAB > &a_residual, int a_iterations)
virtual void restrictResidual (LevelData< EBCellFAB > &a_resCoarse, LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_rhsFine)
virtual void prolongIncrement (LevelData< EBCellFAB > &a_phiThisLevel, const LevelData< EBCellFAB > &a_correctCoarse)
virtual int refToCoarser ()
virtual int refToFiner ()
virtual void AMRResidual (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
virtual void AMRResidualNF (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, const LevelData< EBCellFAB > &a_rhs, bool a_homogeneousBC)
virtual void AMROperator (LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, bool a_homogeneousBC, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
virtual void AMROperatorNF (LevelData< EBCellFAB > &a_LofPhi, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_phiCoarse, bool a_homogeneousBC)
virtual void AMRRestrict (LevelData< EBCellFAB > &a_resCoarse, const LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
virtual void AMRProlong (LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
virtual void AMRUpdateResidual (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_correction, const LevelData< EBCellFAB > &a_coarseCorrection)
void reflux (LevelData< EBCellFAB > &a_residual, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
void gsrbColor (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_lph, const LevelData< EBCellFAB > &a_rhs, const IntVect &a_color)
void getDivFStencil (VoFStencil &a_vofStencil, const VolIndex &a_vof, const DataIndex &a_dit)
void getFluxStencil (VoFStencil &a_fluxStencil, const FaceIndex &a_face, const DataIndex &a_dit)
void getFaceCenteredFluxStencil (VoFStencil &a_fluxStencil, const FaceIndex &a_face, const DataIndex &a_dit)
void incrOpRegularDir (EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const int &a_dir, const DataIndex &a_datInd)
void applyOpIrregular (EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const DataIndex &a_datInd)

Protected Member Functions

void defineStencils ()

Protected Attributes

const IntVect m_ghostCellsPhi
const IntVect m_ghostCellsRHS
RefCountedPtr< EBQuadCFInterpm_quadCFIWithCoar
EBLevelGrid m_eblg
EBLevelGrid m_eblgFine
EBLevelGrid m_eblgCoar
EBLevelGrid m_eblgCoarMG
EBLevelGrid m_eblgCoarsenedFine
RefCountedPtr
< ConductivityBaseDomainBC
m_domainBC
RefCountedPtr
< ConductivityBaseEBBC
m_ebBC
Real m_dxFine
Real m_dx
Real m_dxCoar
RefCountedPtr< LevelData
< EBCellFAB > > 
m_acoef
RefCountedPtr< LevelData
< EBFluxFAB > > 
m_bcoef
RefCountedPtr< LevelData
< BaseIVFAB< Real > > > 
m_bcoIrreg
Real m_alpha
Real m_beta
int m_refToFine
int m_refToCoar
bool m_hasFine
bool m_hasInterpAve
bool m_hasCoar
EBMGAverage m_ebAverage
EBMGInterp m_ebInterp
LayoutData< RefCountedPtr
< EBStencil > > 
m_opEBStencil
LevelData< EBCellFABm_relCoef
LayoutData< VoFIteratorm_vofIterIrreg
LayoutData< VoFIteratorm_vofIterMulti
LayoutData< VoFIteratorm_vofIterDomLo [CH_SPACEDIM]
LayoutData< VoFIteratorm_vofIterDomHi [CH_SPACEDIM]
LayoutData< CFIVSm_loCFIVS [SpaceDim]
LayoutData< CFIVSm_hiCFIVS [SpaceDim]
EBFastFR m_fastFR
bool m_hasMGObjects
bool m_layoutChanged
EBMGAverage m_ebAverageMG
EBMGInterp m_ebInterpMG
DisjointBoxLayout m_dblCoarMG
EBISLayout m_ebislCoarMG
ProblemDomain m_domainCoarMG
Vector< IntVectm_colors

Static Protected Attributes

static bool s_turnOffBCs

Private Member Functions

void incrementFRCoar (EBFastFR &a_fluxReg, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi)
void incrementFRFine (EBFastFR &a_fluxReg, const LevelData< EBCellFAB > &a_phiFine, const LevelData< EBCellFAB > &a_phi, AMRLevelOp< LevelData< EBCellFAB > > *a_finerOp)
void getFlux (FArrayBox &a_flux, const FArrayBox &a_phi, const Box &a_faceBox, const int &a_idir, const Real &a_dx, const DataIndex &a_datInd)
void applyCFBCs (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *const a_phiCoarse, bool a_homogeneousCFBC)
void getOpVoFStencil (VoFStencil &a_stencil, const EBISBox &a_ebisbox, const VolIndex &a_vof)
void getOpVoFStencil (VoFStencil &a_stencil, const int &a_dir, const Vector< VolIndex > &a_allMonotoneVoFs, const EBISBox &a_ebisbox, const VolIndex &a_vof, const bool &a_lowOrder)
void getOpFaceStencil (VoFStencil &a_stencil, const Vector< VolIndex > &a_allMonotoneVofs, const EBISBox &a_ebisbox, const VolIndex &a_vof, int a_dir, const Side::LoHiSide &a_side, const FaceIndex &a_face, const bool &a_lowOrder)
void levelJacobi (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs)
void applyHomogeneousCFBCs (LevelData< EBCellFAB > &a_phi)
void applyHomogeneousCFBCs (EBCellFAB &a_phi, const DataIndex &a_datInd, int a_idir, Side::LoHiSide a_hiorlo)
 EBConductivityOp (const EBConductivityOp &a_opin)
void operator= (const EBConductivityOp &a_opin)

Constructor & Destructor Documentation

virtual EBConductivityOp::~EBConductivityOp (  )  [inline, virtual]

EBConductivityOp::EBConductivityOp (  )  [inline]

EBConductivityOp::EBConductivityOp ( const EBLevelGrid a_eblgFine,
const EBLevelGrid a_eblg,
const EBLevelGrid a_eblgCoar,
const EBLevelGrid a_eblgCoarMG,
const RefCountedPtr< EBQuadCFInterp > &  a_quadCFI,
const RefCountedPtr< ConductivityBaseDomainBC > &  a_domainBC,
const RefCountedPtr< ConductivityBaseEBBC > &  a_ebBC,
const Real a_dx,
const Real a_dxCoar,
const int &  a_refToFine,
const int &  a_refToCoar,
const bool &  a_hasFine,
const bool &  a_hasCoar,
const bool &  a_hasMGObjects,
const bool &  a_layoutChanged,
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 
)

If you are approaching this operator from this interface, consider backing away and using EBConductivityOpFactory to generate these objects. Really. a_eblgFine, : grid at finer level \ a_eblg, : grid at this level \ a_eblgCoar, : grid at coarser level \ a_eblgCoarMG, : grid at intermediate multigrid level \ a_domainBC, : domain boundary conditions at this level \ a_ebBC: eb boundary conditions at this level \ a_dx: grid spacing at this level \ a_origin: offset to lowest corner of the domain \ a_refToFine: refinement ratio to finer level \ a_refToCoar: refinement ratio to coarser level \ a_hasFiner: true if there is a finer AMR level, false otherwise. \ a_hasCoarser: true if there is a coarser AMR level. \ a_hasCoarserMG: true if there is a coarser MultiGrid level. \ a_preCondIters: number of iterations to do for pre-conditioning \ a_relaxType: 0 means point Jacobi, 1 is Gauss-Seidel. \ a_acoef: coefficent of identity \ a_bcoef: coefficient of gradient.\ a_ghostCellsPhi: Number of ghost cells in phi, correction\ a_ghostCellsRhs: Number of ghost cells in RHS, residual, lphi\ Ghost cell arguments are there for caching reasons. Once you set them, an error is thrown if you send in data that does not match.

EBConductivityOp::EBConductivityOp ( const EBConductivityOp a_opin  )  [inline, private]

References MayDay::Error().


Member Function Documentation

virtual void EBConductivityOp::setAlphaAndBeta ( const Real a_alpha,
const Real a_beta 
) [virtual]

set the constants in the equation

Implements TGAHelmOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::diagonalScale ( LevelData< EBCellFAB > &  a_rhs  )  [virtual]

Set the diagonal scaling of the operator. If you are solving rho(x) dphi/dt = L(phi), this would mean multiply by rho. In EB applications, even for constant coefficients, it means to multiply by kappa.

Implements TGAHelmOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::fillGrad ( const LevelData< EBCellFAB > &  a_phi  )  [inline, virtual]

virtual void EBConductivityOp::getFlux ( EBFluxFAB a_flux,
const LevelData< EBCellFAB > &  a_data,
const Box a_grid,
const DataIndex a_dit,
Real  a_scale 
) [virtual]

void EBConductivityOp::getFlux ( EBFaceFAB a_fluxCentroid,
const EBCellFAB a_phi,
const Box a_ghostedBox,
const Box a_fabBox,
const ProblemDomain a_domain,
const EBISBox a_ebisBox,
const Real a_dx,
const DataIndex a_datInd,
const int &  a_idir 
)

void EBConductivityOp::AMRResidualNC ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_rhs,
bool  a_homogeneousBC,
AMRLevelOp< LevelData< EBCellFAB > > *  a_finerOp 
) [virtual]

a_residual = a_rhs - L(a_phiFine, a_phi) no coaser AMR level

Implements AMRLevelOp< LevelData< EBCellFAB > >.

void EBConductivityOp::AMROperatorNC ( LevelData< EBCellFAB > &  a_LofPhi,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
bool  a_homogeneousBC,
AMRLevelOp< LevelData< EBCellFAB > > *  a_finerOp 
) [virtual]

apply AMR operator no coaser AMR level

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::residual ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  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< EBCellFAB > >.

virtual void EBConductivityOp::preCond ( LevelData< EBCellFAB > &  a_cor,
const LevelData< EBCellFAB > &  a_residual 
) [virtual]

Given the current state of the residual the correction, apply your preconditioner to a_cor.

Implements LinearOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::applyOp ( LevelData< EBCellFAB > &  a_opPhi,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > *const   a_phiCoarse,
const bool &  a_homogeneousPhysBC,
const bool &  a_homogeneousCFBC 
) [virtual]

This function assumes that coarse-fine boundary condtions have been dealt with.

virtual void EBConductivityOp::applyOpNoBoundary ( LevelData< EBCellFAB > &  a_opPhi,
const LevelData< EBCellFAB > &  a_phi 
) [virtual]

virtual function called by LevelTGA

Implements TGAHelmOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::applyOp ( LevelData< EBCellFAB > &  a_opPhi,
const LevelData< EBCellFAB > &  a_phi,
bool  a_homogeneousPhysBC 
) [virtual]

this is the linearop function. CFBC is set to homogeneous. phic is null

Implements LinearOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::create ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  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< EBCellFAB > >.

virtual void EBConductivityOp::createCoarsened ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs,
const int &  a_refRat 
) [virtual]

Real EBConductivityOp::AMRNorm ( const LevelData< EBCellFAB > &  a_coarResid,
const LevelData< EBCellFAB > &  a_fineResid,
const int &  a_refRat,
const int &  a_ord 
) [virtual]

virtual void EBConductivityOp::assign ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_rhs 
) [virtual]

Set a_lhs equal to a_rhs.

Implements LinearOp< LevelData< EBCellFAB > >.

virtual Real EBConductivityOp::dotProduct ( const LevelData< EBCellFAB > &  a_1,
const LevelData< EBCellFAB > &  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< EBCellFAB > >.

virtual void EBConductivityOp::incr ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_x,
Real  a_scale 
) [virtual]

Increment by scaled amount (a_lhs += a_scale*a_x).

Implements LinearOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::axby ( LevelData< EBCellFAB > &  a_lhs,
const LevelData< EBCellFAB > &  a_x,
const LevelData< EBCellFAB > &  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< EBCellFAB > >.

virtual void EBConductivityOp::scale ( LevelData< EBCellFAB > &  a_lhs,
const Real a_scale 
) [virtual]

Multiply the input by a given scale (a_lhs *= a_scale).

Implements LinearOp< LevelData< EBCellFAB > >.

virtual Real EBConductivityOp::norm ( const LevelData< EBCellFAB > &  a_rhs,
int  a_ord 
) [virtual]

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< EBCellFAB > >.

virtual void EBConductivityOp::setToZero ( LevelData< EBCellFAB > &  a_lhs  )  [virtual]

Set a_lhs to zero.

Implements LinearOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::setVal ( LevelData< EBCellFAB > &  a_lhs,
const Real a_value 
) [virtual]

virtual void EBConductivityOp::createCoarser ( LevelData< EBCellFAB > &  a_coarse,
const LevelData< EBCellFAB > &  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< EBCellFAB > >.

virtual void EBConductivityOp::relax ( LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_residual,
int  a_iterations 
) [virtual]

Use your relaxtion 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< EBCellFAB > >.

virtual void EBConductivityOp::restrictResidual ( LevelData< EBCellFAB > &  a_resCoarse,
LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_rhsFine 
) [virtual]

Calculate restricted residual: a_resCoarse[2h] = I[h->2h] (a_rhsFine[h] - L[h](a_phiFine[h]))

Implements MGLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::prolongIncrement ( LevelData< EBCellFAB > &  a_phiThisLevel,
const LevelData< EBCellFAB > &  a_correctCoarse 
) [virtual]

Correct the fine solution based on coarse correction: a_phiThisLevel += I[2h->h] (a_correctCoarse)

Implements MGLevelOp< LevelData< EBCellFAB > >.

virtual int EBConductivityOp::refToCoarser (  )  [virtual]

Refinement ratio between this level and coarser level. Returns 1 when there are no coarser AMRLevelOp objects

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual int EBConductivityOp::refToFiner (  )  [virtual]

Refinement ratio between this level and coarser level. Returns 1 when there are no coarser AMRLevelOp objects

virtual void EBConductivityOp::AMRResidual ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoarse,
const LevelData< EBCellFAB > &  a_rhs,
bool  a_homogeneousBC,
AMRLevelOp< LevelData< EBCellFAB > > *  a_finerOp 
) [virtual]

a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse)

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::AMRResidualNF ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoarse,
const LevelData< EBCellFAB > &  a_rhs,
bool  a_homogeneousBC 
) [virtual]

a_residual = a_rhs - L(a_phi, a_phiCoarse)

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::AMROperator ( LevelData< EBCellFAB > &  a_LofPhi,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoarse,
bool  a_homogeneousBC,
AMRLevelOp< LevelData< EBCellFAB > > *  a_finerOp 
) [virtual]

a_residual = a_rhs - L(a_phi, a_phiFine, a_phiCoarse)

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::AMROperatorNF ( LevelData< EBCellFAB > &  a_LofPhi,
const LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_phiCoarse,
bool  a_homogeneousBC 
) [virtual]

a_residual = a_rhs - L(a_phi, a_phiCoarse)

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::AMRRestrict ( LevelData< EBCellFAB > &  a_resCoarse,
const LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_coarseCorrection 
) [virtual]

a_resCoarse = I[h-2h] (a_residual - L(a_correction, a_coarseCorrection))

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::AMRProlong ( LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_coarseCorrection 
) [virtual]

a_correction += I[2h->h](a_coarseCorrection)

Implements AMRLevelOp< LevelData< EBCellFAB > >.

virtual void EBConductivityOp::AMRUpdateResidual ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_correction,
const LevelData< EBCellFAB > &  a_coarseCorrection 
) [virtual]

a_residual = a_residual - L(a_correction, a_coarseCorrection)

Implements AMRLevelOp< LevelData< EBCellFAB > >.

void EBConductivityOp::reflux ( LevelData< EBCellFAB > &  a_residual,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
AMRLevelOp< LevelData< EBCellFAB > > *  a_finerOp 
)

void EBConductivityOp::gsrbColor ( LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_lph,
const LevelData< EBCellFAB > &  a_rhs,
const IntVect a_color 
)

void EBConductivityOp::getDivFStencil ( VoFStencil a_vofStencil,
const VolIndex a_vof,
const DataIndex a_dit 
)

void EBConductivityOp::getFluxStencil ( VoFStencil a_fluxStencil,
const FaceIndex a_face,
const DataIndex a_dit 
)

void EBConductivityOp::getFaceCenteredFluxStencil ( VoFStencil a_fluxStencil,
const FaceIndex a_face,
const DataIndex a_dit 
)

void EBConductivityOp::incrOpRegularDir ( EBCellFAB a_lhs,
const EBCellFAB a_phi,
const bool &  a_homogeneous,
const int &  a_dir,
const DataIndex a_datInd 
)

void EBConductivityOp::applyOpIrregular ( EBCellFAB a_lhs,
const EBCellFAB a_phi,
const bool &  a_homogeneous,
const DataIndex a_datInd 
)

void EBConductivityOp::defineStencils (  )  [protected]

void EBConductivityOp::incrementFRCoar ( EBFastFR a_fluxReg,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi 
) [private]

void EBConductivityOp::incrementFRFine ( EBFastFR a_fluxReg,
const LevelData< EBCellFAB > &  a_phiFine,
const LevelData< EBCellFAB > &  a_phi,
AMRLevelOp< LevelData< EBCellFAB > > *  a_finerOp 
) [private]

void EBConductivityOp::getFlux ( FArrayBox a_flux,
const FArrayBox a_phi,
const Box a_faceBox,
const int &  a_idir,
const Real a_dx,
const DataIndex a_datInd 
) [private]

void EBConductivityOp::applyCFBCs ( LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > *const   a_phiCoarse,
bool  a_homogeneousCFBC 
) [private]

void EBConductivityOp::getOpVoFStencil ( VoFStencil a_stencil,
const EBISBox a_ebisbox,
const VolIndex a_vof 
) [private]

void EBConductivityOp::getOpVoFStencil ( VoFStencil a_stencil,
const int &  a_dir,
const Vector< VolIndex > &  a_allMonotoneVoFs,
const EBISBox a_ebisbox,
const VolIndex a_vof,
const bool &  a_lowOrder 
) [private]

void EBConductivityOp::getOpFaceStencil ( VoFStencil a_stencil,
const Vector< VolIndex > &  a_allMonotoneVofs,
const EBISBox a_ebisbox,
const VolIndex a_vof,
int  a_dir,
const Side::LoHiSide a_side,
const FaceIndex a_face,
const bool &  a_lowOrder 
) [private]

void EBConductivityOp::levelJacobi ( LevelData< EBCellFAB > &  a_phi,
const LevelData< EBCellFAB > &  a_rhs 
) [private]

void EBConductivityOp::applyHomogeneousCFBCs ( LevelData< EBCellFAB > &  a_phi  )  [private]

void EBConductivityOp::applyHomogeneousCFBCs ( EBCellFAB a_phi,
const DataIndex a_datInd,
int  a_idir,
Side::LoHiSide  a_hiorlo 
) [private]

void EBConductivityOp::operator= ( const EBConductivityOp a_opin  )  [inline, private]

References MayDay::Error().


Member Data Documentation

bool EBConductivityOp::s_turnOffBCs [static, protected]

bool EBConductivityOp::m_hasFine [protected]

bool EBConductivityOp::m_hasCoar [protected]


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

Generated on Tue Apr 14 14:23:12 2009 for Chombo + EB by  doxygen 1.5.5