Chombo + EB  3.0
Public Member Functions | Protected Attributes | List of all members
EBLevelTGA Class Reference

#include <EBLevelTGA.H>

Inheritance diagram for EBLevelTGA:
Inheritance graph
[legend]

Public Member Functions

 EBLevelTGA (const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > &a_solver)
 
virtual ~EBLevelTGA ()
 Destructor. More...
 
void computeDiffusion (LevelData< EBCellFAB > &a_DiffusiveTerm, LevelData< EBCellFAB > &a_phiOld, LevelData< EBCellFAB > &a_src, EBFluxRegister *a_FineFluxRegPtr, EBFluxRegister *a_crseFluxRegPtr, const LevelData< EBCellFAB > *a_crsePhiOldPtr, const LevelData< EBCellFAB > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, bool a_rhsAlreadyKappaWeighted=false)
 
void updateSoln (LevelData< EBCellFAB > &a_phiNew, LevelData< EBCellFAB > &a_phiOld, LevelData< EBCellFAB > &a_src, EBFluxRegister *a_fineFluxRegPtr, EBFluxRegister *a_crseFluxRegPtr, const LevelData< EBCellFAB > *a_crsePhiOldPtr, const LevelData< EBCellFAB > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, bool a_rhsAlreadyKappaWeighted=false, int a_fluxStartComponent=0)
 
void setSourceGhostCells (LevelData< EBCellFAB > &a_src, const DisjointBoxLayout &a_grids, int a_lev)
 
void setEBLG (Vector< EBLevelGrid > &a_eblg)
 
- Public Member Functions inherited from BaseLevelTGA< LevelData< EBCellFAB >, EBFluxFAB, EBFluxRegister >
 BaseLevelTGA (const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > &a_solver)
 
virtual ~BaseLevelTGA ()
 Destructor, called after destructors of BaseLevelTGA subclasses. More...
 
void updateSoln (LevelData< EBCellFAB > &a_phiNew, LevelData< EBCellFAB > &a_phiOld, LevelData< EBCellFAB > &a_src, LevelData< EBFluxFAB > &a_flux, EBFluxRegister *a_fineFluxRegPtr, EBFluxRegister *a_crseFluxRegPtr, const LevelData< EBCellFAB > *a_crsePhiOldPtr, const LevelData< EBCellFAB > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, bool a_rhsAlreadyKappaWeighted=false, int a_fluxStartComponent=0)
 
void computeDiffusion (LevelData< EBCellFAB > &a_diffusiveTerm, LevelData< EBCellFAB > &a_phiOld, LevelData< EBCellFAB > &a_src, LevelData< EBFluxFAB > &a_flux, EBFluxRegister *a_fineFluxRegPtr, EBFluxRegister *a_crseFluxRegPtr, const LevelData< EBCellFAB > *a_crsePhiOldPtr, const LevelData< EBCellFAB > *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, bool a_rhsAlreadyKappaWeighted=false)
 
- Public Member Functions inherited from BaseLevelHeatSolver< LevelData< EBCellFAB >, EBFluxFAB, EBFluxRegister >
 BaseLevelHeatSolver (const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > &a_solver)
 
virtual ~BaseLevelHeatSolver ()
 Destructor, called after destructors of BaseLevelHeatSolver subclasses. More...
 
virtual void applyOperator (LevelData< EBCellFAB > &a_ans, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *a_phiC, int a_level, Real a_alpha, Real a_beta, bool a_applyBC)
 
void applyHelm (LevelData< EBCellFAB > &a_ans, const LevelData< EBCellFAB > &a_phi, const LevelData< EBCellFAB > *a_phiC, int a_level, Real a_mu, Real a_dt, bool a_homogeneousBC)
 
void incrementFlux (LevelData< EBFluxFAB > &a_diffusiveFlux, LevelData< EBCellFAB > &a_phi, int a_level, Real a_mu, Real a_dt, Real a_sign, bool a_setToZero)
 
void solveHelm (LevelData< EBCellFAB > &a_phi, LevelData< EBCellFAB > &a_phiC, LevelData< EBCellFAB > &a_rhs, int a_level, Real a_mu, Real a_dt, bool a_zeroPhi=true)
 
void resetSolverAlphaAndBeta (const Real &a_alpha, const Real &a_beta)
 
