26 #include "NamespaceHeader.H" 106 for (
size_t i = 0; i < m_observers.size(); ++i)
116 virtual void createCoarser(T& a_coarse,
const T& a_fine,
bool ghosted) = 0;
125 virtual void relax(T& a_correction,
const T& a_residual,
int a_iterations) = 0 ;
131 virtual void relaxNF(T& a_phi,
const T* a_phiCoarse,
const T& a_rhs,
int a_iterations)
133 this->relax(a_phi, a_rhs, a_iterations);
141 virtual void restrictResidual(T& a_resCoarse, T& a_phiFine,
const T& a_rhsFine) = 0;
144 virtual void restrictResidual(T& a_resCoarse, T& a_phiFine,
const T* a_phiCoarse,
const T& a_rhsFine,
bool homogeneous)
146 restrictResidual(a_resCoarse, a_phiFine, a_rhsFine);
150 virtual void restrictR(T& a_phiCoarse,
const T& a_phiFine)
160 virtual void prolongIncrement(T& a_phiThisLevel,
const T& a_correctCoarse) = 0;
169 if (std::find(m_observers.begin(), m_observers.end(), a_observer) == m_observers.end())
171 m_observers.push_back(a_observer);
172 m_coarseningFactors.push_back(-1);
181 virtual void applyOpMg(T& a_lhs, T& a_phi, T* a_phiCoarse,
bool a_homogeneous)
187 virtual void residualNF(T& a_lhs, T& a_phi,
const T* a_phiCoarse,
const T& a_rhs,
bool a_homogeneous =
false)
189 this->residual(a_lhs, a_phi, a_rhs, a_homogeneous);
197 typename std::vector<MGLevelOpObserver<T>*>::iterator iter =
198 std::find(m_observers.begin(), m_observers.end(), a_observer);
199 if (iter != m_observers.end())
201 (*iter)->clearObservee();
202 m_observers.erase(iter);
204 size_t index = std::distance(m_observers.begin(), iter);
205 m_coarseningFactors.erase(m_coarseningFactors.begin() + index);
220 int a_coarseningFactor)
225 if (std::find(m_observers.begin(), m_observers.end(), a_operator) == m_observers.end())
228 addObserver(a_operator);
231 m_coarseningFactors.back() = a_coarseningFactor;
241 for (
int i = 0; i < m_observers.size(); ++i)
243 if (m_coarseningFactors[i] != -1)
253 m_observers[i]->operatorChanged(*
this);
266 int a_coarseningFactor)
273 return m_observers.size();
315 bool a_homoOnly =
true) = 0;
359 virtual void solve(T& a_phi,
const T& a_rhs,
Real a_tolerance,
int a_maxIterations,
int verbosity= 0);
369 virtual void oneCycle(T& a_e,
const T& a_res);
372 void init(
const T& a_correction,
const T& a_residual);
375 void cycle(
int a_depth, T& a_correction,
const T& a_residual);
390 int m_depth, m_defaultDepth,
m_pre, m_post, m_bottom, m_cycle, m_numMG;
482 if (a_finestOp == NULL)
484 nextOp = a_factory.
MGnewOp(a_domain, 0,
true);
492 while (nextOp != NULL)
499 if ((
m_depth < a_maxDepth) || (a_maxDepth < 0))
508 m_op.
back()->addCoarserObserver(nextOp, 2);
543 for (
int i = 2; i <
m_depth; i++)
556 this->
init(a_phi, a_rhs);
558 T correction, residual;
559 m_op[0]->create(correction, a_phi);
560 m_op[0]->create(residual, a_rhs);
561 m_op[0]->setToZero(a_phi);
562 m_op[0]->residual(residual, a_phi, a_rhs,
false);
564 Real errorno =
m_op[0]->norm(residual, 0);
567 pout() <<
"multigrid::solve initial residual = " << errorno << std::endl;
569 Real compval = a_tolerance*errorno;
570 Real epsilon = 1.0e-16;
571 compval =
Max(compval, epsilon);
572 Real error = errorno;
574 while ((error > compval) && (error > a_tolerance*errorno) && (iter < a_maxIterations))
576 m_op[0]->setToZero(correction);
577 m_op[0]->residual(residual, a_phi, a_rhs,
false);
578 error =
m_op[0]->norm(residual, 0);
581 pout() <<
"multigrid::solve iter = " << iter <<
", residual = " << error << std::endl;
584 this->
cycle(0, correction, residual);
585 m_op[0]->incr(a_phi, correction, 1.0);
591 pout() <<
"multigrid::solve final residual = " << error << std::endl;
603 CH_TIME(
"Multigrid::oneCycle");
608 cycle(0, a_e, a_residual);
612 T correction, residual;
613 m_op[0]->create(correction, a_e);
614 m_op[0]->create(residual, a_residual);
616 m_op[0]->residual(residual, a_e, a_residual);
618 m_op[0]->setToZero(correction);
619 this->
cycle(0, correction, residual);
621 m_op[0]->incr(a_e, correction, 1.0);
637 CH_TIME(
"Multigrid::cycle:bottom-solve");
640 m_op[depth ]->relax(correction, residual, 1);
656 m_op[depth ]->restrictResidual(*(
m_residual[depth+1]), correction, residual);
664 m_op[depth ]->relax(correction, residual,
m_pre);
666 for (
int img = 0; img < cycles; img++)
668 m_op[depth ]->restrictResidual(*(
m_residual[depth+1]), correction, residual);
678 m_op[depth ]->relax(correction, residual,
m_post);
682 m_op[depth ]->relax( correction, residual,
m_pre );
684 m_op[depth ]->restrictResidual(*(
m_residual[depth+1]), correction, residual);
687 for (
int img = 0; img < cycles; img++)
693 m_op[depth ]->relax(correction, residual,
m_post);
703 std::vector<bool>::reverse_iterator it =
m_ownOp.rbegin();
704 for (
int i=
m_op.
size()-1; i>=0; --i, ++it)
706 if (*it)
delete m_op[i];
720 template <
typename T>
729 m_op->removeObserver(
this);
734 #include "NamespaceFooter.H" std::ostream & pout()
Use this in place of std::cout for program output.
void addObserver(MGLevelOpObserver< T > *a_observer)
Definition: MultiGrid.H:167
MGLevelOp()
Constructor.
Definition: MultiGrid.H:94
Definition: LinearSolver.H:28
void clear()
Definition: MultiGrid.H:701
int m_cycle
Definition: MultiGrid.H:390
#define CH_SPACEDIM
Definition: SPACE.H:51
#define CH_assert(cond)
Definition: CHArray.H:37
A class to facilitate interaction with physical boundary conditions.
Definition: ProblemDomain.H:141
Vector< T *> m_residual
Definition: MultiGrid.H:407
int m_defaultDepth
Definition: MultiGrid.H:390
MGLevelOpObserver & operator=(const MGLevelOpObserver &)
MGLevelOpObserver()
Base level Constructor. Called by all subclasses.
Definition: MultiGrid.H:41
Vector< MGLevelOp< T > *> getAllOperators()
Definition: MultiGrid.H:430
one dimensional dynamic array
Definition: Vector.H:53
void removeObserver(MGLevelOpObserver< T > *a_observer)
Definition: MultiGrid.H:195
void setObservee(MGLevelOp< T > *a_observee)
Definition: MultiGrid.H:60
virtual ~MGLevelOpFactory()
Destructor.
Definition: MultiGrid.H:304
void init(const T &a_correction, const T &a_residual)
Definition: MultiGrid.H:534
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.
Definition: MultiGrid.H:461
Vector< MGLevelOp< T > * > m_op
Definition: MultiGrid.H:406
virtual void restrictResidual(T &a_resCoarse, T &a_phiFine, const T *a_phiCoarse, const T &a_rhsFine, bool homogeneous)
full-multigrid version of restrictResidual, useful for FAS-type schemes. Defaults to standard restric...
Definition: MultiGrid.H:144
virtual void oneCycle(T &a_e, const T &a_res)
Execute ONE v-cycle of multigrid.
Definition: MultiGrid.H:601
bool m_homogeneous
Definition: MultiGrid.H:391
int m_bottomCells
Definition: MultiGrid.H:404
Definition: MultiGrid.H:36
std::vector< MGLevelOpObserver< T > * > m_observers
List of observers for the operator.
Definition: MultiGrid.H:279
MultiGrid()
Definition: MultiGrid.H:439
void cycle(int a_depth, T &a_correction, const T &a_residual)
Definition: MultiGrid.H:630
int m_post
Definition: MultiGrid.H:390
std::vector< bool > m_ownOp
Definition: MultiGrid.H:409
void setBottomSolver(LinearSolver< T > *a_bottomSolver)
Definition: MultiGrid.H:525
void resize(unsigned int isize)
Definition: Vector.H:346
ProblemDomain m_topLevelDomain
Definition: MultiGrid.H:400
void push_back(const T &in)
Definition: Vector.H:295
int m_depth
Definition: MultiGrid.H:390
virtual void restrictR(T &a_phiCoarse, const T &a_phiFine)
simple restriction function
Definition: MultiGrid.H:150
void clear()
Definition: Vector.H:180
T & back()
Definition: Vector.H:269
#define CH_TIME(name)
Definition: CH_Timer.H:82
Definition: MultiGrid.H:294
bool m_defined
Definition: MultiGrid.H:402
virtual void residualNF(T &a_lhs, T &a_phi, const T *a_phiCoarse, const T &a_rhs, bool a_homogeneous=false)
Definition: MultiGrid.H:187
double Real
Definition: REAL.H:33
Definition: MultiGrid.H:30
int numObservers() const
Returns the number of objects observing this operator.
Definition: MultiGrid.H:271
size_t size() const
Definition: Vector.H:192
int m_bottom
Definition: MultiGrid.H:390
MultiGrid(const MultiGrid< T > &a_opin)
Definition: MultiGrid.H:412
virtual void finerOperatorChanged(const MGLevelOp< T > &a_operator, int a_coarseningFactor)
Definition: MultiGrid.H:265
void clearObservee()
Definition: MultiGrid.H:68
static void Error(const char *const a_msg=m_nullString, int m_exitCode=CH_DEFAULT_ERROR_CODE)
Print out message to cerr and exit with the specified exit code.
void operator=(const MultiGrid< T > &a_opin)
Definition: MultiGrid.H:417
Definition: MultiGrid.H:331
Vector< T *> m_correction
Definition: MultiGrid.H:408
Definition: LinearSolver.H:156
virtual MGLevelOp< T > * MGnewOp(const ProblemDomain &a_FineindexSpace, int a_depth, bool a_homoOnly=true)=0
void addCoarserObserver(MGLevelOp< T > *a_operator, int a_coarseningFactor)
Definition: MultiGrid.H:219
virtual void applyOpMg(T &a_lhs, T &a_phi, T *a_phiCoarse, bool a_homogeneous)
Apply an operator.
Definition: MultiGrid.H:181
virtual ~MultiGrid()
Definition: MultiGrid.H:453
MGLevelOpFactory()
Base class constructor.
Definition: MultiGrid.H:299
void notifyObserversOfChange()
This should be called whenever the operator's data is updated.
Definition: MultiGrid.H:236
T Max(const T &a_a, const T &a_b)
Definition: Misc.H:39
std::vector< int > m_coarseningFactors
Multigrid coarsening factors.
Definition: MultiGrid.H:282
virtual ~MGLevelOp()
Destructor.
Definition: MultiGrid.H:103
MGLevelOp< T > * m_op
Definition: MultiGrid.H:76
const Box & domainBox() const
Returns the logical computational domain.
Definition: ProblemDomain.H:887
virtual ~MGLevelOpObserver()
Destructor.
Definition: MultiGrid.H:721
virtual void operatorChanged(const MGLevelOp< T > &a_operator)
Definition: MultiGrid.H:53
int m_pre
Definition: MultiGrid.H:390
virtual void relaxNF(T &a_phi, const T *a_phiCoarse, const T &a_rhs, int a_iterations)
specialized no-fine relax function, useful for full-multigrid schemes, defaults to regular relax ...
Definition: MultiGrid.H:131
virtual void solve(T &a_phi, const T &a_rhs, Real a_tolerance, int a_maxIterations, int verbosity=0)
Definition: MultiGrid.H:553
LinearSolver< T > * m_bottomSolver
Definition: MultiGrid.H:392