VCAMRPoissonOp Class Reference

#include <VCAMRPoissonOp.H>

Inheritance diagram for VCAMRPoissonOp:

Inheritance graph
[legend]
Collaboration diagram for VCAMRPoissonOp:

Collaboration graph
[legend]

List of all members.


Detailed Description

Operator for solving variable-coefficient (alpha * aCoef(x) * I - beta * Div(bCoef(x) . Grad)) phi = rho over an AMR hierarchy.

Public Member Functions

void setAlphaAndBeta (const Real &a_alpha, const Real &a_beta)
 for tga stuff
void diagonalScale (LevelData< FArrayBox > &a_rhs)
 for tga stuff
void setCoefs (const RefCountedPtr< LevelData< FArrayBox > > &a_aCoef, const RefCountedPtr< LevelData< FluxBox > > &a_bCoef, const Real &a_alpha, const Real &a_beta)
 also calls reset lambda
void resetLambda ()
 should be called every time coefs are changed.
 VCAMRPoissonOp ()
virtual ~VCAMRPoissonOp ()
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, const Real &a_dxLevel, int a_refRatio, const ProblemDomain &a_domain, BCHolder a_bc)
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)
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)
void computeLambda ()
 Compute lambda once alpha, aCoef, beta, bCoef are defined.
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 applyOpNoBoundary (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi)
 apply operator without any boundary or coarse-fine boundary conditions and no finer level
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 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_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 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 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_homogeneousPhysBC, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp)
virtual void AMROperatorNC (LevelData< FArrayBox > &a_LofPhi, const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, bool a_homogeneousPhysBC, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp)
virtual void AMROperatorNF (LevelData< FArrayBox > &a_LofPhi, const LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_phiCoarse, 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)
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 write (const LevelData< FArrayBox > *a_data, const char *a_filename)

Public Attributes

Real m_alpha
 Identity operator constant coefficient -- if you change this call resetLambda().
RefCountedPtr< LevelData
< FArrayBox > > 
m_aCoef
 Identity operator spatially varying coefficient storage (cell-centered) --- if you change this call resetLambda().
Real m_beta
 Laplacian operator constant coefficient --- if you change this call resetLambda().
RefCountedPtr< LevelData
< FluxBox > > 
m_bCoef
 Laplacian operator spatially varying coefficient storage (face-centered) --- if you change this call resetLambda().
LevelData< FArrayBoxm_lambda
 Reciprocal of the diagonal entry of the operator matrix.
Real m_dxCrse
Real m_dx
ProblemDomain m_domain

Protected Member Functions

void levelGSRB (LevelData< FArrayBox > &a_e, const LevelData< FArrayBox > &a_residual)
void levelJacobi (LevelData< FArrayBox > &a_phi, const LevelData< FArrayBox > &a_rhs)
void homogeneousCFInterp (LevelData< FArrayBox > &a_phif)
void homogeneousCFInterp (LevelData< FArrayBox > &a_phif, 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 getFlux (FArrayBox &a_flux, const FArrayBox &a_data, const FluxBox &a_bCoef, const Box &a_facebox, int a_dir, int a_ref=1) const
 computes flux over face-centered a_facebox.
void reflux (const LevelData< FArrayBox > &a_phiFine, const LevelData< FArrayBox > &a_phi, LevelData< FArrayBox > &a_residual, AMRLevelOp< LevelData< FArrayBox > > *a_finerOp)
void computeOperatorNoBCs (LevelData< FArrayBox > &a_lhs, const LevelData< FArrayBox > &a_phi)
 utility function which computes operator after all bc's have been set

Protected Attributes

LevelDataOps< FArrayBoxm_levelOps
BCHolder m_bc
LayoutData< CFIVSm_loCFIVS [SpaceDim]
LayoutData< CFIVSm_hiCFIVS [SpaceDim]
Copier m_exchangeCopier
QuadCFInterp m_interpWithCoarser
int m_refToCoarser
int m_refToFiner

Constructor & Destructor Documentation

VCAMRPoissonOp::VCAMRPoissonOp (  )  [inline]

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


Member Function Documentation

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

for tga stuff

Implements TGAHelmOp< LevelData< FArrayBox > >.

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

void VCAMRPoissonOp::setCoefs ( const RefCountedPtr< LevelData< FArrayBox > > &  a_aCoef,
const RefCountedPtr< LevelData< FluxBox > > &  a_bCoef,
const Real a_alpha,
const Real a_beta 
)

also calls reset lambda

void VCAMRPoissonOp::resetLambda (  ) 

should be called every time coefs are changed.

void VCAMRPoissonOp::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 VCAMRPoissonOp::define ( const DisjointBoxLayout a_grids,
const DisjointBoxLayout a_baseBAPtr,
const Real a_dxLevel,
int  a_refRatio,
const ProblemDomain a_domain,
BCHolder  a_bc 
)

define function for AMRLevelOp which has no finer AMR level

void VCAMRPoissonOp::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 
)

full define function for AMRLevelOp with both coarser and finer levels

void VCAMRPoissonOp::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 
)

