NWOViscousTensorOp Class Reference

#include <NWOViscousTensorOp.H>

Inheritance diagram for NWOViscousTensorOp:

Inheritance graph
[legend]

List of all members.


Detailed Description

operator is alpha a I + (divF) = alpha*acoef I + beta*div(eta(grad B + grad BT) + lambda I div B ) beta, lambda, and eta are incorporated into the flux F

AMRLevelOp functions

static int s_prolongType
 prolongation type
RefCountedPtr< LevelData
< FluxBox > > 
m_eta
RefCountedPtr< LevelData
< FluxBox > > 
m_lambda
RefCountedPtr< LevelData
< FArrayBox > > 
m_acoef
Real m_alpha
Real m_beta
int m_refToCoar
int m_refToFine
BCHolder m_bc
DisjointBoxLayout m_grids
DisjointBoxLayout m_gridsCoar
DisjointBoxLayout m_gridsFine
int m_ncomp
Real m_dx
Real m_dxCrse
ProblemDomain m_domain
Real m_safety
Real m_relaxTolerance
int m_relaxMinIter
LevelDataOps< FArrayBoxm_levelOps
Copier m_exchangeCopier
NWOQuadCFInterp m_interpWithCoarser
Vector< IntVectm_colors
LevelData< FArrayBoxm_relaxCoef
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 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 AMRProlong (LevelData< FArrayBox > &a_correction, const LevelData< FArrayBox > &a_coarseCorrection)
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 outputLevel (LevelData< FArrayBox > &a_rhs, string &a_name)
virtual void outputAMR (Vector< LevelData< FArrayBox > * > &a_rhs, string &a_name)
void homogeneousCFInterp (LevelData< FArrayBox > &a_phif)
void homogeneousCFInterpPhi (LevelData< FArrayBox > &a_phif, const DataIndex &a_datInd, int a_idir, Side::LoHiSide a_hiorlo)
void homogeneousCFInterpTanGrad (LevelData< FArrayBox > &a_tanGrad, const LevelData< FArrayBox > &a_phi, const DataIndex &a_datInd, int a_idir, Side::LoHiSide a_hiorlo)
void interpOnIVSHomo (LevelData< FArrayBox > &a_phif, const DataIndex &a_datInd, const int a_idir, const Side::LoHiSide a_hiorlo, const IntVectSet &a_interpIVS)
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)
void AMROperatorNF (LevelData< FArrayBox > &a_LofPhi, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_phiCoarse, bool a_homogeneousBC)
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)
void cfinterp (const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phiCoarse)
RefCountedPtr< LevelData
< FluxBox > > 
getEta () const
 access function
RefCountedPtr< LevelData
< FluxBox > > 
getLambda () const
RefCountedPtr< LevelData
< FArrayBox > > 
getACoef () const
Real getAlpha () const
Real getBeta () const
void defineRelCoef ()
 NWOViscousTensorOp ()
 weak construction is bad
 NWOViscousTensorOp (const NWOViscousTensorOp &a_opin)
void operator= (const NWOViscousTensorOp &a_opin)

Public Types

enum  prolongationType { piecewiseConstant = 0, linearInterp, NUM_INTERP_TYPE }

Public Member Functions

virtual void diagonalScale (LevelData< FArrayBox > &a_rhs, bool a_kappaWeighted)
virtual void setAlphaAndBeta (const Real &a_alpha, const Real &a_beta)
virtual void diagonalScale (LevelData< FArrayBox > &a_rhs)
virtual void divideByIdentityCoef (LevelData< FArrayBox > &a_rhs)
virtual ~NWOViscousTensorOp ()
 NWOViscousTensorOp (const DisjointBoxLayout &a_grids, const DisjointBoxLayout &a_gridsFine, const DisjointBoxLayout &a_gridsCoar, const RefCountedPtr< LevelData< FluxBox > > &a_eta, const RefCountedPtr< LevelData< FluxBox > > &a_lambda, const RefCountedPtr< LevelData< FArrayBox > > &a_acoef, Real a_alpha, Real a_beta, int a_refToFine, int a_refToCoar, const ProblemDomain &a_domain, const Real &a_dxLevel, const Real &a_dxCoar, BCFunc a_bc, Real a_safety=VTOP_DEFAULT_SAFETY, Real a_relaxTolerance=0.0, int a_relaxMinIter=1000)
 NWOViscousTensorOp (const DisjointBoxLayout &a_grids, const DisjointBoxLayout &a_gridsFine, const DisjointBoxLayout &a_gridsCoar, const RefCountedPtr< LevelData< FluxBox > > &a_eta, const RefCountedPtr< LevelData< FluxBox > > &a_lambda, const RefCountedPtr< LevelData< FArrayBox > > &a_acoef, Real a_alpha, Real a_beta, int a_refToFine, int a_refToCoar, const ProblemDomain &a_domain, const Real &a_dxLevel, const Real &a_dxCoar, BCHolder &a_bc, Real a_safety=VTOP_DEFAULT_SAFETY, Real a_relaxTolerance=0.0, int a_relaxMinIter=1000)
