Chombo + EB
3.2
|
#include <MultiGrid.H>
Public Member Functions | |
MultiGrid () | |
virtual | ~MultiGrid () |
virtual void | define (MGLevelOpFactory< T > &a_factory, LinearSolver< T > *a_bottomSolver, const ProblemDomain &a_domain, int a_maxDepth=-1, MGLevelOp< T > *a_finestLevelOp=NULL) |
Function to define a MultiGrid object. More... | |
virtual void | solve (T &a_phi, const T &a_rhs, Real a_tolerance, int a_maxIterations, int verbosity=0) |
virtual void | oneCycle (T &a_e, const T &a_res) |
Execute ONE v-cycle of multigrid. More... | |
void | init (const T &a_correction, const T &a_residual) |
void | cycle (int a_depth, T &a_correction, const T &a_residual) |
void | clear () |
void | setBottomSolver (LinearSolver< T > *a_bottomSolver) |
Vector< MGLevelOp< T > *> | getAllOperators () |
Public Attributes | |
int | m_depth |
int | m_defaultDepth |
int | m_pre |
int | m_post |
int | m_bottom |
int | m_cycle |
int | m_numMG |
bool | m_homogeneous |
LinearSolver< T > * | m_bottomSolver |
ProblemDomain | m_topLevelDomain |
Protected Attributes | |
bool | m_defined |
int | m_bottomCells |
Vector< MGLevelOp< T > * > | m_op |
Vector< T *> | m_residual |
Vector< T *> | m_correction |
std::vector< bool > | m_ownOp |
Private Member Functions | |
MultiGrid (const MultiGrid< T > &a_opin) | |
void | operator= (const MultiGrid< T > &a_opin) |
Class to execute geometric multigrid. This class is not meant to be particularly user-friendly and a good option for people who want something a tad less raw is to use AMRMultiGrid instead.
References MultiGrid< T >::clear().
References MayDay::Error().
|
virtual |
Function to define a MultiGrid object.
a_factory is the factory for generating operators. a_bottomSolver is called at the bottom of v-cycle. a_domain is the problem domain at the top of the vcycle. maxDepth defines the location of the bottom of the v-cycle. The vycle will terminate (hit bottom) when the factory returns NULL for a paticular depth if maxdepth = -1. Otherwise the vycle terminates at maxdepth.
Reimplemented in FASMultiGrid< T >.
References Vector< T >::back(), CH_SPACEDIM, CH_TIME, MultiGrid< T >::clear(), ProblemDomain::domainBox(), MultiGrid< T >::m_bottomCells, MultiGrid< T >::m_bottomSolver, MultiGrid< T >::m_correction, MultiGrid< T >::m_defaultDepth, MultiGrid< T >::m_defined, MultiGrid< T >::m_depth, MultiGrid< T >::m_op, MultiGrid< T >::m_ownOp, MultiGrid< T >::m_residual, MultiGrid< T >::m_topLevelDomain, MGLevelOpFactory< T >::MGnewOp(), Box::numPts(), Vector< T >::push_back(), and Vector< T >::resize().
Referenced by FASMultiGrid< T >::define().
|
virtual |
solve L(a_phi) = a_rhs. Tolerance is how much you want the norm of the error reduced. verbosity is how chatty you want the function to be. maxIterations is the maximum number of v-cycles. This does the whole residual correction switcharoo and calls oneCycle up to maxIterations times, evaluating the residual as it goes.
Reimplemented in FASMultiGrid< T >.
References CH_TIME, Vector< T >::clear(), MultiGrid< T >::cycle(), MultiGrid< T >::init(), MultiGrid< T >::m_op, Max(), and pout().
|
virtual |
Execute ONE v-cycle of multigrid.
If you want the solution to converge, you need to iterate this. See solve() or AMRMultiGrid::solve for a more automatic solve() function. This operates residual-correction form of solution so all boundary conditions are assumed to be homogeneous. L(a_e) = a_res
Reimplemented in FASMultiGrid< T >.
References CH_TIME, Vector< T >::clear(), MultiGrid< T >::cycle(), MultiGrid< T >::m_homogeneous, and MultiGrid< T >::m_op.
void MultiGrid< T >::init | ( | const T & | a_correction, |
const T & | a_residual | ||
) |
References CH_TIME, MultiGrid< T >::m_correction, MultiGrid< T >::m_depth, MultiGrid< T >::m_op, and MultiGrid< T >::m_residual.
Referenced by FASMultiGrid< T >::solve(), and MultiGrid< T >::solve().
void MultiGrid< T >::cycle | ( | int | a_depth, |
T & | a_correction, | ||
const T & | a_residual | ||
) |
References CH_TIME, MultiGrid< T >::m_bottom, MultiGrid< T >::m_bottomCells, MultiGrid< T >::m_bottomSolver, MultiGrid< T >::m_correction, MultiGrid< T >::m_cycle, MultiGrid< T >::m_depth, MultiGrid< T >::m_op, MultiGrid< T >::m_post, MultiGrid< T >::m_pre, and MultiGrid< T >::m_residual.
Referenced by MultiGrid< T >::oneCycle(), and MultiGrid< T >::solve().
void MultiGrid< T >::clear | ( | ) |
void MultiGrid< T >::setBottomSolver | ( | LinearSolver< T > * | a_bottomSolver | ) |
References MultiGrid< T >::m_bottomSolver, MultiGrid< T >::m_depth, and MultiGrid< T >::m_op.
for changing coefficients — not for the faint of heart.
References MGLevelOpObserver< T >::m_op.
References MayDay::Error().
int MultiGrid< T >::m_depth |
Public solver parameters. m_pre and m_post are the ones that usually get set and are the number of relaxations performed before and after multigrid recursion. See AMRMultiGrid for a more user-friendly interface.
Referenced by FASMultiGrid< T >::cycle(), MultiGrid< T >::cycle(), MultiGrid< T >::define(), MultiGrid< T >::init(), FASMultiGrid< T >::oneCycle(), and MultiGrid< T >::setBottomSolver().
int MultiGrid< T >::m_defaultDepth |
Referenced by MultiGrid< T >::define().
int MultiGrid< T >::m_pre |
Referenced by FASMultiGrid< T >::cycle(), and MultiGrid< T >::cycle().
int MultiGrid< T >::m_post |
Referenced by FASMultiGrid< T >::cycle(), and MultiGrid< T >::cycle().
int MultiGrid< T >::m_bottom |
Referenced by FASMultiGrid< T >::cycle(), and MultiGrid< T >::cycle().
int MultiGrid< T >::m_cycle |
Referenced by FASMultiGrid< T >::cycle(), and MultiGrid< T >::cycle().
int MultiGrid< T >::m_numMG |
bool MultiGrid< T >::m_homogeneous |
Referenced by FASMultiGrid< T >::define(), and MultiGrid< T >::oneCycle().
LinearSolver<T>* MultiGrid< T >::m_bottomSolver |
Referenced by MultiGrid< T >::cycle(), MultiGrid< T >::define(), and MultiGrid< T >::setBottomSolver().
ProblemDomain MultiGrid< T >::m_topLevelDomain |
Referenced by MultiGrid< T >::define().
|
protected |
Referenced by MultiGrid< T >::clear(), and MultiGrid< T >::define().
|
protected |
Referenced by MultiGrid< T >::cycle(), and MultiGrid< T >::define().
Referenced by MultiGrid< T >::clear(), FASMultiGrid< T >::cycle(), MultiGrid< T >::cycle(), MultiGrid< T >::define(), MultiGrid< T >::init(), FASMultiGrid< T >::oneCycle(), MultiGrid< T >::oneCycle(), MultiGrid< T >::setBottomSolver(), FASMultiGrid< T >::solve(), MultiGrid< T >::solve(), and MGLevelOpObserver< LevelData< NodeFArrayBox > >::~MGLevelOpObserver().
|
protected |
Referenced by MultiGrid< T >::clear(), and MultiGrid< T >::define().