00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011
00012
00013
00014
00015 #ifndef _MULTIGRID_H_
00016 #define _MULTIGRID_H_
00017
00018 #include "CH_Timer.H"
00019 #include "LinearSolver.H"
00020 #include "ProblemDomain.H"
00021 #include "REAL.H"
00022 #include "Box.H"
00023 #include "parstream.H"
00024 #include "RefCountedPtr.H"
00025 #include <vector>
00026 #include "NamespaceHeader.H"
00027
00028
00029 template <typename T>
00030 class MGLevelOp;
00031
00032
00033
00034
00035 template <typename T>
00036 class MGLevelOpObserver
00037 {
00038 public:
00039
00040
00041 MGLevelOpObserver()
00042 :m_op(NULL)
00043 {
00044 }
00045
00046
00047 virtual ~MGLevelOpObserver();
00048
00049
00050
00051
00052
00053 virtual void operatorChanged(const MGLevelOp<T>& a_operator)
00054 {
00055 }
00056
00057
00058
00059
00060 void setObservee(MGLevelOp<T>* a_observee)
00061 {
00062 CH_assert(m_op == NULL);
00063 m_op = a_observee;
00064 }
00065
00066
00067
00068 void clearObservee()
00069 {
00070 m_op = NULL;
00071 }
00072
00073 private:
00074
00075
00076 MGLevelOp<T>* m_op;
00077
00078
00079 MGLevelOpObserver(const MGLevelOpObserver&);
00080 MGLevelOpObserver& operator=(const MGLevelOpObserver&);
00081 };
00082
00083
00084
00085
00086
00087
00088 template <typename T>
00089 class MGLevelOp : public LinearOp<T>, public MGLevelOpObserver<T>
00090 {
00091 public:
00092
00093
00094 MGLevelOp()
00095 :LinearOp<T>(),
00096 MGLevelOpObserver<T>(),
00097 m_observers(),
00098 m_coarseningFactors()
00099 {
00100 }
00101
00102
00103 virtual ~MGLevelOp()
00104 {
00105
00106 for (size_t i = 0; i < m_observers.size(); ++i)
00107 m_observers[i]->clearObservee();
00108 }
00109
00110
00111
00112
00113
00114
00115
00116 virtual void createCoarser(T& a_coarse, const T& a_fine, bool ghosted) = 0;
00117
00118
00119
00120
00121
00122
00123
00124
00125 virtual void relax(T& a_correction, const T& a_residual, int a_iterations) = 0 ;
00126
00127
00128
00129
00130
00131 virtual void relaxNF(T& a_phi, const T* a_phiCoarse, const T& a_rhs, int a_iterations)
00132 {
00133 this->relax(a_phi, a_rhs, a_iterations);
00134 }
00135
00136
00137
00138
00139
00140
00141 virtual void restrictResidual(T& a_resCoarse, T& a_phiFine, const T& a_rhsFine) = 0;
00142
00143
00144 virtual void restrictResidual(T& a_resCoarse, T& a_phiFine, const T* a_phiCoarse, const T& a_rhsFine, bool homogeneous)
00145 {
00146 restrictResidual(a_resCoarse, a_phiFine, a_rhsFine);
00147 }
00148
00149
00150 virtual void restrictR(T& a_phiCoarse, const T& a_phiFine)
00151 {
00152
00153 }
00154
00155
00156
00157
00158
00159
00160 virtual void prolongIncrement(T& a_phiThisLevel, const T& a_correctCoarse) = 0;
00161
00162
00163
00164
00165
00166
00167 void addObserver(MGLevelOpObserver<T>* a_observer)
00168 {
00169 if (std::find(m_observers.begin(), m_observers.end(), a_observer) == m_observers.end())
00170 {
00171 m_observers.push_back(a_observer);
00172 m_coarseningFactors.push_back(-1);
00173 a_observer->setObservee(this);
00174 }
00175 }
00176
00177
00178
00179
00180
00181 virtual void applyOpMg(T& a_lhs, T& a_phi, T* a_phiCoarse, bool a_homogeneous)
00182 {
00183
00184 }
00185
00186
00187 virtual void residualNF(T& a_lhs, T& a_phi, const T* a_phiCoarse, const T& a_rhs, bool a_homogeneous = false)
00188 {
00189 this->residual(a_lhs, a_phi, a_rhs, a_homogeneous);
00190 }
00191
00192
00193
00194
00195 void removeObserver(MGLevelOpObserver<T>* a_observer)
00196 {
00197 typename std::vector<MGLevelOpObserver<T>*>::iterator iter =
00198 std::find(m_observers.begin(), m_observers.end(), a_observer);
00199 if (iter != m_observers.end())
00200 {
00201 (*iter)->clearObservee();
00202 m_observers.erase(iter);
00203
00204 size_t index = std::distance(m_observers.begin(), iter);
00205 m_coarseningFactors.erase(m_coarseningFactors.begin() + index);
00206 }
00207 }
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219 void addCoarserObserver(MGLevelOp<T>* a_operator,
00220 int a_coarseningFactor)
00221 {
00222 CH_assert(a_operator != NULL);
00223 CH_assert(a_operator != this);
00224
00225 if (std::find(m_observers.begin(), m_observers.end(), a_operator) == m_observers.end())
00226 {
00227
00228 addObserver(a_operator);
00229
00230
00231 m_coarseningFactors.back() = a_coarseningFactor;
00232 }
00233 }
00234
00235
00236 void notifyObserversOfChange()
00237 {
00238
00239
00240
00241 for (int i = 0; i < m_observers.size(); ++i)
00242 {
00243 if (m_coarseningFactors[i] != -1)
00244 {
00245 MGLevelOp<T>* op = dynamic_cast<MGLevelOp<T>*>(m_observers[i]);
00246 CH_assert(op != NULL);
00247
00248 op->finerOperatorChanged(*this, m_coarseningFactors[i]);
00249 }
00250 else
00251 {
00252
00253 m_observers[i]->operatorChanged(*this);
00254 }
00255 }
00256 }
00257
00258
00259
00260
00261
00262
00263
00264
00265 virtual void finerOperatorChanged(const MGLevelOp<T>& a_operator,
00266 int a_coarseningFactor)
00267 {
00268 }
00269
00270
00271 int numObservers() const
00272 {
00273 return m_observers.size();
00274 }
00275
00276 private:
00277
00278
00279 std::vector<MGLevelOpObserver<T>*> m_observers;
00280
00281
00282 std::vector<int> m_coarseningFactors;
00283
00284
00285 MGLevelOp(const MGLevelOp&);
00286 MGLevelOp& operator=(const MGLevelOp&);
00287 };
00288
00289
00290
00291
00292
00293 template <class T>
00294 class MGLevelOpFactory
00295 {
00296 public:
00297
00298
00299 MGLevelOpFactory()
00300 {
00301 }
00302
00303
00304 virtual ~MGLevelOpFactory()
00305 {
00306 }
00307
00308
00309
00310
00311
00312
00313 virtual MGLevelOp<T>* MGnewOp(const ProblemDomain& a_FineindexSpace,
00314 int a_depth,
00315 bool a_homoOnly = true) = 0;
00316
00317 private:
00318
00319
00320 MGLevelOpFactory(const MGLevelOpFactory&);
00321 MGLevelOpFactory& operator=(const MGLevelOpFactory&);
00322 };
00323
00324
00325
00326
00327
00328
00329
00330 template <class T>
00331 class MultiGrid
00332 {
00333 public:
00334 MultiGrid();
00335 virtual ~MultiGrid();
00336
00337
00338
00339
00340
00341
00342
00343
00344
00345
00346 virtual void define(MGLevelOpFactory<T>& a_factory,
00347 LinearSolver<T>* a_bottomSolver,
00348 const ProblemDomain& a_domain,
00349 int a_maxDepth = -1,
00350 MGLevelOp<T>* a_finestLevelOp = NULL);
00351
00352
00353
00354
00355
00356
00357
00358
00359 virtual void solve(T& a_phi, const T& a_rhs, Real a_tolerance, int a_maxIterations, int verbosity= 0);
00360
00361
00362
00363
00364
00365
00366
00367
00368
00369 virtual void oneCycle(T& a_e, const T& a_res);
00370
00371
00372 void init(const T& a_correction, const T& a_residual);
00373
00374
00375 void cycle(int a_depth, T& a_correction, const T& a_residual);
00376
00377
00378 void clear();
00379
00380
00381 void setBottomSolver(LinearSolver<T>* a_bottomSolver);
00382
00383
00384
00385
00386
00387
00388
00389
00390 int m_depth, m_defaultDepth, m_pre, m_post, m_bottom, m_cycle, m_numMG;
00391 bool m_homogeneous;
00392 LinearSolver<T>* m_bottomSolver;
00393
00394
00395
00396
00397
00398 Vector< MGLevelOp<T>* > getAllOperators();
00399
00400 ProblemDomain m_topLevelDomain;
00401 protected:
00402 bool m_defined;
00403
00404 int m_bottomCells;
00405
00406 Vector<MGLevelOp<T>*> m_op;
00407 Vector< T* > m_residual;
00408 Vector< T* > m_correction;
00409 std::vector<bool> m_ownOp;
00410
00411 private:
00412 MultiGrid(const MultiGrid<T>& a_opin)
00413 {
00414 MayDay::Error("invalid operator");
00415 }
00416
00417 void operator=(const MultiGrid<T>& a_opin)
00418 {
00419 MayDay::Error("invalid operator");
00420 }
00421 };
00422
00423
00424
00425
00426
00427
00428 template <class T>
00429 Vector< MGLevelOp<T>* >
00430 MultiGrid<T>::getAllOperators()
00431 {
00432 Vector<MGLevelOp<T>*> retval = m_op;
00433 return retval;
00434 }
00435
00436
00437
00438 template <class T>
00439 MultiGrid<T>::MultiGrid()
00440 :m_pre(3),
00441 m_post(3),
00442 m_bottom(0),
00443 m_cycle(1),
00444 m_numMG(1),
00445 m_homogeneous(true),
00446 m_defined(false)
00447 {
00448 }
00449
00450
00451
00452 template <class T>
00453 MultiGrid<T>::~MultiGrid()
00454 {
00455 clear();
00456 }
00457
00458
00459
00460 template <class T>
00461 void MultiGrid<T>::define(MGLevelOpFactory<T>& a_factory,
00462 LinearSolver<T>* a_bottomSolver,
00463 const ProblemDomain& a_domain,
00464 int a_maxDepth,
00465 MGLevelOp<T>* a_finestOp)
00466 {
00467
00468 CH_TIME("MultiGrid::define");
00469 clear();
00470 m_depth = 0;
00471 m_bottomSolver = a_bottomSolver;
00472 m_residual.resize(m_depth);
00473 m_correction.resize(m_depth);
00474 m_op.resize(m_depth);
00475 m_ownOp.resize(0);
00476 m_topLevelDomain = a_domain;
00477
00478 m_bottomCells = a_domain.domainBox().numPts();
00479
00480 MGLevelOp<T>* nextOp;
00481 bool ownOp;
00482 if (a_finestOp == NULL)
00483 {
00484 nextOp = a_factory.MGnewOp(a_domain, 0, true);
00485 ownOp = true;
00486 }
00487 else
00488 {
00489 nextOp = a_finestOp;
00490 ownOp = false;
00491 }
00492 while (nextOp != NULL)
00493 {
00494 m_residual.push_back(new T());
00495 m_correction.push_back(new T());
00496 m_op.push_back(nextOp);
00497 m_ownOp.push_back(ownOp);
00498 m_depth++;
00499 if ((m_depth < a_maxDepth) || (a_maxDepth < 0))
00500 {
00501
00502 nextOp = a_factory.MGnewOp(a_domain , m_depth, true);
00503 ownOp = true;
00504 if (nextOp != NULL)
00505 {
00506
00507
00508 m_op.back()->addCoarserObserver(nextOp, 2);
00509 for (int i=0; i<CH_SPACEDIM; ++i) m_bottomCells /= 2;
00510 }
00511 }
00512 else
00513 {
00514 nextOp = NULL;
00515 }
00516 }
00517 m_defaultDepth = m_depth;
00518 m_bottomSolver->define(m_op[m_depth-1], true);
00519 m_defined = true;
00520 }
00521
00522
00523
00524 template <class T>
00525 void MultiGrid<T>::setBottomSolver(LinearSolver<T>* a_bottomSolver)
00526 {
00527 m_bottomSolver = a_bottomSolver;
00528 m_bottomSolver->define(m_op[m_depth-1], true);
00529 }
00530
00531
00532
00533 template <class T>
00534 void MultiGrid<T>::init(const T& a_e, const T& a_residual)
00535 {
00536 CH_TIME("MutliGrid::init");
00537 if (m_depth > 1)
00538 {
00539 m_op[0]->createCoarser(*(m_residual[1]), a_residual, false);
00540 m_op[0]->createCoarser(*(m_correction[1]), a_e, true);
00541 }
00542
00543 for (int i = 2; i < m_depth; i++)
00544 {
00545 m_op[i-1]->createCoarser(*(m_residual[i]), *(m_residual[i-1]), false);
00546 m_op[i-1]->createCoarser(*(m_correction[i]), *(m_correction[i-1]), true);
00547 }
00548 }
00549
00550
00551
00552 template <class T>
00553 void MultiGrid<T>::solve(T& a_phi, const T& a_rhs, Real a_tolerance, int a_maxIterations, int verbosity)
00554 {
00555 CH_TIME("MultiGrid::solve");
00556 this->init(a_phi, a_rhs);
00557
00558 T correction, residual;
00559 m_op[0]->create(correction, a_phi);
00560 m_op[0]->create(residual, a_rhs);
00561 m_op[0]->setToZero(a_phi);
00562 m_op[0]->residual(residual, a_phi, a_rhs, false);
00563
00564 Real errorno = m_op[0]->norm(residual, 0);
00565 if (verbosity > 2)
00566 {
00567 pout() << "multigrid::solve initial residual = " << errorno << std::endl;
00568 }
00569 Real compval = a_tolerance*errorno;
00570 Real epsilon = 1.0e-16;
00571 compval = Max(compval, epsilon);
00572 Real error = errorno;
00573 int iter = 0;
00574 while ((error > compval) && (error > a_tolerance*errorno) && (iter < a_maxIterations))
00575 {
00576 m_op[0]->setToZero(correction);
00577 m_op[0]->residual(residual, a_phi, a_rhs, false);
00578 error = m_op[0]->norm(residual, 0);
00579 if (verbosity > 3)
00580 {
00581 pout() << "multigrid::solve iter = " << iter << ", residual = " << error << std::endl;
00582 }
00583
00584 this->cycle(0, correction, residual);
00585 m_op[0]->incr(a_phi, correction, 1.0);
00586
00587 iter++;
00588 }
00589 if (verbosity > 2)
00590 {
00591 pout() << "multigrid::solve final residual = " << error << std::endl;
00592 }
00593
00594 m_op[0]->clear(correction);
00595 m_op[0]->clear(residual);
00596 }
00597
00598
00599
00600 template <class T>
00601 void MultiGrid<T>::oneCycle(T& a_e, const T& a_residual)
00602 {
00603 CH_TIME("Multigrid::oneCycle");
00604
00605
00606 if (m_homogeneous)
00607 {
00608 cycle(0, a_e, a_residual);
00609 }
00610 else
00611 {
00612 T correction, residual;
00613 m_op[0]->create(correction, a_e);
00614 m_op[0]->create(residual, a_residual);
00615
00616 m_op[0]->residual(residual, a_e, a_residual);
00617
00618 m_op[0]->setToZero(correction);
00619 this->cycle(0, correction, residual);
00620
00621 m_op[0]->incr(a_e, correction, 1.0);
00622 m_op[0]->clear(correction);
00623 m_op[0]->clear(residual);
00624 }
00625 }
00626
00627
00628
00629 template <class T>
00630 void MultiGrid<T>::cycle(int depth, T& correction, const T& residual)
00631 {
00632 CH_TIME("Multigrid::cycle");
00633
00634
00635 if (depth == m_depth-1)
00636 {
00637 CH_TIME("Multigrid::cycle:bottom-solve");
00638 if (m_bottomCells == 1)
00639 {
00640 m_op[depth ]->relax(correction, residual, 1);
00641 }
00642 else
00643 {
00644 m_op[depth ]->relax(correction, residual, m_bottom);
00645 m_bottomSolver->solve(correction, residual);
00646 }
00647
00648 }
00649 else
00650 {
00651 int cycles = m_cycle;
00652 if ( cycles < 0 )
00653 {
00654 cycles = -cycles;
00655
00656 m_op[depth ]->restrictResidual(*(m_residual[depth+1]), correction, residual);
00657 m_op[depth ]->setToZero(*(m_correction[depth+1]));
00658
00659
00660 cycle(depth+1, *(m_correction[depth+1]), *(m_residual[depth+1]));
00661
00662 m_op[depth ]->prolongIncrement(correction, *(m_correction[depth+1]));
00663
00664 m_op[depth ]->relax(correction, residual, m_pre);
00665
00666 for (int img = 0; img < cycles; img++)
00667 {
00668 m_op[depth ]->restrictResidual(*(m_residual[depth+1]), correction, residual);
00669 m_op[depth ]->setToZero(*(m_correction[depth+1]));
00670
00671 m_cycle = 1;
00672 cycle(depth+1, *(m_correction[depth+1]), *(m_residual[depth+1]));
00673 m_cycle = -cycles;
00674
00675 m_op[depth ]->prolongIncrement(correction, *(m_correction[depth+1]));
00676 }
00677
00678 m_op[depth ]->relax(correction, residual, m_post);
00679 }
00680 else
00681 {
00682 m_op[depth ]->relax( correction, residual, m_pre );
00683
00684 m_op[depth ]->restrictResidual(*(m_residual[depth+1]), correction, residual);
00685 m_op[depth ]->setToZero(*(m_correction[depth+1]));
00686
00687 for (int img = 0; img < cycles; img++)
00688 {
00689 cycle(depth+1, *(m_correction[depth+1]), *(m_residual[depth+1]));
00690 }
00691 m_op[depth ]->prolongIncrement(correction, *(m_correction[depth+1]));
00692
00693 m_op[depth ]->relax(correction, residual, m_post);
00694 }
00695 }
00696 }
00697
00698
00699
00700 template <class T>
00701 void MultiGrid<T>::clear()
00702 {
00703 std::vector<bool>::reverse_iterator it = m_ownOp.rbegin();
00704 for (int i=m_op.size()-1; i>=0; --i, ++it)
00705 {
00706 if (*it) delete m_op[i];
00707 delete m_correction[i];
00708 delete m_residual[i];
00709 }
00710 m_op.resize(0);
00711 m_residual.resize(0);
00712 m_correction.resize(0);
00713 m_defined = false;
00714 }
00715
00716
00717
00718
00719
00720 template <typename T>
00721 MGLevelOpObserver<T>::~MGLevelOpObserver()
00722 {
00723 if (m_op != NULL)
00724 {
00725
00726
00727
00728
00729 m_op->removeObserver(this);
00730 }
00731 }
00732
00733
00734 #include "NamespaceFooter.H"
00735 #endif