void fillGrad (const LevelData< FArrayBox > &a_phiFine)
void cellGrad (FArrayBox &a_gradPhi, const FArrayBox &a_phi, const Box &a_grid)
void define (const DisjointBoxLayout &a_grids, const DisjointBoxLayout &a_gridsFine, const DisjointBoxLayout &a_gridsCoar, const RefCountedPtr< LevelData< FluxBox > > &a_eta, const RefCountedPtr< LevelData< FluxBox > > &a_lambda, const RefCountedPtr< LevelData< FArrayBox > > &a_acoef, Real a_alpha, Real a_beta, int a_refToFine, int a_refToCoar, const ProblemDomain &a_domain, const Real &a_dxLevel, const Real &a_dxCoar, BCHolder &a_bc, Real a_safety=VTOP_DEFAULT_SAFETY, Real a_relaxTolerance=0.0, int a_relaxMinIter=1000)
virtual void residual (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs, bool a_homogeneous=false)
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 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)
void reflux (const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, LevelData< FArrayBox > &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)
void applyOpNoBoundary (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi)
void getFlux (FArrayBox &a_flux, const FArrayBox &a_data, const FArrayBox &a_etaFace, const FArrayBox &a_lambdaFace, const Box &a_facebox, int a_dir, int ref=1)
void computeOperatorNoBCs (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi)
 utility function which computes operator after all bc's have been set
virtual void assign (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_rhs)
virtual Real dotProduct (const LevelData< FArrayBox > &a_1, const LevelData< FArrayBox > &a_2)
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, Real 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 dx () const
virtual Real dxCrse () const
virtual void setToZero (LevelData< FArrayBox > &a_x)
MGLevelOp functions
virtual void relax (LevelData< FArrayBox > &a_e, const LevelData< FArrayBox > &a_residual, int iterations)
virtual void createCoarser (LevelData< FArrayBox > &a_coarse, const LevelData< FArrayBox > &a_fine, bool 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)

Static Public Member Functions

static void loHiCenterFace (Box &a_loBox, int &a_hasLo, Box &a_hiBox, int &a_hasHi, Box &a_centerBox, const ProblemDomain &a_eblg, const Box &a_inBox, const int &a_dir)
static void getFaceDivAndGrad (FArrayBox &a_faceDiv, FArrayBox &a_faceGrad, const FArrayBox &a_data, const FArrayBox &a_gradData, const ProblemDomain &a_domain, const Box &a_faceBox, const int &a_faceDir, const Real a_dx)
static void getFlux (FArrayBox &a_flux, const FArrayBox &a_data, const FArrayBox &a_etaFace, const FArrayBox &a_lamFace, const Box &a_faceBox, const ProblemDomain &a_domain, const Real &a_dx, const Real &a_beta, int a_dir, int a_ref)

Member Enumeration Documentation

Enumerator:
piecewiseConstant 
linearInterp 
NUM_INTERP_TYPE 


Constructor & Destructor Documentation

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