virtual LevelTGAHelmOp< LevelData< EBCellFAB >, EBFluxFAB > * newOp (const ProblemDomain &a_indexSpace, RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &a_opFact)
 
int size () const
 Returns the number of grid levels on which this integrator operates. More...
 

Protected Attributes

bool m_isEBLGSet
 
Vector< EBLevelGridm_eblg
 
- Protected Attributes inherited from BaseLevelTGA< LevelData< EBCellFAB >, EBFluxFAB, EBFluxRegister >
Real m_mu1
 
Real m_mu2
 
Real m_mu3
 
Real m_mu4
 
Real m_r1
 
- Protected Attributes inherited from BaseLevelHeatSolver< LevelData< EBCellFAB >, EBFluxFAB, EBFluxRegister >
Vector< DisjointBoxLayoutm_grids
 The disjoint box layouts at every AMR grid level. More...
 
Vector< int > m_refRat
 The refinement ratios between AMR grid levels. More...
 
ProblemDomain m_level0Domain
 The coarsest domain on which the Helmholtz equation is integrated. More...
 
Vector< LevelTGAHelmOp< LevelData< EBCellFAB >, EBFluxFAB > * > m_ops
 
RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > m_solver
 The multigrid solver used to solve the Helmholtz equation. More...
 

Additional Inherited Members

- Protected Member Functions inherited from BaseLevelHeatSolver< LevelData< EBCellFAB >, EBFluxFAB, EBFluxRegister >
void timeInterp (LevelData< EBCellFAB > &a_data, const LevelData< EBCellFAB > &a_oldData, const LevelData< EBCellFAB > &a_newData, Real a_time, Real a_oldTime, Real a_newTime, int a_level)
 

Detailed Description

This class implements an imlplicit integrator for a diffusion equation using an approach devised by Twizzell, Gumel, and Arigu.

Constructor & Destructor Documentation

◆ EBLevelTGA()

EBLevelTGA::EBLevelTGA ( const Vector< DisjointBoxLayout > &  a_grids,
const Vector< int > &  a_refRat,
const ProblemDomain a_level0Domain,
RefCountedPtr< AMRLevelOpFactory< LevelData< EBCellFAB > > > &  a_opFact,
const RefCountedPtr< AMRMultiGrid< LevelData< EBCellFAB > > > &  a_solver 
)
inline

Create a new TGA level integrator.

Parameters
a_gridsThe disjoint box layout on which the level integrator is defined.
a_refRatThe refinement ratios for the boxes.
a_level0DomainThe coarsest grid level defining the problem domain.
a_opFactA factory that generates level diffusion operators.
a_solverA multigrid solver.

References m_isEBLGSet.

◆ ~EBLevelTGA()

virtual EBLevelTGA::~EBLevelTGA ( )
inlinevirtual

Destructor.

References computeDiffusion(), setSourceGhostCells(), and updateSoln().

Member Function Documentation

◆ computeDiffusion()

void EBLevelTGA::computeDiffusion ( LevelData< EBCellFAB > &  a_DiffusiveTerm,
LevelData< EBCellFAB > &  a_phiOld,
LevelData< EBCellFAB > &  a_src,
EBFluxRegister a_FineFluxRegPtr,
EBFluxRegister a_crseFluxRegPtr,
const LevelData< EBCellFAB > *  a_crsePhiOldPtr,
const LevelData< EBCellFAB > *  a_crsePhiNewPtr,
Real  a_oldTime,
Real  a_crseOldTime,
Real  a_crseNewTime,
Real  a_dt,
int  a_level,
bool  a_zeroPhi = true,
bool  a_rhsAlreadyKappaWeighted = false 
)

Computes the time-centered diffusive term for explicit updating. This integrates the diffusion equation for a new value of a quantity phi and then computes the diffusion term using a finite difference stencil in time.

