#include <NWOEBConductivityOp.H>
Public Member Functions | |
NWOEBConductivityOp (const EBLevelGrid &a_eblgFine, const EBLevelGrid &a_eblg, const EBLevelGrid &a_eblgCoar, const EBLevelGrid &a_eblgCoarMG, const RefCountedPtr< NWOEBQuadCFInterp > &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, const int &a_relaxType) | |
~NWOEBConductivityOp () | |
virtual void | resetACoefficient (RefCountedPtr< LevelData< EBCellFAB > > &a_acoef) |
Real | dx () const |
virtual void | kappaScale (LevelData< EBCellFAB > &a_rhs) |
for eb only. kappa weight the rhs but do not multiply by the identity coefficient | |
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) |
void | setAlphaAndBeta (const Real &a_alpha, const Real &a_beta) |
void | diagonalScale (LevelData< EBCellFAB > &a_rhs, bool a_kappaWeighted) |
void | divideByIdentityCoef (LevelData< EBCellFAB > &a_rhs) |
void | fillGrad (const LevelData< EBCellFAB > &a_phi) |
These functions are part of the LevelTGA interface...... | |
void | fillPhiGhost (const EBCellFAB &a_phi, const DataIndex &a_datInd, bool a_homog) const |
void | fillPhiGhost (const LevelData< EBCellFAB > &a_phi, bool a_homog) const |
void | getFlux (EBFluxFAB &a_flux, const LevelData< EBCellFAB > &a_data, const Box &a_grid, const DataIndex &a_dit, Real a_scale) |
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 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 | |
void | applyOpIrregular (EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const DataIndex &a_datInd) |
void | applyOpRegular (EBCellFAB &a_lhs, const EBCellFAB &a_phi, const bool &a_homogeneous, const DataIndex &a_datInd) |
virtual void | applyOp (LevelData< EBCellFAB > &a_opPhi, const LevelData< EBCellFAB > &a_phi, bool a_homogeneousPhysBC=false) |
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 Real | localMaxNorm (const LevelData< EBCellFAB > &a_rhs) |
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, bool a_skip_res=false) |
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 | gsrbColor (LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > &a_rhs, const IntVect &a_color) |
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 | homogeneousCFInterp (LevelData< EBCellFAB > &a_phif) |
const EBLevelGrid & | getEBLG () const |
void | defineStencils () |
do not call this one unless you really know what you are doing | |
void | getAlphaDiagWeight (LayoutData< BaseIVFAB< Real > > const *&a_alphaDiagWeight) |
void | getAlphaBeta (Real &a_alpha, Real &a_beta) |
const RefCountedPtr< BaseDomainBC > | getDomainBC () |
const RefCountedPtr< LevelData < EBCellFAB > > | getAScalingCoefficients () |
const RefCountedPtr< LevelData < EBFluxFAB > > | getBScalingCoefficients () |
void | getEBBCFluxStencil (LayoutData< BaseIVFAB< VoFStencil > > const *&a_ebbcFluxStencil) |
Static Public Member Functions | |
static void | setForceNoEBCF (bool a_forceNoEBCF) |
static void | doLazyRelax (bool a_doLazyRelax) |
Protected Member Functions | |
Real | getSafety () |
void | incrOpRegularAllDirs (Box *a_loBox, Box *a_hiBox, int *a_hasLo, int *a_hasHi, Box &a_curDblBox, Box &a_curPhiBox, int a_nComps, BaseFab< Real > &a_curOpPhiFAB, const BaseFab< Real > &a_curPhiFAB, bool a_homogeneousPhysBC, const DataIndex &a_dit) |
void | applyDomainFlux (Box *a_loBox, Box *a_hiBox, int *a_hasLo, int *a_hasHi, Box &a_dblBox, int a_nComps, BaseFab< Real > &a_phiFAB, bool a_homogeneousPhysBC, const DataIndex &a_dit) |
void | GSColorAllIrregular (EBCellFAB &a_phi, const EBCellFAB &a_rhs, const int &a_icolor, const DataIndex &a_dit) |
virtual void | calculateAlphaWeight () |
virtual void | calculateRelaxationCoefficient () |
void | defineEBCFStencils () |
void | getFluxEBCF (EBFaceFAB &a_flux, const EBCellFAB &a_phi, const Box &a_ghostedBox, Vector< FaceIndex > &a_faceitEBCF, Vector< VoFStencil > &a_stenEBCF) |
void | getFluxRegOnly (EBFaceFAB &a_fluxCentroid, const EBCellFAB &a_phi, const Box &a_ghostedBox, const Real &a_dx, const DataIndex &a_datInd, const int &a_idir) |
Protected Attributes | |
LayoutData< Vector< FaceIndex > > | m_faceitCoar [2 *SpaceDim] |
LayoutData< Vector< VoFStencil > > | m_stencilCoar [2 *SpaceDim] |
int | m_relaxType |
const IntVect | m_ghostPhi |
const IntVect | m_ghostRHS |
RefCountedPtr< NWOEBQuadCFInterp > | m_interpWithCoarser |
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 |
LayoutData< BaseIVFAB< Real > > | m_alphaDiagWeight |
LayoutData< BaseIVFAB< Real > > | m_betaDiagWeight |
int | m_refToFine |
int | m_refToCoar |
bool | m_hasEBCF |
bool | m_hasFine |
bool | m_hasInterpAve |
bool | m_hasCoar |
EBMGAverage | m_ebAverage |
EBMGInterp | m_ebInterp |
LayoutData< RefCountedPtr < VCAggStencil > > | m_opEBStencil |
LevelData< EBCellFAB > | m_relCoef |
Multigrid relaxation coefficient. | |
LayoutData< VoFIterator > | m_vofIterIrreg |
LayoutData< VoFIterator > | m_vofIterMulti |
LayoutData< VoFIterator > | m_vofIterDomLo [CH_SPACEDIM] |
LayoutData< VoFIterator > | m_vofIterDomHi [CH_SPACEDIM] |
LayoutData< CFIVS > | m_loCFIVS [SpaceDim] |
LayoutData< CFIVS > | m_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< IntVect > | m_colors |
Copier | m_exchangeCopier |
Static Protected Attributes | |
static bool | s_turnOffBCs |
static bool | s_forceNoEBCF |
static bool | s_doLazyRelax |
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 | 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) |
NWOEBConductivityOp () | |
Default constructor. Creates an undefined conductivity operator. | |
NWOEBConductivityOp (const NWOEBConductivityOp &a_opin) | |
void | operator= (const NWOEBConductivityOp &a_opin) |
NWOEBConductivityOp::NWOEBConductivityOp | ( | const EBLevelGrid & | a_eblgFine, | |
const EBLevelGrid & | a_eblg, | |||
const EBLevelGrid & | a_eblgCoar, | |||
const EBLevelGrid & | a_eblgCoarMG, | |||
const RefCountedPtr< NWOEBQuadCFInterp > & | 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, | |||
const int & | a_relaxType | |||
) |
Constructs a conductivity operator using the given data. This constructor is for time-independent a and b coefficients. If you are approaching this operator from this interface, consider backing away and using NWOEBConductivityOpFactory to generate these objects. Really. 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.
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 |
NWOEBConductivityOp::~NWOEBConductivityOp | ( | ) | [inline] |
NWOEBConductivityOp::NWOEBConductivityOp | ( | ) | [private] |
Default constructor. Creates an undefined conductivity operator.
NWOEBConductivityOp::NWOEBConductivityOp | ( | const NWOEBConductivityOp & | a_opin | ) | [inline, private] |
References MayDay::Error().
virtual void NWOEBConductivityOp::resetACoefficient | ( | RefCountedPtr< LevelData< EBCellFAB > > & | a_acoef | ) | [inline, virtual] |
This sets the data storage for the a coefficient to a different object and recalculates the stuff it depends on. Use this only if you know what you're doing.
References calculateAlphaWeight(), calculateRelaxationCoefficient(), and m_acoef.
Real NWOEBConductivityOp::dx | ( | ) | const [inline, virtual] |
Return dx at this level of refinement
Reimplemented from LinearOp< LevelData< EBCellFAB > >.
References m_dx.
for eb only. kappa weight the rhs but do not multiply by the identity coefficient
Reimplemented from TGAHelmOp< LevelData< EBCellFAB > >.
void NWOEBConductivityOp::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 NWOEBConductivityOp::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 > >.
Sets the scaling constants in the Helmholtz equation.
a_alpha | The scaling constant that multiplies the identity term. | |
a_beta | The scaling constant that multiplies the derivative term. |
Implements TGAHelmOp< LevelData< EBCellFAB > >.
void NWOEBConductivityOp::diagonalScale | ( | LevelData< EBCellFAB > & | a_rhs, | |
bool | a_kappaWeighted | |||
) | [virtual] |
Scales the right hand side of the Helmholtz equation by the identity term in the operator. If you are solving rho(x) dphi/dt = L(phi), this means multiply by rho (or kappa * rho in embedded boundary calculations.
a_rhs | The right hand side of the equation to be scaled. | |
a_kappaWeighted | If set to true, a_rhs will be scaled by the volume fraction in addition to the identity term. |
Reimplemented from TGAHelmOp< LevelData< EBCellFAB > >.
Divides the right hand side of the Helmholtz equation by the identity coefficient rho(x) in the equation rho(x) dphi/dt = L(phi).
a_rhs | The right hand side of the equation to be scaled. |
Implements TGAHelmOp< LevelData< EBCellFAB > >.
These functions are part of the LevelTGA interface......
Implements LevelTGAHelmOp< LevelData< EBCellFAB >, EBFluxFAB >.
void NWOEBConductivityOp::fillPhiGhost | ( | const EBCellFAB & | a_phi, | |
const DataIndex & | a_datInd, | |||
bool | a_homog | |||
) | const |
void NWOEBConductivityOp::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 void NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::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.
Referenced by applyOpNoBoundary().
virtual void NWOEBConductivityOp::applyOpNoBoundary | ( | LevelData< EBCellFAB > & | a_opPhi, | |
const LevelData< EBCellFAB > & | a_phi | |||
) | [inline, virtual] |
virtual function called by LevelTGA
Implements TGAHelmOp< LevelData< EBCellFAB > >.
References applyOp(), and s_turnOffBCs.
void NWOEBConductivityOp::applyOpIrregular | ( | EBCellFAB & | a_lhs, | |
const EBCellFAB & | a_phi, | |||
const bool & | a_homogeneous, | |||
const DataIndex & | a_datInd | |||
) |
void NWOEBConductivityOp::applyOpRegular | ( | EBCellFAB & | a_lhs, | |
const EBCellFAB & | a_phi, | |||
const bool & | a_homogeneous, | |||
const DataIndex & | a_datInd | |||
) |
virtual void NWOEBConductivityOp::applyOp | ( | LevelData< EBCellFAB > & | a_opPhi, | |
const LevelData< EBCellFAB > & | a_phi, | |||
bool | a_homogeneousPhysBC = false | |||
) | [virtual] |
this is the linearop function. CFBC is set to homogeneous. phic is null
Implements LinearOp< LevelData< EBCellFAB > >.
virtual void NWOEBConductivityOp::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 NWOEBConductivityOp::createCoarsened | ( | LevelData< EBCellFAB > & | a_lhs, | |
const LevelData< EBCellFAB > & | a_rhs, | |||
const int & | a_refRat | |||
) | [virtual] |
Implements AMRLevelOp< LevelData< EBCellFAB > >.
Real NWOEBConductivityOp::AMRNorm | ( | const LevelData< EBCellFAB > & | a_coarResid, | |
const LevelData< EBCellFAB > & | a_fineResid, | |||
const int & | a_refRat, | |||
const int & | a_ord | |||
) | [virtual] |
Reimplemented from AMRLevelOp< LevelData< EBCellFAB > >.
virtual void NWOEBConductivityOp::assign | ( | LevelData< EBCellFAB > & | a_lhs, | |
const LevelData< EBCellFAB > & | a_rhs | |||
) | [virtual] |
Set a_lhs equal to a_rhs.
Implements LinearOp< LevelData< EBCellFAB > >.
virtual Real NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::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 > >.
Reimplemented from AMRLevelOp< LevelData< EBCellFAB > >.
Set a_lhs to zero.
Implements LinearOp< LevelData< EBCellFAB > >.
virtual void NWOEBConductivityOp::setVal | ( | LevelData< EBCellFAB > & | a_lhs, | |
const Real & | a_value | |||
) | [virtual] |
virtual void NWOEBConductivityOp::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 NWOEBConductivityOp::relax | ( | LevelData< EBCellFAB > & | a_correction, | |
const LevelData< EBCellFAB > & | 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< EBCellFAB > >.
virtual void NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::refToFiner | ( | ) | [virtual] |
Refinement ratio between this level and coarser level. Returns 1 when there are no coarser AMRLevelOp objects
virtual void NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::AMRRestrict | ( | LevelData< EBCellFAB > & | a_resCoarse, | |
const LevelData< EBCellFAB > & | a_residual, | |||
const LevelData< EBCellFAB > & | a_correction, | |||
const LevelData< EBCellFAB > & | a_coarseCorrection, | |||
bool | a_skip_res = false | |||
) | [virtual] |
a_resCoarse = I[h-2h] (a_residual - L(a_correction, a_coarseCorrection))
Implements AMRLevelOp< LevelData< EBCellFAB > >.
virtual void NWOEBConductivityOp::AMRProlong | ( | LevelData< EBCellFAB > & | a_correction, | |
const LevelData< EBCellFAB > & | a_coarseCorrection | |||
) | [virtual] |
a_correction += I[2h->h](a_coarseCorrection)
Implements AMRLevelOp< LevelData< EBCellFAB > >.
virtual void NWOEBConductivityOp::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 NWOEBConductivityOp::gsrbColor | ( | LevelData< EBCellFAB > & | a_phi, | |
const LevelData< EBCellFAB > & | a_rhs, | |||
const IntVect & | a_color | |||
) |
void NWOEBConductivityOp::reflux | ( | LevelData< EBCellFAB > & | a_residual, | |
const LevelData< EBCellFAB > & | a_phiFine, | |||
const LevelData< EBCellFAB > & | a_phi, | |||
AMRLevelOp< LevelData< EBCellFAB > > * | a_finerOp | |||
) |
void NWOEBConductivityOp::gsrbColor | ( | LevelData< EBCellFAB > & | a_phi, | |
const LevelData< EBCellFAB > & | a_lph, | |||
const LevelData< EBCellFAB > & | a_rhs, | |||
const IntVect & | a_color | |||
) |
void NWOEBConductivityOp::getDivFStencil | ( | VoFStencil & | a_vofStencil, | |
const VolIndex & | a_vof, | |||
const DataIndex & | a_dit | |||
) |
void NWOEBConductivityOp::getFluxStencil | ( | VoFStencil & | a_fluxStencil, | |
const FaceIndex & | a_face, | |||
const DataIndex & | a_dit | |||
) |
void NWOEBConductivityOp::getFaceCenteredFluxStencil | ( | VoFStencil & | a_fluxStencil, | |
const FaceIndex & | a_face, | |||
const DataIndex & | a_dit | |||
) |
const EBLevelGrid& NWOEBConductivityOp::getEBLG | ( | ) | const [inline] |
References m_eblg.
static void NWOEBConductivityOp::setForceNoEBCF | ( | bool | a_forceNoEBCF | ) | [inline, static] |
References s_forceNoEBCF.
static void NWOEBConductivityOp::doLazyRelax | ( | bool | a_doLazyRelax | ) | [inline, static] |
References s_doLazyRelax.
void NWOEBConductivityOp::defineStencils | ( | ) |
do not call this one unless you really know what you are doing
void NWOEBConductivityOp::getAlphaDiagWeight | ( | LayoutData< BaseIVFAB< Real > > const *& | a_alphaDiagWeight | ) | [inline] |
References m_alphaDiagWeight.
const RefCountedPtr<BaseDomainBC> NWOEBConductivityOp::getDomainBC | ( | ) | [inline] |
References m_domainBC.
const RefCountedPtr<LevelData<EBCellFAB> > NWOEBConductivityOp::getAScalingCoefficients | ( | ) | [inline] |
References m_acoef.
const RefCountedPtr<LevelData<EBFluxFAB> > NWOEBConductivityOp::getBScalingCoefficients | ( | ) | [inline] |
References m_bcoef.
void NWOEBConductivityOp::getEBBCFluxStencil | ( | LayoutData< BaseIVFAB< VoFStencil > > const *& | a_ebbcFluxStencil | ) | [inline] |
References m_ebBC.
Real NWOEBConductivityOp::getSafety | ( | ) | [protected] |
void NWOEBConductivityOp::incrOpRegularAllDirs | ( | Box * | a_loBox, | |
Box * | a_hiBox, | |||
int * | a_hasLo, | |||
int * | a_hasHi, | |||
Box & | a_curDblBox, | |||
Box & | a_curPhiBox, | |||
int | a_nComps, | |||
BaseFab< Real > & | a_curOpPhiFAB, | |||
const BaseFab< Real > & | a_curPhiFAB, | |||
bool | a_homogeneousPhysBC, | |||
const DataIndex & | a_dit | |||
) | [protected] |
void NWOEBConductivityOp::applyDomainFlux | ( | Box * | a_loBox, | |
Box * | a_hiBox, | |||
int * | a_hasLo, | |||
int * | a_hasHi, | |||
Box & | a_dblBox, | |||
int | a_nComps, | |||
BaseFab< Real > & | a_phiFAB, | |||
bool | a_homogeneousPhysBC, | |||
const DataIndex & | a_dit | |||
) | [protected] |
void NWOEBConductivityOp::GSColorAllIrregular | ( | EBCellFAB & | a_phi, | |
const EBCellFAB & | a_rhs, | |||
const int & | a_icolor, | |||
const DataIndex & | a_dit | |||
) | [protected] |
virtual void NWOEBConductivityOp::calculateAlphaWeight | ( | ) | [protected, virtual] |
Referenced by resetACoefficient().
virtual void NWOEBConductivityOp::calculateRelaxationCoefficient | ( | ) | [protected, virtual] |
Referenced by resetACoefficient().
void NWOEBConductivityOp::defineEBCFStencils | ( | ) | [protected] |
void NWOEBConductivityOp::getFluxEBCF | ( | EBFaceFAB & | a_flux, | |
const EBCellFAB & | a_phi, | |||
const Box & | a_ghostedBox, | |||
Vector< FaceIndex > & | a_faceitEBCF, | |||
Vector< VoFStencil > & | a_stenEBCF | |||
) | [protected] |
void NWOEBConductivityOp::getFluxRegOnly | ( | EBFaceFAB & | a_fluxCentroid, | |
const EBCellFAB & | a_phi, | |||
const Box & | a_ghostedBox, | |||
const Real & | a_dx, | |||
const DataIndex & | a_datInd, | |||
const int & | a_idir | |||
) | [protected] |
void NWOEBConductivityOp::incrementFRCoar | ( | EBFastFR & | a_fluxReg, | |
const LevelData< EBCellFAB > & | a_phiFine, | |||
const LevelData< EBCellFAB > & | a_phi | |||
) | [private] |
void NWOEBConductivityOp::incrementFRFine | ( | EBFastFR & | a_fluxReg, | |
const LevelData< EBCellFAB > & | a_phiFine, | |||
const LevelData< EBCellFAB > & | a_phi, | |||
AMRLevelOp< LevelData< EBCellFAB > > * | a_finerOp | |||
) | [private] |
void NWOEBConductivityOp::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 NWOEBConductivityOp::getOpVoFStencil | ( | VoFStencil & | a_stencil, | |
const EBISBox & | a_ebisbox, | |||
const VolIndex & | a_vof | |||
) | [private] |
void NWOEBConductivityOp::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 NWOEBConductivityOp::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 NWOEBConductivityOp::operator= | ( | const NWOEBConductivityOp & | a_opin | ) | [inline, private] |
References MayDay::Error().
bool NWOEBConductivityOp::s_turnOffBCs [static, protected] |
Referenced by applyOpNoBoundary().
bool NWOEBConductivityOp::s_forceNoEBCF [static, protected] |
Referenced by setForceNoEBCF().
bool NWOEBConductivityOp::s_doLazyRelax [static, protected] |
Referenced by doLazyRelax().
LayoutData< Vector<FaceIndex> > NWOEBConductivityOp::m_faceitCoar[2 *SpaceDim] [protected] |
LayoutData< Vector<VoFStencil> > NWOEBConductivityOp::m_stencilCoar[2 *SpaceDim] [protected] |
int NWOEBConductivityOp::m_relaxType [protected] |
const IntVect NWOEBConductivityOp::m_ghostPhi [protected] |
const IntVect NWOEBConductivityOp::m_ghostRHS [protected] |
EBLevelGrid NWOEBConductivityOp::m_eblg [protected] |
Referenced by getEBLG().
EBLevelGrid NWOEBConductivityOp::m_eblgFine [protected] |
EBLevelGrid NWOEBConductivityOp::m_eblgCoar [protected] |
EBLevelGrid NWOEBConductivityOp::m_eblgCoarMG [protected] |
EBLevelGrid NWOEBConductivityOp::m_eblgCoarsenedFine [protected] |
Referenced by getDomainBC().
RefCountedPtr<ConductivityBaseEBBC> NWOEBConductivityOp::m_ebBC [protected] |
Referenced by getEBBCFluxStencil().
Real NWOEBConductivityOp::m_dxFine [protected] |
Real NWOEBConductivityOp::m_dx [protected] |
Referenced by dx().
Real NWOEBConductivityOp::m_dxCoar [protected] |
RefCountedPtr<LevelData<EBCellFAB> > NWOEBConductivityOp::m_acoef [protected] |
Referenced by getAScalingCoefficients(), and resetACoefficient().
RefCountedPtr<LevelData<EBFluxFAB> > NWOEBConductivityOp::m_bcoef [protected] |
Referenced by getBScalingCoefficients().
RefCountedPtr<LevelData<BaseIVFAB<Real> > > NWOEBConductivityOp::m_bcoIrreg [protected] |
Real NWOEBConductivityOp::m_alpha [protected] |
Referenced by getAlphaBeta().
Real NWOEBConductivityOp::m_beta [protected] |
Referenced by getAlphaBeta().
LayoutData<BaseIVFAB<Real> > NWOEBConductivityOp::m_alphaDiagWeight [protected] |
Referenced by getAlphaDiagWeight().
LayoutData<BaseIVFAB<Real> > NWOEBConductivityOp::m_betaDiagWeight [protected] |
int NWOEBConductivityOp::m_refToFine [protected] |
int NWOEBConductivityOp::m_refToCoar [protected] |
bool NWOEBConductivityOp::m_hasEBCF [protected] |
bool NWOEBConductivityOp::m_hasFine [protected] |
bool NWOEBConductivityOp::m_hasInterpAve [protected] |
bool NWOEBConductivityOp::m_hasCoar [protected] |
EBMGAverage NWOEBConductivityOp::m_ebAverage [protected] |
EBMGInterp NWOEBConductivityOp::m_ebInterp [protected] |
LayoutData<RefCountedPtr<VCAggStencil> > NWOEBConductivityOp::m_opEBStencil [protected] |
LevelData<EBCellFAB> NWOEBConductivityOp::m_relCoef [protected] |
Multigrid relaxation coefficient.
LayoutData<VoFIterator > NWOEBConductivityOp::m_vofIterIrreg [protected] |
LayoutData<VoFIterator > NWOEBConductivityOp::m_vofIterMulti [protected] |
LayoutData<VoFIterator > NWOEBConductivityOp::m_vofIterDomLo[CH_SPACEDIM] [protected] |
LayoutData<VoFIterator > NWOEBConductivityOp::m_vofIterDomHi[CH_SPACEDIM] [protected] |
LayoutData<CFIVS> NWOEBConductivityOp::m_loCFIVS[SpaceDim] [protected] |
LayoutData<CFIVS> NWOEBConductivityOp::m_hiCFIVS[SpaceDim] [protected] |
EBFastFR NWOEBConductivityOp::m_fastFR [protected] |
bool NWOEBConductivityOp::m_hasMGObjects [protected] |
bool NWOEBConductivityOp::m_layoutChanged [protected] |
EBMGAverage NWOEBConductivityOp::m_ebAverageMG [protected] |
EBMGInterp NWOEBConductivityOp::m_ebInterpMG [protected] |
DisjointBoxLayout NWOEBConductivityOp::m_dblCoarMG [protected] |
EBISLayout NWOEBConductivityOp::m_ebislCoarMG [protected] |
ProblemDomain NWOEBConductivityOp::m_domainCoarMG [protected] |
Vector<IntVect> NWOEBConductivityOp::m_colors [protected] |
Copier NWOEBConductivityOp::m_exchangeCopier [protected] |