00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011
00012
00013
00014
00015 #ifndef _AMRMULTIGRID_H_
00016 #define _AMRMULTIGRID_H_
00017
00018 #include "MultiGrid.H"
00019 #include "REAL.H"
00020 #include "Box.H"
00021 #include "NoOpSolver.H"
00022 #include "parstream.H"
00023 #include "CH_Timer.H"
00024 #include "Copier.H"
00025 #include "SPMD.H"
00026 #include "Misc.H"
00027
00028 #include "NamespaceHeader.H"
00029
00030
00031
00032
00033
00034 template <typename T>
00035 class AMRLevelOp : public MGLevelOp<T>
00036 {
00037 public:
00038
00039
00040 virtual void dumpAMR(Vector<T*>& a_data, string name)
00041 {
00042 }
00043
00044 virtual void dumpLevel(T& a_data, string name)
00045 {
00046 }
00047
00048
00049 AMRLevelOp()
00050 :MGLevelOp<T>()
00051 {
00052 }
00053
00054 virtual void dumpStuff(Vector<T*> data, string filename)
00055 {
00056 }
00057
00058
00059 virtual ~AMRLevelOp()
00060 {
00061 }
00062
00063 virtual Real AMRNorm(const T& a_coarResid,
00064 const T& a_fineResid,
00065 const int& a_refRat,
00066 const int& a_ord)
00067 {
00068 return this->norm(a_coarResid, 0);
00069 }
00070
00071
00072 virtual void outputLevel(T& a_rhs, string& a_name)
00073 {
00074 }
00075
00076
00077 virtual void outputAMR(Vector<T*>& a_rhs, string& a_name)
00078 {
00079 }
00080
00081
00082
00083
00084
00085 virtual int refToCoarser() = 0;
00086
00087
00088
00089
00090
00091 virtual void AMRResidual(T& a_residual, const T& a_phiFine, const T& a_phi,
00092 const T& a_phiCoarse, const T& a_rhs,
00093 bool a_homogeneousDomBC,
00094 AMRLevelOp<T>* a_finerOp) = 0;
00095
00096
00097
00098
00099
00100
00101 virtual void AMRResidualNF(T& a_residual, const T& a_phi, const T& a_phiCoarse,
00102 const T& a_rhs, bool a_homogeneousBC) = 0;
00103
00104
00105
00106
00107
00108
00109 virtual void AMRResidualNC(T& a_residual, const T& a_phiFine, const T& a_phi,
00110 const T& a_rhs, bool a_homogeneousBC,
00111 AMRLevelOp<T>* a_finerOp) = 0;
00112
00113
00114
00115
00116
00117
00118
00119 virtual void AMROperator(T& a_LofPhi,
00120 const T& a_phiFine, const T& a_phi,
00121 const T& a_phiCoarse,
00122 bool a_homogeneousDomBC,
00123 AMRLevelOp<T>* a_finerOp) = 0;
00124
00125
00126
00127
00128
00129
00130 virtual void AMROperatorNF(T& a_LofPhi,
00131 const T& a_phi,
00132 const T& a_phiCoarse,
00133 bool a_homogeneousBC) = 0;
00134
00135
00136
00137
00138
00139
00140 virtual void AMROperatorNC(T& a_LofPhi,
00141 const T& a_phiFine,
00142 const T& a_phi,
00143 bool a_homogeneousBC,
00144 AMRLevelOp<T>* a_finerOp) = 0;
00145
00146
00147
00148
00149 virtual void AMRRestrict(T& a_resCoarse, const T& a_residual, const T& a_correction,
00150 const T& a_coarseCorrection) = 0;
00151
00152
00153
00154
00155
00156 virtual void AMRProlong(T& a_correction, const T& a_coarseCorrection) = 0;
00157
00158
00159
00160
00161
00162 virtual void AMRUpdateResidual(T& a_residual, const T& a_correction,
00163 const T& a_coarseCorrection) = 0;
00164
00165
00166
00167
00168
00169 virtual void createCoarsened(T& a_lhs,
00170 const T& a_rhs,
00171 const int& a_refRat) = 0;
00172
00173
00174
00175
00176
00177
00178
00179
00180
00181
00182
00183
00184 virtual void buildCopier(Copier& a_copier, const T& a_lhs, const T& a_rhs)
00185 {
00186 }
00187
00188 virtual void assignCopier(T& a_lhs, const T& a_rhs, const Copier& a_copier)
00189 {
00190 this->assign(a_lhs, a_rhs);
00191 }
00192
00193 virtual void zeroCovered(T& a_lhs, T& a_rhs, const Copier& a_copier)
00194 {
00195 this->setToZero(a_rhs);
00196 this->assignCopier(a_lhs, a_rhs, a_copier);
00197 }
00198 virtual Real localMaxNorm(const T& a_phi)
00199 {
00200 return this->norm(a_phi, 0);
00201 }
00202
00203
00204 virtual void AMRProlongS(T& a_correction, const T& a_coarseCorrection,
00205 T& a_temp, const Copier& a_copier)
00206 {
00207 AMRProlong(a_correction, a_coarseCorrection);
00208 }
00209 virtual void AMRRestrictS(T& a_resCoarse, const T& a_residual, const T& a_correction,
00210 const T& a_coarseCorrection, T& scratch)
00211 {
00212 AMRRestrict(a_resCoarse, a_residual, a_correction, a_coarseCorrection);
00213 }
00214
00215 virtual unsigned int orderOfAccuracy(void) const
00216 {
00217 return 2;
00218 }
00219
00220
00221 virtual void enforceCFConsistency(T& a_coarseCorrection, const T& a_correction)
00222 {
00223 }
00224 };
00225
00226
00227
00228
00229
00230 template <class T>
00231 class AMRLevelOpFactory : public MGLevelOpFactory<T>
00232 {
00233 public:
00234 virtual ~AMRLevelOpFactory()
00235 {
00236 }
00237
00238
00239
00240
00241
00242
00243 virtual AMRLevelOp<T>* AMRnewOp(const ProblemDomain& a_indexSpace)=0;
00244
00245
00246
00247
00248
00249 virtual int refToFiner(const ProblemDomain& a_indexSpace) const =0;
00250
00251 };
00252
00253
00254
00255
00256 template <class T>
00257 class AMRMultiGridInspector
00258 {
00259 public:
00260
00261
00262 AMRMultiGridInspector()
00263 {
00264 }
00265
00266
00267 virtual ~AMRMultiGridInspector()
00268 {
00269 }
00270
00271
00272
00273
00274
00275
00276
00277
00278 virtual void recordResiduals(const Vector<T*>& a_residuals,
00279 int a_minLevel,
00280 int a_maxLevel,
00281 int a_iter) = 0;
00282
00283
00284
00285
00286
00287
00288
00289
00290 virtual void recordCorrections(const Vector<T*>& a_corrections,
00291 int a_minLevel,
00292 int a_maxLevel,
00293 int a_iter) = 0;
00294
00295 private:
00296
00297 AMRMultiGridInspector(const AMRMultiGridInspector&);
00298 AMRMultiGridInspector& operator=(const AMRMultiGridInspector&);
00299 };
00300
00301
00302
00303
00304
00305 template <class T>
00306 class AMRMultiGrid
00307 {
00308 public:
00309
00310 AMRMultiGrid();
00311
00312 virtual ~AMRMultiGrid();
00313
00314 void outputAMR(Vector<T*>& a_data, string& a_name, int a_lmax, int a_lbase);
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325 virtual void define(const ProblemDomain& a_coarseDomain,
00326 AMRLevelOpFactory<T>& a_factory,
00327 LinearSolver<T>* a_bottomSolver,
00328 int a_numLevels);
00329
00330
00331
00332 void addInspector(RefCountedPtr<AMRMultiGridInspector<T> >& a_inspector)
00333 {
00334 CH_assert(!a_inspector.isNull());
00335 m_inspectors.push_back(a_inspector);
00336 }
00337
00338
00339
00340
00341
00342
00343 virtual void solve(Vector<T*>& a_phi, const Vector<T*>& a_rhs,
00344 int l_max, int l_base, bool a_zeroPhi=true,
00345 bool forceHomogeneous = false);
00346
00347
00348
00349
00350
00351 virtual void solveNoInit(Vector<T*>& a_phi, const Vector<T*>& a_rhs,
00352 int l_max, int l_base, bool a_zeroPhi=true,
00353 bool forceHomogeneous = false);
00354
00355
00356 virtual void solveNoInitResid(Vector<T*>& a_phi, Vector<T*>& a_finalResid,
00357 const Vector<T*>& a_rhs,
00358 int l_max, int l_base, bool a_zeroPhi=true,
00359 bool forceHomogeneous = false);
00360
00361 void relaxOnlyHomogeneous(Vector<T*>& a_phi, const Vector<T*>& a_rhs,
00362 int l_max, int l_base);
00363
00364 virtual void AMRVCycle(Vector<T*>& a_correction,
00365 Vector<T*>& a_residual,
00366 int l, int l_max, int l_base);
00367
00368 void setMGCycle(int a_numMG);
00369
00370 void init(const Vector<T*>& a_phi, const Vector<T*>& a_rhs,
00371 int l_max, int l_base);
00372
00373 void revert(const Vector<T*>& a_phi, const Vector<T*>& a_rhs,
00374 int l_max, int l_base);
00375
00376
00377 Real m_eps, m_hang, m_normThresh;
00378 bool m_solverParamsSet;
00379 int m_imin, m_iterMax, m_verbosity, m_exitStatus;
00380 int m_pre, m_post, m_bottom, m_numMG;
00381
00382
00383
00384 int m_maxDepth;
00385 AMRLevelOp<T>& levelOp(int level);
00386
00387
00388
00389
00390 Real m_convergenceMetric;
00391
00392
00393 Real m_bottomSolverEpsCushion;
00394
00395
00396
00397
00398
00399 Real computeAMRResidual(Vector<T*>& a_resid,
00400 Vector<T*>& a_phi,
00401 const Vector<T*>& a_rhs,
00402 int l_max,
00403 int l_base,
00404 bool a_homogeneousBC=false,
00405 bool a_computeNorm=true);
00406
00407
00408
00409 Real computeAMRResidual(Vector<T*>& a_phi,
00410 const Vector<T*>& a_rhs,
00411 int l_max,
00412 int l_min);
00413
00414
00415
00416
00417
00418 void computeAMROperator(Vector<T*>& a_lph,
00419 Vector<T*>& a_phi,
00420 int l_max,
00421 int l_base,
00422 bool a_homogeneousBC=false);
00423
00424
00425
00426
00427
00428 Vector< MGLevelOp<T>* > getAllOperators();
00429 Vector< MGLevelOp<T>* > getOperatorsOp();
00430 Vector< Vector< MGLevelOp<T>* > > getOperatorsMG();
00431 Vector< AMRLevelOp<T> * >& getAMROperators()
00432 {
00433 return m_op;
00434 }
00435
00436
00437
00438
00439
00440
00441
00442
00443
00444
00445
00446
00447
00448 void setSolverParameters(const int& a_pre,
00449 const int& a_post,
00450 const int& a_bottom,
00451 const int& a_numMG,
00452 const int& a_iterMax,
00453 const Real& a_eps,
00454 const Real& a_hang,
00455 const Real& a_normThresh);
00456
00457
00458
00459
00460
00461
00462
00463
00464 void setBottomSolver(int l_max, int l_base);
00465
00466 void setBottomSolverEpsCushion(Real a_bottomSolverEpsCushion);
00467
00468 protected:
00469
00470 void relax(T& phi, T& R, int depth, int nRelax = 2);
00471
00472 void computeAMRResidualLevel(Vector<T*>& a_resid,
00473 Vector<T*>& a_phi,
00474 const Vector<T*>& a_rhs,
00475 int l_max, int l_base, int ilev,
00476 bool a_homogeneousBC);
00477
00478 Vector<AMRLevelOp<T>*> m_op;
00479 Vector<MultiGrid<T> *> m_mg;
00480 Vector<T*> m_correction;
00481 Vector<T*> m_residual;
00482 Vector<T*> m_resC;
00483 Vector<Copier> m_resCopier;
00484 Vector<Copier> m_reverseCopier;
00485
00486 NoOpSolver<T> m_nosolve;
00487
00488 LinearSolver<T>* m_bottomSolver;
00489
00490 Vector<char> m_hasInitBeenCalled;
00491
00492 void clear();
00493
00494 private:
00495
00496
00497 Vector<RefCountedPtr<AMRMultiGridInspector<T> > > m_inspectors;
00498
00499
00500 AMRMultiGrid(const AMRMultiGrid<T>&);
00501 AMRMultiGrid& operator=(const AMRMultiGrid<T>&);
00502 };
00503
00504
00505
00506
00507
00508
00509
00510 template <class T>
00511 void
00512 AMRMultiGrid<T>::outputAMR(Vector<T*>& a_data, string& a_name, int a_lmax, int a_lbase)
00513 {
00514 Vector<T*> outputData;
00515 for (int ilev = a_lbase; ilev <= a_lmax; ilev++)
00516 {
00517 outputData.push_back(a_data[ilev]);
00518 }
00519 m_op[a_lbase]->outputAMR(outputData, a_name);
00520 }
00521
00522 template <class T>
00523 void
00524 AMRMultiGrid<T>::setSolverParameters(const int& a_pre,
00525 const int& a_post,
00526 const int& a_bottom,
00527 const int& a_numMG,
00528 const int& a_iterMax,
00529 const Real& a_eps,
00530 const Real& a_hang,
00531 const Real& a_normThresh)
00532 {
00533 m_solverParamsSet = true;
00534 m_pre = a_pre;
00535 m_post = a_post;
00536 m_bottom = a_bottom;
00537 m_eps = a_eps;
00538 m_hang = a_hang;
00539 m_normThresh = a_normThresh;
00540 m_iterMax = a_iterMax;
00541 for (int img = 0; img < m_mg.size(); img++)
00542 {
00543 m_mg[img]->m_pre = a_pre;
00544 m_mg[img]->m_post = a_post;
00545 m_mg[img]->m_bottom = a_bottom;
00546 }
00547 setMGCycle(a_numMG);
00548 m_bottomSolverEpsCushion = 1.0;
00549 }
00550 template <class T>
00551 Vector< MGLevelOp<T> * >
00552 AMRMultiGrid<T>::getAllOperators()
00553 {
00554 Vector< MGLevelOp<T>* > retval;
00555 for (int iop = 0; iop < m_op.size(); iop++)
00556 {
00557 MGLevelOp<T>* operPtr = (MGLevelOp<T>*) m_op[iop];
00558 retval.push_back(operPtr);
00559 }
00560
00561 for (int img = 0; img < m_mg.size(); img++)
00562 {
00563 Vector< MGLevelOp<T>* > mgOps = m_mg[img]->getAllOperators();
00564 retval.append(mgOps);
00565 }
00566 return retval;
00567 }
00568
00569 template <class T>
00570 Vector< MGLevelOp<T> * >
00571 AMRMultiGrid<T>::getOperatorsOp()
00572 {
00573 Vector< MGLevelOp<T>* > retval;
00574 for (int iop = 0; iop < m_op.size(); iop++)
00575 {
00576 MGLevelOp<T>* operPtr = (MGLevelOp<T>*) m_op[iop];
00577 retval.push_back(operPtr);
00578 }
00579 return retval;
00580 }
00581
00582 template <class T>
00583 Vector<Vector< MGLevelOp<T> * > >
00584 AMRMultiGrid<T>::getOperatorsMG()
00585 {
00586 Vector< Vector< MGLevelOp<T>* > > retval(m_mg.size());
00587
00588 for (int img = 0; img < m_mg.size(); img++)
00589 {
00590 retval[img] = m_mg[img]->getAllOperators();
00591 }
00592 return retval;
00593 }
00594
00595 template <class T>
00596 AMRMultiGrid<T>::AMRMultiGrid()
00597 :m_eps(1E-6),
00598 m_hang(1E-15),
00599 m_normThresh(1E-30),
00600 m_imin(5),
00601 m_iterMax(20),
00602 m_verbosity(3),
00603 m_pre(2),
00604 m_post(2),
00605 m_bottom(2),
00606 m_numMG(1),
00607 m_maxDepth(-1),
00608 m_convergenceMetric(0.),
00609 m_bottomSolverEpsCushion(1.0),
00610 m_bottomSolver(NULL),
00611 m_inspectors()
00612 {
00613 m_solverParamsSet = false;
00614
00615
00616
00617 }
00618 template <class T>
00619 void AMRMultiGrid<T>::setMGCycle(int a_numMG)
00620 {
00621 for (int ilev = 0; ilev < m_op.size(); ilev++)
00622 {
00623 m_mg[ilev]->m_numMG = a_numMG;
00624 m_mg[ilev]->m_cycle = a_numMG;
00625 }
00626 m_numMG = a_numMG;
00627 }
00628 template <class T>
00629 void AMRMultiGrid<T>::relax(T& a_correction, T& a_residual, int depth, int a_numSmooth)
00630 {
00631 CH_TIME("AMRMultiGrid::relax");
00632
00633
00634 if (m_op[depth]->refToCoarser() > 2)
00635 {
00636
00637
00638
00639 int intermediateDepth = 0;
00640 int r = m_op[depth]->refToCoarser();
00641 while (r>2)
00642 {
00643 r/=2;
00644 intermediateDepth++;
00645 }
00646 int tmp = m_mg[depth]->m_depth;
00647 m_mg[depth]->m_depth = intermediateDepth;
00648
00649
00650
00651
00652 m_mg[depth]->cycle(0, a_correction, a_residual);
00653 m_mg[depth]->m_depth = tmp;
00654 }
00655 else
00656 {
00657
00658
00659 m_op[depth]->relax(a_correction, a_residual, a_numSmooth);
00660 }
00661
00662 }
00663
00664 template <class T>
00665 AMRMultiGrid<T>::~AMRMultiGrid()
00666 {
00667 CH_TIME("~AMRMultiGrid");
00668 clear();
00669 }
00670
00671 template <class T>
00672 AMRLevelOp<T>& AMRMultiGrid<T>::levelOp(int level)
00673 {
00674 return *(m_op[level]);
00675 }
00676
00677
00678 template <class T>
00679 Real AMRMultiGrid<T>::computeAMRResidual(Vector<T*>& a_resid,
00680 Vector<T*>& a_phi,
00681 const Vector<T*>& a_rhs,
00682 int l_max,
00683 int l_base,
00684 bool a_homogeneousBC,
00685 bool a_computeNorm)
00686 {
00687 CH_TIME("AMRMultiGrid::computeAMRResidual");
00688
00689 Real rnorm = 0;
00690 Real localNorm = 0;
00691 for (int ilev = l_base; ilev <= l_max; ilev++)
00692 {
00693
00694 computeAMRResidualLevel(a_resid,
00695 a_phi,
00696 a_rhs,
00697 l_max, l_base, ilev, a_homogeneousBC);
00698 if (a_computeNorm)
00699 {
00700 if (ilev == l_max)
00701 {
00702 localNorm =m_op[ilev]->localMaxNorm(*a_resid[ilev]);
00703 }
00704 else
00705 {
00706 m_op[ilev]->zeroCovered(*a_resid[ilev], *m_resC[ilev+1], m_resCopier[ilev+1]);
00707 localNorm = m_op[ilev]->localMaxNorm(*a_resid[ilev]);
00708 }
00709 const int normType = 0;
00710 if (normType==2)
00711 localNorm = m_op[ilev]->norm(*a_resid[ilev],normType);
00712 rnorm = Max(localNorm, rnorm);
00713 }
00714 }
00715 #ifdef CH_MPI
00716 if (a_computeNorm)
00717 {
00718 CH_TIME("MPI_Allreduce");
00719 Real recv;
00720 int result = MPI_Allreduce(&rnorm, &recv, 1, MPI_CH_REAL,
00721 MPI_MAX, Chombo_MPI::comm);
00722 if (result != MPI_SUCCESS)
00723 {
00724
00725 MayDay::Error("sorry, but I had a communcation error on norm");
00726 }
00727 rnorm = recv;
00728 }
00729 #endif
00730
00731 return rnorm;
00732 }
00733
00734 template <class T>
00735 Real AMRMultiGrid<T>::computeAMRResidual(Vector<T*>& a_phi,
00736 const Vector<T*>& a_rhs,
00737 int l_max,
00738 int l_min)
00739 {
00740 return computeAMRResidual(m_residual,
00741 a_phi,
00742 a_rhs,
00743 l_max,
00744 l_min,
00745 false,
00746 true);
00747
00748
00749 }
00750
00751 template <class T>
00752 void AMRMultiGrid<T>::computeAMROperator(Vector<T*>& a_lph,
00753 Vector<T*>& a_phi,
00754 int l_max,
00755 int l_base,
00756 bool a_homogeneousBC)
00757 {
00758 for (int ilev = l_base; ilev <= l_max; ilev++)
00759 {
00760 m_op[ilev]->create(*m_residual[ilev], *a_lph[ilev]);
00761 m_op[ilev]->setToZero(*m_residual[ilev]);
00762 }
00763 computeAMRResidual(a_lph, a_phi, m_residual, l_max, l_base, a_homogeneousBC, false);
00764
00765
00766 for (int ilev = l_base; ilev <= l_max; ilev++)
00767 {
00768 m_op[ilev]->scale(*a_lph[ilev], -1.);
00769 }
00770 }
00771
00772
00773
00774 template <class T>
00775 void AMRMultiGrid<T>::computeAMRResidualLevel(Vector<T*>& a_resid,
00776 Vector<T*>& a_phi,
00777 const Vector<T*>& a_rhs,
00778 int l_max, int l_base, int ilev,
00779 bool a_homogeneousBC)
00780 {
00781 CH_TIME("AMRMultiGrid<T>::computeAMRResidualLevel");
00782 CH_assert(m_hasInitBeenCalled[ilev]=='t');
00783
00784
00785 if (l_max != l_base)
00786 {
00787 if (ilev == l_max)
00788 {
00789 m_op[l_max]->AMRResidualNF(*(a_resid[l_max]), *(a_phi[l_max]),
00790 *(a_phi[l_max-1]), *(a_rhs[l_max]),
00791 a_homogeneousBC);
00792 }
00793 else if (ilev == l_base)
00794 {
00795 if (l_base == 0)
00796 {
00797 m_op[l_base]->AMRResidualNC(*(a_resid[l_base]), *(a_phi[l_base+1]),
00798 *(a_phi[l_base]), *(a_rhs[l_base]),
00799 a_homogeneousBC, m_op[l_base+1]);
00800 }
00801 else
00802 {
00803 m_op[l_base]->AMRResidual(*a_resid[l_base], *a_phi[l_base+1], *a_phi[l_base],
00804 *a_phi[l_base-1], *a_rhs[l_base],
00805 a_homogeneousBC, m_op[l_base+1]);
00806 }
00807 }
00808 else
00809 {
00810 m_op[ilev]->AMRResidual(*a_resid[ilev], *a_phi[ilev+1], *a_phi[ilev],
00811 *a_phi[ilev-1], *a_rhs[ilev],
00812 a_homogeneousBC, m_op[ilev+1]);
00813 }
00814 }
00815 else
00816 {
00817 CH_assert(ilev == l_base);
00818 if (l_base == 0)
00819 {
00820 m_op[l_max]->residual(*a_resid[l_max], *a_phi[l_max], *a_rhs[l_max],a_homogeneousBC);
00821 }
00822 else
00823 {
00824 m_op[l_max]->AMRResidualNF(*(a_resid[l_max]), *(a_phi[l_max]),
00825 *(a_phi[l_max-1]), *(a_rhs[l_max]),
00826 a_homogeneousBC);
00827 }
00828 }
00829
00830 }
00831
00832 template<class T>
00833 void AMRMultiGrid<T>::solve(Vector<T*>& a_phi, const Vector<T*>& a_rhs,
00834 int l_max, int l_base, bool a_zeroPhi,
00835 bool a_forceHomogeneous)
00836 {
00837 CH_TIME("AMRMultiGrid::solve");
00838 init(a_phi, a_rhs, l_max, l_base);
00839 solveNoInit(a_phi, a_rhs, l_max, l_base, a_zeroPhi, a_forceHomogeneous);
00840
00841 revert(a_phi, a_rhs, l_max, l_base);
00842 }
00843
00844 template<class T>
00845 void AMRMultiGrid<T>::solveNoInit(Vector<T*>& a_phi,
00846 const Vector<T*>& a_rhs,
00847 int l_max, int l_base, bool a_zeroPhi,
00848 bool a_forceHomogeneous)
00849 {
00850 Vector<T*> uberResidual(a_rhs.size());
00851 int lowlim = l_base;
00852 if (l_base > 0)
00853 lowlim--;
00854
00855 for (int ilev = lowlim; ilev <= l_max; ilev++)
00856 {
00857 uberResidual[ilev] = new T();
00858 if (ilev >= l_base)
00859 {
00860 m_op[ilev]->create(*uberResidual[ilev], *a_rhs[ilev]);
00861 }
00862 }
00863 solveNoInitResid(a_phi, uberResidual, a_rhs, l_max, l_base, a_zeroPhi, a_forceHomogeneous);
00864 for (int i = lowlim; i <= l_max; i++)
00865 {
00866 m_op[i]->clear(*uberResidual[i]);
00867 delete uberResidual[i];
00868 }
00869 }
00870
00871 template<class T>
00872 void AMRMultiGrid<T>::solveNoInitResid(Vector<T*>& a_phi, Vector<T*>& uberResidual,
00873 const Vector<T*>& a_rhs,
00874 int l_max, int l_base, bool a_zeroPhi,
00875 bool a_forceHomogeneous)
00876 {
00877 CH_TIMERS("AMRMultiGrid::solveNoInit");
00878 CH_TIMER("AMRMultiGrid::AMRVcycle", vtimer);
00879 setBottomSolver(l_max, l_base);
00880
00881 CH_assert(l_base <= l_max);
00882 CH_assert(a_rhs.size() == a_phi.size());
00883
00884
00885
00886 Vector<T*> uberCorrection(a_rhs.size());
00887
00888 int lowlim = l_base;
00889 if (l_base > 0)
00890 lowlim--;
00891
00892 bool outputIntermediates = false;
00893
00894 for (int ilev = lowlim; ilev <= l_max; ilev++)
00895 {
00896 uberCorrection[ilev] = new T();
00897 m_op[ilev]->create(*uberCorrection[ilev], *a_phi[ilev]);
00898 if (ilev >= l_base)
00899 {
00900 m_op[ilev]->create(*uberResidual[ilev], *a_rhs[ilev]);
00901 m_op[ilev]->setToZero(*(uberResidual[ilev]));
00902 }
00903 m_op[ilev]->setToZero(*(uberCorrection[ilev]));
00904 }
00905
00906 if (a_zeroPhi)
00907 for (int ilev = l_base; ilev <=l_max; ++ilev)
00908 {
00909 m_op[ilev]->setToZero(*(a_phi[ilev]));
00910 }
00911
00912
00913 Real initial_rnorm = 0;
00914 {
00915 CH_TIME("Initial AMR Residual");
00916 initial_rnorm=computeAMRResidual(uberResidual, a_phi, a_rhs, l_max, l_base, a_forceHomogeneous, true);
00917 }
00918
00919 if (m_convergenceMetric != 0.)
00920 {
00921 initial_rnorm = m_convergenceMetric;
00922 }
00923
00924 Real rnorm = initial_rnorm;
00925 Real norm_last = 2*initial_rnorm;
00926
00927
00928 m_bottomSolver->setConvergenceMetrics(initial_rnorm, m_bottomSolverEpsCushion*m_eps);
00929
00930 int iter=0;
00931 if (m_verbosity >= 2)
00932 {
00933 pout() << " AMRMultiGrid:: iteration = " << iter << ", residual norm = " << rnorm << std::endl;
00934 }
00935
00936 bool goNorm = rnorm > m_normThresh;
00937 bool goRedu = rnorm > m_eps*initial_rnorm;
00938 bool goIter = iter < m_iterMax;
00939 bool goHang = iter < m_imin || rnorm <(1-m_hang)*norm_last;
00940 while (goIter && goRedu && goHang && goNorm)
00941 {
00942
00943
00944 for (int i = 0; i < m_inspectors.size(); ++i)
00945 m_inspectors[i]->recordResiduals(uberResidual, l_base, l_max, iter);
00946
00947 if (outputIntermediates)
00948 {
00949 char strresname[100];
00950 sprintf(strresname, "amrmg.res.iter.%03d", iter);
00951 string nameres(strresname);
00952 outputAMR(uberResidual, nameres, l_max, l_base);
00953 }
00954
00955 norm_last = rnorm;
00956
00957
00958 CH_START(vtimer);
00959 AMRVCycle(uberCorrection, uberResidual, l_max, l_max, l_base);
00960 CH_STOP(vtimer);
00961
00962
00963 char charname[100];
00964 sprintf(charname, "resid_iter%03d.%dd.hdf5", iter, SpaceDim);
00965 string sname(charname);
00966
00967
00968 for (int i = 0; i < m_inspectors.size(); ++i)
00969 m_inspectors[i]->recordCorrections(uberCorrection, l_base, l_max, iter);
00970
00971
00972 for (int ilev = l_base; ilev <= l_max; ilev++)
00973 {
00974 if (outputIntermediates)
00975 {
00976 char strcorname[100];
00977 sprintf(strcorname, "amrmg.phi.iter.%03d", iter);
00978 string namecor(strcorname);
00979 outputAMR(a_phi, namecor, l_max, l_base);
00980 }
00981
00982 m_op[ilev]->incr(*(a_phi[ilev]), *(uberCorrection[ilev]), 1.0);
00983
00984 if (outputIntermediates)
00985 {
00986 char strcorname[100];
00987 sprintf(strcorname, "amrmg.cor.iter.%03d", iter);
00988 string namecor(strcorname);
00989 outputAMR(uberCorrection, namecor, l_max, l_base);
00990 }
00991
00992 m_op[ilev]->setToZero(*(uberCorrection[ilev]));
00993
00994 }
00995
00996
00997
00998
00999
01000 if (m_op[0]->orderOfAccuracy()>2)
01001 {
01002 for (int ilev=l_max; ilev>l_base; ilev--)
01003 {
01004 m_op[ilev]->enforceCFConsistency(*a_phi[ilev-1], *a_phi[ilev]);
01005 }
01006 }
01007
01008
01009
01010 rnorm = computeAMRResidual(uberResidual, a_phi, a_rhs, l_max, l_base, a_forceHomogeneous, true);
01011 iter++;
01012 if (m_verbosity >= 3)
01013 {
01014 pout() << " AMRMultiGrid:: iteration = " << iter << ", residual norm = " << rnorm;
01015 if (rnorm > 0.0)
01016 {
01017 pout() << ", rate = " << norm_last/rnorm;
01018 }
01019 pout() << std::endl;
01020 }
01021
01022 goNorm = rnorm > m_normThresh;
01023 goRedu = rnorm > m_eps*initial_rnorm;
01024 goIter = iter < m_iterMax;
01025 goHang = iter < m_imin || rnorm <(1-m_hang)*norm_last;
01026 }
01027 if ((rnorm > 10.*initial_rnorm) && (rnorm > 10.*m_eps))
01028 {
01029 pout() << "solver seems to have blown up" << endl;
01030 MayDay::Error("kaboom");
01031 }
01032 m_exitStatus = int(!goRedu) + int(!goIter)*2 + int(!goHang)*4 + int(!goNorm)*8;
01033 if (m_verbosity >= 2)
01034 {
01035 pout() << " AMRMultiGrid:: iteration = " << iter << ", residual norm = " << rnorm << std::endl;
01036 }
01037 if (m_verbosity > 1)
01038 {
01039 if (!goIter && goRedu && goNorm)
01040 {
01041 pout() << " AMRMultiGrid:: WARNING: Exit because max iteration count exceeded" << std::endl;
01042 }
01043 if (!goHang && goRedu && goNorm)
01044 {
01045 pout() << " AMRMultiGrid:: WARNING: Exit because of solver hang" << std::endl;
01046 }
01047 if (m_verbosity > 4)
01048 {
01049 pout() << " AMRMultiGrid:: exitStatus = " << m_exitStatus << std::endl;
01050 }
01051 }
01052 for (int i = lowlim; i <= l_max; i++)
01053 {
01054 m_op[i]->clear(*uberCorrection[i]);
01055 delete uberCorrection[i];
01056 }
01057 }
01058
01059
01060 template<class T>
01061 void AMRMultiGrid<T>::relaxOnlyHomogeneous(Vector<T*>& a_phi, const Vector<T*>& a_rhs,
01062 int l_max, int l_base)
01063 {
01064 CH_TIME("AMRMultiGrid::relaxOnly");
01065 CH_assert(l_max == l_base);
01066 CH_assert(l_max == 0);
01067 init(a_phi, a_rhs, l_max, l_base);
01068
01069 CH_assert(a_rhs.size() == a_phi.size());
01070
01071 Vector<T*> uberResidual(a_rhs.size());
01072
01073 for (int ilev = 0; ilev < m_op.size(); ilev++)
01074 {
01075 uberResidual[ilev] = new T();
01076 m_op[ilev]->create(*uberResidual[ilev], *a_rhs[ilev]);
01077 }
01078
01079
01080 m_op[0]->residual(*uberResidual[0], *a_phi[0], *a_rhs[0], true);
01081 Real initial_rnorm = m_op[0]->norm(*uberResidual[0], 0);
01082 Real rnorm = initial_rnorm;
01083 Real norm_last = 2*initial_rnorm;
01084
01085
01086 int iter=0;
01087 if (m_verbosity >= 3)
01088 {
01089 pout() << " AMRMultiGrid::relaxOnly iteration = " << iter << ", residual norm = " << rnorm << std::endl;
01090 }
01091 bool goNorm = rnorm > m_normThresh;
01092 bool goRedu = rnorm > m_eps*initial_rnorm;
01093 bool goIter = iter < m_iterMax;
01094 bool goHang = iter < m_imin || rnorm <(1-m_hang)*norm_last;
01095 while (goIter && goRedu && goHang && goNorm)
01096 {
01097 norm_last = rnorm;
01098 m_op[0]->relax(*a_phi[0], *a_rhs[0], 1);
01099
01100 iter++;
01101
01102 m_op[0]->residual(*uberResidual[0], *a_phi[0], *a_rhs[0], true);
01103 rnorm = m_op[0]->norm(*uberResidual[0], 2);
01104 if (m_verbosity >= 4)
01105 {
01106 pout() << " AMRMultiGrid::relaxOnly iteration = " << iter << ", residual norm = " << rnorm
01107 << ", rate = " << norm_last/rnorm << std::endl;
01108 }
01109 goNorm = rnorm > m_normThresh;
01110 goRedu = rnorm > m_eps*initial_rnorm;
01111 goIter = iter < m_iterMax;
01112 goHang = iter < m_imin || rnorm <(1-m_hang)*norm_last;
01113 }
01114 m_exitStatus = 0;
01115 for (int i = 0; i < m_op.size(); i++)
01116 {
01117 m_op[i]->clear(*uberResidual[i]);
01118 delete uberResidual[i];
01119 }
01120 }
01121
01122 template <class T>
01123 void AMRMultiGrid<T>::clear()
01124 {
01125 for (int i = 0; i < m_op.size(); i++)
01126 {
01127 m_op[i]->clear(*m_correction[i]);
01128 m_op[i]->clear(*m_residual[i]);
01129 m_op[i]->clear(*m_resC[i]);
01130
01131 delete m_correction[i];
01132 delete m_residual[i];
01133 delete m_resC[i];
01134 delete m_op[i];
01135 delete m_mg[i];
01136
01137 m_correction[i] = NULL;
01138 m_residual[i] = NULL;
01139 m_resC[i] = NULL;
01140 m_op[i] = NULL;
01141 m_mg[i] = NULL;
01142 }
01143 m_hasInitBeenCalled.resize(0);
01144 }
01145
01146 template <class T>
01147 void AMRMultiGrid<T>::init(const Vector<T*>& a_phi, const Vector<T*>& a_rhs,
01148 int l_max, int l_base)
01149 {
01150
01151 CH_TIME("AMRMultiGrid::init");
01152
01153
01154 for (int i = l_base; i <= l_max; i++)
01155 {
01156 m_hasInitBeenCalled[i] = 't';
01157 AMRLevelOp<T>& op = *(m_op[i]);
01158
01159 op.create(*m_correction[i], *a_phi[i]);
01160 op.create(*m_residual[i], *a_rhs[i]);
01161 int r = op.refToCoarser();
01162 if (i!= l_base)
01163 {
01164 int intermediateDepth = 0;
01165
01166 while (r>2)
01167 {
01168 r/=2;
01169 intermediateDepth++;
01170 }
01171 m_mg[i]->m_depth = intermediateDepth;
01172 }
01173 m_mg[i]->init(*a_phi[i], *a_rhs[i]);
01174
01175 if (i != l_base)
01176 {
01177 r = op.refToCoarser();
01178 op.createCoarsened(*m_resC[i], *a_rhs[i], r);
01179 op.buildCopier(m_resCopier[i], *a_rhs[i-1], *m_resC[i]);
01180 m_reverseCopier[i] = m_resCopier[i];
01181 m_reverseCopier[i].reverse();
01182 }
01183 }
01184 }
01185 template <class T>
01186 void AMRMultiGrid<T>::revert(const Vector<T*>& a_phi, const Vector<T*>& a_rhs,
01187 int l_max, int l_base)
01188 {
01189 CH_TIME("AMRMultiGrid::revert");
01190 for (int i = l_base; i <= l_max; i++)
01191 {
01192 if (i!= l_base)
01193 {
01194 m_mg[i]->m_depth = m_mg[i]->m_defaultDepth;
01195 }
01196 }
01197 }
01198
01199 template <class T>
01200 void AMRMultiGrid<T>::setBottomSolver(int l_max, int l_base)
01201 {
01202
01203 for (int ilev = l_base; ilev <= l_max; ilev++)
01204 {
01205
01206 if (!m_solverParamsSet)
01207 {
01208 m_mg[ilev]->m_pre = m_pre;
01209 m_mg[ilev]->m_post = m_post;
01210 m_mg[ilev]->m_bottom = m_bottom;
01211 }
01212 m_mg[ilev]->m_bottomSolver = &m_nosolve;
01213 }
01214
01215 m_mg[l_base]->setBottomSolver(m_bottomSolver);
01216 }
01217
01218 template <class T>
01219 void AMRMultiGrid<T>::setBottomSolverEpsCushion(Real a_bottomSolverEpsCushion)
01220 {
01221 CH_assert((a_bottomSolverEpsCushion > 0.) &&
01222 (a_bottomSolverEpsCushion <= 1.));
01223 m_bottomSolverEpsCushion = a_bottomSolverEpsCushion;
01224 }
01225
01226 template <class T>
01227 void AMRMultiGrid<T>::define(const ProblemDomain& a_coarseDomain,
01228 AMRLevelOpFactory<T>& a_factory,
01229 LinearSolver<T>* a_bottomSolver,
01230 int a_maxAMRLevels)
01231 {
01232 CH_TIME("AMRMultiGrid::define");
01233 this->clear();
01234 m_op.resize( a_maxAMRLevels, NULL);
01235 m_mg.resize( a_maxAMRLevels, NULL);
01236 m_hasInitBeenCalled.resize(a_maxAMRLevels, 'f');
01237
01238 m_correction.resize( a_maxAMRLevels, NULL);
01239 m_residual. resize( a_maxAMRLevels, NULL);
01240 m_resC. resize( a_maxAMRLevels, NULL);
01241 m_resCopier. resize( a_maxAMRLevels);
01242 m_reverseCopier.resize( a_maxAMRLevels);
01243 m_bottomSolver = a_bottomSolver;
01244
01245 ProblemDomain current = a_coarseDomain;
01246 for (int i = 0; i < a_maxAMRLevels; i++)
01247 {
01248
01249 m_correction[i] = new T();
01250 m_residual[i] = new T();
01251 m_resC[i] = new T();
01252 m_mg[i]= new MultiGrid<T>();
01253 m_op[i] = a_factory.AMRnewOp(current);
01254 m_mg[i]->define(a_factory, &m_nosolve, current, m_maxDepth, m_op[i]);
01255
01256
01257
01258 if (i < a_maxAMRLevels-1)
01259 {
01260 current.refine(a_factory.refToFiner(current));
01261 }
01262 }
01263 }
01264
01265 template<class T>
01266 void AMRMultiGrid<T>::AMRVCycle(Vector<T*>& a_uberCorrection,
01267 Vector<T*>& a_uberResidual,
01268 int ilev, int l_max, int l_base)
01269 {
01270 if (ilev == l_max)
01271 {
01272 for (int level = l_base; level <= l_max; level++)
01273 {
01274 m_op[level]->assignLocal(*m_residual[level], *a_uberResidual[level]);
01275 m_op[level]->setToZero(*m_correction[level]);
01276 }
01277 }
01278
01279 if (l_max == l_base)
01280 {
01281 CH_assert(ilev == l_base);
01282 m_mg[l_base]->oneCycle(*(a_uberCorrection[ilev]), *(a_uberResidual[ilev]));
01283 }
01284 else if (ilev == l_base)
01285 {
01286 m_mg[l_base]->oneCycle(*(m_correction[ilev]), *(m_residual[ilev]));
01287
01288 m_op[ilev]->incr(*(a_uberCorrection[ilev]), *(m_correction[ilev]), 1.0);
01289 }
01290 else
01291 {
01292
01293
01294 this->relax(*(m_correction[ilev]), *(m_residual[ilev]), ilev, m_pre);
01295 m_op[ilev]->incr(*(a_uberCorrection[ilev]), *(m_correction[ilev]), 1.0);
01296
01297
01298 m_op[ilev-1]->setToZero(*(m_correction[ilev-1]));
01299
01300
01301
01302 computeAMRResidualLevel(m_residual,
01303 a_uberCorrection,
01304 a_uberResidual,
01305 l_max, l_base, ilev-1,
01306 true);
01307
01308
01309 m_op[ilev]->AMRRestrictS(*(m_resC[ilev]),
01310 *(m_residual[ilev]),
01311 *(m_correction[ilev]),
01312 *(m_correction[ilev-1]),
01313 *(a_uberCorrection[ilev]));
01314
01315
01316
01317 m_op[ilev-1]->assignCopier(*m_residual[ilev-1], *(m_resC[ilev]), m_resCopier[ilev]);
01318
01319
01320
01321 for (int img = 0; img < m_numMG; img++)
01322 {
01323 AMRVCycle(a_uberCorrection, a_uberResidual, ilev-1, l_max, l_base);
01324 }
01325
01326
01327
01328
01329 m_op[ilev]->AMRProlongS(*(m_correction[ilev]), *(m_correction[ilev-1]),
01330 *m_resC[ilev], m_reverseCopier[ilev]);
01331
01332 m_op[ilev]->AMRUpdateResidual(*(m_residual[ilev]), *(m_correction[ilev]), *(m_correction[ilev-1]));
01333
01334
01335
01336 T& dCorr = *(a_uberCorrection[ilev]);
01337 m_op[ilev]->setToZero(dCorr);
01338 this->relax(dCorr, *(m_residual[ilev]), ilev, m_post);
01339
01340
01341 m_op[ilev]->incr(*(m_correction[ilev]), dCorr, 1.0);
01342
01343
01344 m_op[ilev]->assignLocal(*(a_uberCorrection[ilev]), *(m_correction[ilev]));
01345
01346 }
01347 }
01348
01349 #include "NamespaceFooter.H"
01350 #endif