NWOViscousTensorOp::NWOViscousTensorOp ( const DisjointBoxLayout a_grids,
const DisjointBoxLayout a_gridsFine,
const DisjointBoxLayout a_gridsCoar,
const RefCountedPtr< LevelData< FluxBox > > &  a_eta,
const RefCountedPtr< LevelData< FluxBox > > &  a_lambda,
const RefCountedPtr< LevelData< FArrayBox > > &  a_acoef,
Real  a_alpha,
Real  a_beta,
int  a_refToFine,
int  a_refToCoar,
const ProblemDomain a_domain,
const Real a_dxLevel,
const Real a_dxCoar,
BCFunc  a_bc,
Real  a_safety = VTOP_DEFAULT_SAFETY,
Real  a_relaxTolerance = 0.0,
int  a_relaxMinIter = 1000 
)

NWOViscousTensorOp::NWOViscousTensorOp ( const DisjointBoxLayout a_grids,
const DisjointBoxLayout a_gridsFine,
const DisjointBoxLayout a_gridsCoar,
const RefCountedPtr< LevelData< FluxBox > > &  a_eta,
const RefCountedPtr< LevelData< FluxBox > > &  a_lambda,
const RefCountedPtr< LevelData< FArrayBox > > &  a_acoef,
Real  a_alpha,
Real  a_beta,
int  a_refToFine,
int  a_refToCoar,
const ProblemDomain a_domain,
const Real a_dxLevel,
const Real a_dxCoar,
BCHolder a_bc,
Real  a_safety = VTOP_DEFAULT_SAFETY,
Real  a_relaxTolerance = 0.0,
int  a_relaxMinIter = 1000 
)

NWOViscousTensorOp::NWOViscousTensorOp (  )  [inline, private]

weak construction is bad

References MayDay::Error().

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

References MayDay::Error().


Member Function Documentation

