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