Chombo + EB + MF  3.2
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType > Class Template Referenceabstract

#include <BaseLevelTGA.H>

Inheritance diagram for BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >:
Inheritance graph
[legend]

Public Member Functions

 BaseLevelTGA (const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelDataType > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelDataType > > &a_solver)
 
virtual ~BaseLevelTGA ()
 Destructor, called after destructors of BaseLevelTGA subclasses. More...
 
void updateSoln (LevelDataType &a_phiNew, LevelDataType &a_phiOld, LevelDataType &a_src, LevelData< FluxDataType > &a_flux, FluxRegisterType *a_fineFluxRegPtr, FluxRegisterType *a_crseFluxRegPtr, const LevelDataType *a_crsePhiOldPtr, const LevelDataType *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 (LevelDataType &a_diffusiveTerm, LevelDataType &a_phiOld, LevelDataType &a_src, LevelData< FluxDataType > &a_flux, FluxRegisterType *a_fineFluxRegPtr, FluxRegisterType *a_crseFluxRegPtr, const LevelDataType *a_crsePhiOldPtr, const LevelDataType *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< LevelDataType, FluxDataType, FluxRegisterType >
 BaseLevelHeatSolver (const Vector< DisjointBoxLayout > &a_grids, const Vector< int > &a_refRat, const ProblemDomain &a_level0Domain, RefCountedPtr< AMRLevelOpFactory< LevelDataType > > &a_opFact, const RefCountedPtr< AMRMultiGrid< LevelDataType > > &a_solver)
 
virtual ~BaseLevelHeatSolver ()
 Destructor, called after destructors of BaseLevelHeatSolver subclasses. More...
 
virtual void updateSoln (LevelDataType &a_phiNew, LevelDataType &a_phiOld, LevelDataType &a_src, FluxRegisterType *a_fineFluxRegPtr, FluxRegisterType *a_crseFluxRegPtr, const LevelDataType *a_crsePhiOldPtr, const LevelDataType *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 void applyOperator (LevelDataType &a_ans, const LevelDataType &a_phi, const LevelDataType *a_phiC, int a_level, Real a_alpha, Real a_beta, bool a_applyBC)
 
void applyHelm (LevelDataType &a_ans, const LevelDataType &a_phi, const LevelDataType *a_phiC, int a_level, Real a_mu, Real a_dt, bool a_homogeneousBC)
 
void incrementFlux (LevelData< FluxDataType > &a_diffusiveFlux, LevelDataType &a_phi, int a_level, Real a_mu, Real a_dt, Real a_sign, bool a_setToZero)
 
void solveHelm (LevelDataType &a_phi, LevelDataType &a_phiC, LevelDataType &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< LevelDataType, FluxDataType > * newOp (const ProblemDomain &a_indexSpace, RefCountedPtr< AMRLevelOpFactory< LevelDataType > > &a_opFact)
 
int size () const
 Returns the number of grid levels on which this integrator operates. More...
 

Protected Member Functions

virtual void setSourceGhostCells (LevelDataType &a_src, const DisjointBoxLayout &a_dbl, int a_level)=0
 
- Protected Member Functions inherited from BaseLevelHeatSolver< LevelDataType, FluxDataType, FluxRegisterType >
void timeInterp (LevelDataType &a_data, const LevelDataType &a_oldData, const LevelDataType &a_newData, Real a_time, Real a_oldTime, Real a_newTime, int a_level)
 

Protected Attributes

Real m_mu1
 
Real m_mu2
 
Real m_mu3
 
Real m_mu4
 
Real m_r1
 
- Protected Attributes inherited from BaseLevelHeatSolver< LevelDataType, FluxDataType, FluxRegisterType >
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< LevelDataType, FluxDataType > *> m_ops
 
RefCountedPtr< AMRMultiGrid< LevelDataType > > m_solver
 The multigrid solver used to solve the Helmholtz equation. More...
 

Private Member Functions

void updateSolnWithTimeIndependentOp (LevelDataType &a_phiNew, LevelDataType &a_phiOld, LevelDataType &a_src, LevelData< FluxDataType > &a_flux, FluxRegisterType *a_fineFluxRegPtr, FluxRegisterType *a_crseFluxRegPtr, const LevelDataType *a_crsePhiOldPtr, const LevelDataType *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, int a_fluxStartComponent=0)
 
void updateSolnWithTimeDependentOp (LevelDataType &a_phiNew, LevelDataType &a_phiOld, LevelDataType &a_src, LevelData< FluxDataType > &a_flux, FluxRegisterType *a_fineFluxRegPtr, FluxRegisterType *a_crseFluxRegPtr, const LevelDataType *a_crsePhiOldPtr, const LevelDataType *a_crsePhiNewPtr, Real a_oldTime, Real a_crseOldTime, Real a_crseNewTime, Real a_dt, int a_level, bool a_zeroPhi=true, int a_fluxStartComponent=0)
 
BaseLevelTGAoperator= (const BaseLevelTGA &)
 
 BaseLevelTGA (const BaseLevelTGA &a_opin)
 
 BaseLevelTGA ()
 

Detailed Description

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
class BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >

This base class implements the 2nd-order implicit L0-stable time integration algorithm developed by Twizell, Gumel, and Arigu for solving elliptic equations. It relies upon linear algebraic operations defined in the underlying Helmholtz operators.

Template Parameters
LevelDataTypeThe type used to store data at a grid level. This is usually LevelData<T>, where T is some cell-centered FArrayBox type.
FluxDataTypeThe type used to store flux data at a grid level. This is usually an array box clas that stores fluxes.
FluxRegisterTypeThe type used to store flux register data for interactions between different grid levels. This is usually a flux register class.

Constructor & Destructor Documentation

◆ BaseLevelTGA() [1/3]

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::BaseLevelTGA ( const Vector< DisjointBoxLayout > &  a_grids,
const Vector< int > &  a_refRat,
const ProblemDomain a_level0Domain,
RefCountedPtr< AMRLevelOpFactory< LevelDataType > > &  a_opFact,
const RefCountedPtr< AMRMultiGrid< LevelDataType > > &  a_solver 
)
inline

Initializes the base class of a TGA time integrator. This must be called by any subclass of BaseLevelTGA.

Parameters
a_gridsThe DisjointBoxLayout on which the TGA scheme is to operate.
a_refRatioAn array containing the refinement ratios between the hierarchical AMR grid levels for the domain.
a_level0DomainThe domain at the coarsest AMR grid level.
a_opFactA factory typename LevelDataTypehat is used to generate Helmholtz operators to be used by the scheme.
a_solverAn AMR Multigrid solver for solving the linear systems at each stage of the TGA integration scheme.

◆ ~BaseLevelTGA()

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
virtual BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::~BaseLevelTGA ( )
inlinevirtual

Destructor, called after destructors of BaseLevelTGA subclasses.

◆ BaseLevelTGA() [2/3]

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::BaseLevelTGA ( const BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType > &  a_opin)
private

◆ BaseLevelTGA() [3/3]

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::BaseLevelTGA ( )
private

Member Function Documentation

◆ updateSoln()

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
void BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::updateSoln ( LevelDataType &  a_phiNew,
LevelDataType &  a_phiOld,
LevelDataType &  a_src,
LevelData< FluxDataType > &  a_flux,
FluxRegisterType *  a_fineFluxRegPtr,
FluxRegisterType *  a_crseFluxRegPtr,
const LevelDataType *  a_crsePhiOldPtr,
const LevelDataType *  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 
)
inlinevirtual

Integrates the helmholtz equation represented by this object, placing the new solution in a_phiNew.

Parameters
a_phiNewThe new solution (the value of phi at time n + 1) will be stored here.
a_phiOldThe old solution (the value of phi at time n).
a_srcThe source term on the right hand side of the Helmholtz equation.
a_fluxThis will store the flux computed at the current grid level during the solution of the Helmholtz equation.
a_fineFluxRegPtrA pointer to the flux register representing the finer grid level adjacent to this one, or NULL if there is no finer grid level.
a_crseFluxRegPtrA pointer to the flux register representing the coarser grid level adjacent to this one, or NULL if there is no coarser grid level.
a_oldTimeThe time at the beginning of the integration step at the current grid level.
a_crseOldTimeThe time at the beginning of the integration step at the coarser adjacent grid level. This parameter is ignored if there is no coarser grid level.
a_crseNewTimeThe time at the end of the integration step at the coarser adjacent grid level. This parameter is ignored if there is no coarser grid level.
a_dtThe size of the integration step at the current grid level.
a_levelThe current grid level.
a_zeroPhiIf set to true, a_phiNew will be set to zero before the integration takes place. Otherwise, a_phiNew is assumed to be an initial estimate for the solution in the iterative linear solve.
a_fluxStartComponentAn index identifying the component at which flux data begins within a_fineFluxRegPtr and a_crseFluxRegPtr.

Implements BaseLevelHeatSolver< LevelDataType, FluxDataType, FluxRegisterType >.

Referenced by BaseLevelTGA< LevelData< FArrayBox >, FluxBox, LevelFluxRegister >::computeDiffusion().

◆ computeDiffusion()

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
void BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::computeDiffusion ( LevelDataType &  a_diffusiveTerm,
LevelDataType &  a_phiOld,
LevelDataType &  a_src,
LevelData< FluxDataType > &  a_flux,
FluxRegisterType *  a_fineFluxRegPtr,
FluxRegisterType *  a_crseFluxRegPtr,
const LevelDataType *  a_crsePhiOldPtr,
const LevelDataType *  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 
)
inlinevirtual

Computes the time-centered diffusion term L(phi). This can be used to find contributions to the solution from diffusion. The diffusion term is computed by computing a finite difference approximation for d phi/dt using the updated and original values of phi and the time step. Most of the arguments given here are passed along to updateSoln and retain their significance therein.

Parameters
a_diffusiveTermThe diffusion term L(phi) will be stored here.
a_phiOldThe old solution (the value of phi at time n).
a_srcThe source term on the right hand side of the Helmholtz equation (used to fine the new value of phi).
a_fineFluxRegPtrA pointer to the flux register representing the finer grid level adjacent to this one, or NULL if there is no finer grid level.
a_crseFluxRegPtrA pointer to the flux register representing the coarser grid level adjacent to this one, or NULL if there is no coarser grid level.
a_crsePhiOldTimeA pointer to the value of phi at the beginning of the integration step at the coarser grid level adjacent to this one, or NULL if there is no coarser grid level.
a_crsePhiNewTimeA pointer to the value of phi at the end of the integration step at the coarser grid level adjacent to this one, or NULL if there is no coarser grid level.
a_oldTimeThe time at the beginning of the integration step at the current grid level.
a_crseOldTimeThe time at the beginning of the integration step at the coarser adjacent grid level. This parameter is ignored if there is no coarser grid level.
a_crseNewTimeThe time at the end of the integration step at the coarser adjacent grid level. This parameter is ignored if there is no coarser grid level.
a_dtThe size of the integration step at the current grid level.
a_levelThe current grid level.
a_zeroPhiIf set to true, the new value of phi will be set to zero before the integration takes place. Otherwise, it will be set to the value in a_diffusiveTerm.

Reimplemented from BaseLevelHeatSolver< LevelDataType, FluxDataType, FluxRegisterType >.

◆ setSourceGhostCells()

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
virtual void BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::setSourceGhostCells ( LevelDataType &  a_src,
const DisjointBoxLayout a_dbl,
int  a_level 
)
protectedpure 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.

Implemented in EBLevelTGA, and LevelTGA.

Referenced by BaseLevelTGA< LevelData< FArrayBox >, FluxBox, LevelFluxRegister >::computeDiffusion(), BaseLevelTGA< LevelData< FArrayBox >, FluxBox, LevelFluxRegister >::updateSolnWithTimeDependentOp(), and BaseLevelTGA< LevelData< FArrayBox >, FluxBox, LevelFluxRegister >::updateSolnWithTimeIndependentOp().

◆ updateSolnWithTimeIndependentOp()

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
void BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::updateSolnWithTimeIndependentOp ( LevelDataType &  a_phiNew,
LevelDataType &  a_phiOld,
LevelDataType &  a_src,
LevelData< FluxDataType > &  a_flux,
FluxRegisterType *  a_fineFluxRegPtr,
FluxRegisterType *  a_crseFluxRegPtr,
const LevelDataType *  a_crsePhiOldPtr,
const LevelDataType *  a_crsePhiNewPtr,
Real  a_oldTime,
Real  a_crseOldTime,
Real  a_crseNewTime,
Real  a_dt,
int  a_level,
bool  a_zeroPhi = true,
int  a_fluxStartComponent = 0 
)
inlineprivate

◆ updateSolnWithTimeDependentOp()

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
void BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::updateSolnWithTimeDependentOp ( LevelDataType &  a_phiNew,
LevelDataType &  a_phiOld,
LevelDataType &  a_src,
LevelData< FluxDataType > &  a_flux,
FluxRegisterType *  a_fineFluxRegPtr,
FluxRegisterType *  a_crseFluxRegPtr,
const LevelDataType *  a_crsePhiOldPtr,
const LevelDataType *  a_crsePhiNewPtr,
Real  a_oldTime,
Real  a_crseOldTime,
Real  a_crseNewTime,
Real  a_dt,
int  a_level,
bool  a_zeroPhi = true,
int  a_fluxStartComponent = 0 
)
inlineprivate

◆ operator=()

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
BaseLevelTGA& BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::operator= ( const BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType > &  )
private

Member Data Documentation

◆ m_mu1

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
Real BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::m_mu1
protected

The times within the integration step at which the operators are evaluated.

Referenced by BaseLevelTGA< LevelData< FArrayBox >, FluxBox, LevelFluxRegister >::BaseLevelTGA().

◆ m_mu2

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
Real BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::m_mu2
protected

◆ m_mu3

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
Real BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::m_mu3
protected

◆ m_mu4

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
Real BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::m_mu4
protected

◆ m_r1

template<class LevelDataType, class FluxDataType, class FluxRegisterType>
Real BaseLevelTGA< LevelDataType, FluxDataType, FluxRegisterType >::m_r1
protected

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