Parameters
a_diffusionTermThe term computed by this method.
a_phiOldThe value of phi at the beginning of the time step, at the current grid level.
a_srcThe source term in the diffusion equation at the current grid level.
a_fineFluxRegPtrA pointer to the flux register at the finer adjacent grid level (or NULL if there is none).
a_crseFluxRegPtrA pointer to the flux register at the coarser adjacent grid level (or NULL if there is none).
a_crsePhiOldPtrA pointer to the LevelData object holding the old value of phi on the coarser adjacent grid level (or NULL if there is none).
a_crsePhiNewPtrA pointer to the LevelData object holding the updated value of phi on the coarser adjacent grid level (or NULL if there is none).
a_oldTimeThe time at the beginning of the integration step on the current grid level.
a_crseOldTimeThe time at the beginning of the integration step on the coarser adjacent grid level.
a_crseNewTimeThe time at the end of the integration step on the coarser adjacent grid level.
a_dtThe time step on the current grid level.
a_levelThe current grid level.
a_zeroPhiIf this flag is true, the initial estimate of phi will be set to zero. Otherwise, a_phiOld will be used.

Referenced by EBLevelBackwardEuler::~EBLevelBackwardEuler(), EBLevelCrankNicolson::~EBLevelCrankNicolson(), and ~EBLevelTGA().

◆ updateSoln()

void EBLevelTGA::updateSoln ( LevelData< EBCellFAB > &  a_phiNew,
LevelData< EBCellFAB > &  a_phiOld,
LevelData< EBCellFAB > &  a_src,
EBFluxRegister a_fineFluxRegPtr,
EBFluxRegister a_crseFluxRegPtr,
const LevelData< EBCellFAB > *  a_crsePhiOldPtr,
const LevelData< EBCellFAB > *  a_crsePhiNewPtr,
Real  a_oldTime,
Real  a_crseOldTime,
Real  a_crseNewTime,
Real  a_dt,
int  a_level,
bool  a_zeroPhi = true,
bool  a_rhsAlreadyKappaWeighted = false,
int  a_fluxStartComponent = 0 
)
virtual

Integrates the diffusion equation, storing the result in a_phiNew.

Parameters
a_phiNewThe new value of phi at the current grid level.
a_phiOldThe value of phi at the beginning of the time step, at the current grid level.
a_srcThe source term in the diffusion equation at the current grid level.
a_fineFluxRegPtrA pointer to the flux register at the finer adjacent grid level (or NULL if there is none).
a_crseFluxRegPtrA pointer to the flux register at the coarser adjacent grid level (or NULL if there is none).
a_crsePhiOldPtrA pointer to the LevelData object holding the old value of phi on the coarser adjacent grid level (or NULL if there is none).
a_crsePhiNewPtrA pointer to the LevelData object holding the updated value of phi on the coarser adjacent grid level (or NULL if there is none).
a_oldTimeThe time at the beginning of the integration step on the current grid level.
a_crseOldTimeThe time at the beginning of the integration step on the coarser adjacent grid level.
a_crseNewTimeThe time at the end of the integration step on the coarser adjacent grid level.
a_dtThe time step on the current grid level.
a_levelThe current grid level.
a_zeroPhiIf this flag is true, the initial estimate of phi will be set to zero. Otherwise, a_phiOld will be used.

Reimplemented from BaseLevelHeatSolver< LevelData< EBCellFAB >, EBFluxFAB, EBFluxRegister >.

Referenced by EBLevelBackwardEuler::~EBLevelBackwardEuler(), EBLevelCrankNicolson::~EBLevelCrankNicolson(), and ~EBLevelTGA().

◆ setSourceGhostCells()

void EBLevelTGA::setSourceGhostCells ( LevelData< EBCellFAB > &  a_src,
const DisjointBoxLayout a_dbl,
int  a_level 
)
virtual

Sets the value of the source term in the Helmholtz equation on ghost cells. This method should be overridden by subclasses of BaseLevelTGA.

Parameters
a_srcThe source term in the Helmholtz equation whose ghost cell values are to be set.
a_dblThe disjoint box layout that indexes a_src.
a_levelThe grid level at which the ghost cell values are to be set.

Implements BaseLevelTGA< LevelData< EBCellFAB >, EBFluxFAB, EBFluxRegister >.

Referenced by ~EBLevelTGA().

◆ setEBLG()

void EBLevelTGA::setEBLG ( Vector< EBLevelGrid > &  a_eblg)
inline

References m_eblg, and m_isEBLGSet.

Member Data Documentation

◆ m_isEBLGSet

bool EBLevelTGA::m_isEBLGSet
protected

◆ m_eblg

Vector<EBLevelGrid> EBLevelTGA::m_eblg
protected

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