full define function for AMRLevelOp with finer levels, but no coarser

void VCAMRPoissonOp::computeLambda (  ) 

Compute lambda once alpha, aCoef, beta, bCoef are defined.

virtual void VCAMRPoissonOp::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 VCAMRPoissonOp::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 VCAMRPoissonOp::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 VCAMRPoissonOp::applyOpNoBoundary ( LevelData< FArrayBox > &  a_ans,
const LevelData< FArrayBox > &  a_phi 
) [virtual]

apply operator without any boundary or coarse-fine boundary conditions and no finer level

Implements TGAHelmOp< LevelData< FArrayBox > >.

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

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

Set a_lhs equal to a_rhs.

Implements LinearOp< LevelData< FArrayBox > >.

virtual Real VCAMRPoissonOp::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 VCAMRPoissonOp::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 VCAMRPoissonOp::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 VCAMRPoissonOp::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 VCAMRPoissonOp::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 void VCAMRPoissonOp::setToZero ( LevelData< FArrayBox > &  a_lhs  )  [virtual]

Set a_lhs to zero.

Implements LinearOp< LevelData< FArrayBox > >.

virtual void VCAMRPoissonOp::relax ( LevelData< FArrayBox > &  a_correction,
const LevelData< FArrayBox > &  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< FArrayBox > >.

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

returns 1 when there are no coarser AMRLevelOp objects

Implements AMRLevelOp< LevelData< FArrayBox > >.

References m_refToCoarser.

virtual void VCAMRPoissonOp::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 VCAMRPoissonOp::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 VCAMRPoissonOp::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 VCAMRPoissonOp::AMROperator ( LevelData< FArrayBox > &  a_LofPhi,
const LevelData< FArrayBox > &  a_phiFine,
const LevelData< FArrayBox > &  a_phi,
const LevelData< FArrayBox > &  a_phiCoarse,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< FArrayBox > > *  a_finerOp 
) [virtual]

apply AMR operator, including coarse-fine matching

Implements AMRLevelOp< LevelData< FArrayBox > >.

virtual void VCAMRPoissonOp::AMROperatorNC ( LevelData< FArrayBox > &  a_LofPhi,
const LevelData< FArrayBox > &  a_phiFine,
const LevelData< FArrayBox > &  a_phi,
bool  a_homogeneousPhysBC,
AMRLevelOp< LevelData< FArrayBox > > *  a_finerOp 
) [virtual]

AMR operator assuming no more coarser AMR levels

Implements AMRLevelOp< LevelData< FArrayBox > >.

virtual void VCAMRPoissonOp::AMROperatorNF ( LevelData< FArrayBox > &  a_LofPhi,
const LevelData< FArrayBox > &  a_phi,
const LevelData< FArrayBox > &  a_phiCoarse,
bool  a_homogeneousPhysBC 
) [virtual]

AMR operator assuming no finer level

Implements AMRLevelOp< LevelData< FArrayBox > >.

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

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

Implements AMRLevelOp< LevelData< FArrayBox > >.

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

void VCAMRPoissonOp::levelGSRB ( LevelData< FArrayBox > &  a_e,
const LevelData< FArrayBox > &  a_residual 
) [protected]

void VCAMRPoissonOp::levelJacobi ( LevelData< FArrayBox > &  a_phi,
const LevelData< FArrayBox > &  a_rhs 
) [protected]

void VCAMRPoissonOp::homogeneousCFInterp ( LevelData< FArrayBox > &  a_phif  )  [protected]

void VCAMRPoissonOp::homogeneousCFInterp ( LevelData< FArrayBox > &  a_phif,
const DataIndex a_datInd,
int  a_idir,
Side::LoHiSide  a_hiorlo 
) [protected]

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

void VCAMRPoissonOp::getFlux ( FArrayBox a_flux,
const FArrayBox a_data,
const FluxBox a_bCoef,
const Box a_facebox,
int  a_dir,
int  a_ref = 1 
) const [protected]

computes flux over face-centered a_facebox.

void VCAMRPoissonOp::reflux ( const LevelData< FArrayBox > &  a_phiFine,
const LevelData< FArrayBox > &  a_phi,
LevelData< FArrayBox > &  a_residual,
AMRLevelOp< LevelData< FArrayBox > > *  a_finerOp 
) [protected]

void VCAMRPoissonOp::computeOperatorNoBCs ( LevelData< FArrayBox > &  a_lhs,
const LevelData< FArrayBox > &  a_phi 
) [protected]

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


Member Data Documentation

Identity operator constant coefficient -- if you change this call resetLambda().

Identity operator spatially varying coefficient storage (cell-centered) --- if you change this call resetLambda().

Laplacian operator constant coefficient --- if you change this call resetLambda().

Laplacian operator spatially varying coefficient storage (face-centered) --- if you change this call resetLambda().

Reciprocal of the diagonal entry of the operator matrix.

Referenced by refToCoarser().


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

Generated on Tue Apr 14 14:24:02 2009 for Chombo + EB by  doxygen 1.5.5