virtual void NWOViscousTensorOp::diagonalScale ( LevelData< FArrayBox > &  a_rhs,
bool  a_kappaWeighted 
) [inline, 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.

Parameters:
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< FArrayBox > >.

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

Sets the scaling constants in the Helmholtz equation.

Parameters:
a_alpha The scaling constant that multiplies the identity term.
a_beta The scaling constant that multiplies the derivative term.

Implements TGAHelmOp< LevelData< FArrayBox > >.

References defineRelCoef(), m_alpha, and m_beta.

virtual void NWOViscousTensorOp::diagonalScale ( LevelData< FArrayBox > &  a_rhs  )  [inline, virtual]

Scales the right hand side of the Helmholtz equation by the identity term in the operator. This version assumes volume fraction weighting.

Parameters:
a_rhs The right hand side of the equation to be scaled.

Reimplemented from TGAHelmOp< LevelData< FArrayBox > >.

References BoxLayout::dataIterator(), LevelData< T >::disjointBoxLayout(), m_acoef, LayoutIterator::ok(), and SpaceDim.

virtual void NWOViscousTensorOp::divideByIdentityCoef ( LevelData< FArrayBox > &  a_rhs  )  [inline, virtual]

Divides the right hand side of the Helmholtz equation by the identity coefficient rho(x) in the equation rho(x) dphi/dt = L(phi).

Parameters:
a_rhs The right hand side of the equation to be scaled.

Implements TGAHelmOp< LevelData< FArrayBox > >.

References BoxLayout::dataIterator(), LevelData< T >::disjointBoxLayout(), m_acoef, LayoutIterator::ok(), and SpaceDim.

void NWOViscousTensorOp::fillGrad ( const LevelData< FArrayBox > &  a_phiFine  )  [inline, virtual]

Fillgrad here is a placeholder that does not do anything. The whole fillGrad paradigm is an artifact of the how old ViscousTensorOp worked. Sometimes performance tuning eventually forces you to do the right thing.

Implements LevelTGAHelmOp< LevelData< FArrayBox >, FluxBox >.

static void NWOViscousTensorOp::loHiCenterFace ( Box a_loBox,
int &  a_hasLo,
Box a_hiBox,
int &  a_hasHi,
Box a_centerBox,
const ProblemDomain a_eblg,
const Box a_inBox,
const int &  a_dir 
) [static]

void NWOViscousTensorOp::cellGrad ( FArrayBox a_gradPhi,
const FArrayBox a_phi,
const Box a_grid 
)

static void NWOViscousTensorOp::getFaceDivAndGrad ( FArrayBox a_faceDiv,
FArrayBox a_faceGrad,
const FArrayBox a_data,
const FArrayBox a_gradData,
const ProblemDomain a_domain,
const Box a_faceBox,
const int &  a_faceDir,
const Real  a_dx 
) [static]

void NWOViscousTensorOp::define ( const DisjointBoxLayout a_grids,
const DisjointBoxLayout a_gridsFine,
const DisjointBoxLayout a_gridsCoar,
const RefCountedPtr< LevelData< FluxBox > > &  a_eta,
const RefCountedPtr< LevelData< FluxBox > > &  a_lambda,
const RefCountedPtr< LevelData< FArrayBox > > &  a_acoef,
Real  a_alpha,
Real  a_beta,
int  a_refToFine,
int  a_refToCoar,
const ProblemDomain a_domain,
const Real a_dxLevel,
const Real a_dxCoar,
BCHolder a_bc,
Real  a_safety = VTOP_DEFAULT_SAFETY,
Real  a_relaxTolerance = 0.0,
int  a_relaxMinIter = 1000 
)

virtual void NWOViscousTensorOp::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 NWOViscousTensorOp::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 > >.

virtual void NWOViscousTensorOp::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 NWOViscousTensorOp::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 NWOViscousTensorOp::createCoarsened ( LevelData< FArrayBox > &  a_lhs,
const LevelData< FArrayBox > &  a_rhs,
const int &  a_refRat 
) [virtual]

void NWOViscousTensorOp::reflux ( const LevelData< FArrayBox > &  a_phiFine,
const LevelData< FArrayBox > &  a_phi,
LevelData< FArrayBox > &  residual,
AMRLevelOp< LevelData< FArrayBox > > *  a_finerOp 
)

virtual void NWOViscousTensorOp::getFlux ( FluxBox a_flux,
const LevelData< FArrayBox > &  a_data,
const Box a_grid,
const DataIndex a_dit,
Real  a_scale 
) [inline, virtual]

void NWOViscousTensorOp::applyOpNoBoundary ( LevelData< FArrayBox > &  a_ans,
const LevelData< FArrayBox > &  a_phi 
) [inline, virtual]

Apply the differential operator without any boundary or coarse-fine boundary conditions and no finer level

Parameters:
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 > >.

References computeOperatorNoBCs().

void NWOViscousTensorOp::getFlux ( FArrayBox a_flux,
const FArrayBox a_data,
const FArrayBox a_etaFace,
const FArrayBox a_lambdaFace,
const Box a_facebox,
int  a_dir,
int  ref = 1 
)

static void NWOViscousTensorOp::getFlux ( FArrayBox a_flux,
const FArrayBox a_data,
const FArrayBox a_etaFace,
const FArrayBox a_lamFace,
const Box a_faceBox,
const ProblemDomain a_domain,
const Real a_dx,
const Real a_beta,
int  a_dir,
int  a_ref 
) [static]

void NWOViscousTensorOp::computeOperatorNoBCs ( LevelData< FArrayBox > &  a_lhs,
const LevelData< FArrayBox > &  a_phi 
)

utility function which computes operator after all bc's have been set

Referenced by applyOpNoBoundary().

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

Set a_lhs equal to a_rhs.

Implements LinearOp< LevelData< FArrayBox > >.

virtual Real NWOViscousTensorOp::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 NWOViscousTensorOp::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 NWOViscousTensorOp::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 > >.

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

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

Implements LinearOp< LevelData< FArrayBox > >.

virtual Real NWOViscousTensorOp::norm ( const LevelData< FArrayBox > &  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< FArrayBox > >.

virtual Real NWOViscousTensorOp::dx (  )  const [inline, virtual]

Return dx at this level of refinement

Reimplemented from LinearOp< LevelData< FArrayBox > >.

References m_dx.

virtual Real NWOViscousTensorOp::dxCrse (  )  const [inline, virtual]

References m_dxCrse.

virtual void NWOViscousTensorOp::setToZero ( LevelData< FArrayBox > &  a_lhs  )  [virtual]

Set a_lhs to zero.

Implements LinearOp< LevelData< FArrayBox > >.

virtual void NWOViscousTensorOp::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 NWOViscousTensorOp::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 NWOViscousTensorOp::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 > >.

virtual void NWOViscousTensorOp::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 NWOViscousTensorOp::refToCoarser (  )  [inline, virtual]

returns 1 when there are no coarser AMRLevelOp objects

Implements AMRLevelOp< LevelData< FArrayBox > >.

References m_refToCoar.

virtual void NWOViscousTensorOp::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 NWOViscousTensorOp::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 NWOViscousTensorOp::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 NWOViscousTensorOp::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 NWOViscousTensorOp::AMRProlong ( LevelData< FArrayBox > &  a_correction,
const LevelData< FArrayBox > &  a_coarseCorrection 
) [virtual]

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

Implements AMRLevelOp< LevelData< FArrayBox > >.

virtual void NWOViscousTensorOp::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 NWOViscousTensorOp::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 > >.

virtual void NWOViscousTensorOp::outputLevel ( LevelData< FArrayBox > &  a_rhs,
string &  a_name 
) [inline, virtual]

Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.

References writeLevelname().

Referenced by outputAMR().

virtual void NWOViscousTensorOp::outputAMR ( Vector< LevelData< FArrayBox > * > &  a_rhs,
string &  a_name 
) [inline, virtual]

Reimplemented from AMRLevelOp< LevelData< FArrayBox > >.

References outputLevel().

void NWOViscousTensorOp::homogeneousCFInterp ( LevelData< FArrayBox > &  a_phif  ) 

void NWOViscousTensorOp::homogeneousCFInterpPhi ( LevelData< FArrayBox > &  a_phif,
const DataIndex a_datInd,
int  a_idir,
Side::LoHiSide  a_hiorlo 
)

does homogeneous coarse/fine interpolation for phi

void NWOViscousTensorOp::homogeneousCFInterpTanGrad ( LevelData< FArrayBox > &  a_tanGrad,
const LevelData< FArrayBox > &  a_phi,
const DataIndex a_datInd,
int  a_idir,
Side::LoHiSide  a_hiorlo 
)

does homogeneous coarse/fine interpolation for tangential gradient (needs phi ghost cells to be filled in first, so should call homogeneousCFInterpPhi first)

void NWOViscousTensorOp::interpOnIVSHomo ( LevelData< FArrayBox > &  a_phif,
const DataIndex a_datInd,
const int  a_idir,
const Side::LoHiSide  a_hiorlo,
const IntVectSet a_interpIVS 
)

void NWOViscousTensorOp::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 > >.

void NWOViscousTensorOp::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 NWOViscousTensorOp::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 > >.

void NWOViscousTensorOp::cfinterp ( const LevelData< FArrayBox > &  a_phiFine,
const LevelData< FArrayBox > &  a_phiCoarse 
)

RefCountedPtr<LevelData<FluxBox> > NWOViscousTensorOp::getEta (  )  const [inline]

access function

References m_eta.

RefCountedPtr<LevelData<FluxBox> > NWOViscousTensorOp::getLambda (  )  const [inline]

References m_lambda.

RefCountedPtr<LevelData<FArrayBox> > NWOViscousTensorOp::getACoef (  )  const [inline]

References m_acoef.

Real NWOViscousTensorOp::getAlpha (  )  const [inline]

References m_alpha.

Real NWOViscousTensorOp::getBeta (  )  const [inline]

References m_beta.

void NWOViscousTensorOp::defineRelCoef (  )  [protected]

Referenced by setAlphaAndBeta().

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

References MayDay::Error().


Member Data Documentation

prolongation type

make this public/static to make it easy to change, while also ensuring that it remains constent across all instances

Referenced by getEta().

Referenced by getLambda().

Referenced by getAlpha(), and setAlphaAndBeta().

Referenced by getBeta(), and setAlphaAndBeta().

Referenced by refToCoarser().

int NWOViscousTensorOp::m_ncomp [protected]

Referenced by dx().

Referenced by dxCrse().

Referenced by getFlux().


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

Generated on Fri Apr 5 04:25:14 2019 for Chombo + EB by  doxygen 1.5.5