00001 #ifdef CH_LANG_CC
00002
00003
00004
00005
00006
00007
00008
00009 #endif
00010
00011 #ifndef _MULTILEVELLINEAROPI_H_
00012 #define _MULTILEVELLINEAROPI_H_
00013
00014 #include "AMRIO.H"
00015
00016 #include "NamespaceHeader.H"
00017
00018
00019
00020 template<class T>
00021 MultilevelLinearOp<T>::MultilevelLinearOp()
00022 {
00023
00024 m_use_multigrid_preconditioner = true;
00025 m_num_mg_iterations = 1;
00026 m_num_mg_smooth = 4;
00027
00028 m_preCondSolverDepth = -1;
00029 m_precondBottomSolverPtr = NULL;
00030 }
00031
00032
00033 template<class T>
00034 void
00035 MultilevelLinearOp<T>::define(const Vector<DisjointBoxLayout>& a_vectGrids,
00036 const Vector<int>& a_refRatios,
00037 const Vector<ProblemDomain>& a_domains,
00038 const Vector<RealVect>& a_vectDx,
00039 RefCountedPtr<AMRLevelOpFactory<LevelData<T> > >& a_opFactory,
00040 int a_lBase)
00041 {
00042 CH_TIME("MultilevelLinearOp::define");
00043
00044 int numLevels = a_vectGrids.size();
00045
00046
00047 CH_assert (a_lBase >= 0);
00048 CH_assert (a_lBase < numLevels);
00049 CH_assert ( a_domains.size() == numLevels);
00050
00051 CH_assert (a_refRatios.size() >= numLevels -1);
00052
00053 m_lBase = a_lBase;
00054
00055 m_vectGrids = a_vectGrids;
00056 m_refRatios = a_refRatios;
00057 m_domains = a_domains;
00058 m_vectDx = a_vectDx;
00059 m_vectOperators.resize(numLevels);
00060
00061
00062 int coarsestLevel = Max(0, m_lBase-1);
00063 for (int level=coarsestLevel; level<numLevels; level++)
00064 {
00065 RefCountedPtr<AMRLevelOp<LevelData<T> > > op = RefCountedPtr<AMRLevelOp<LevelData<T> > >(a_opFactory->AMRnewOp(m_domains[level]));
00066 m_vectOperators[level] = op;
00067 }
00068
00069
00070 if (m_use_multigrid_preconditioner)
00071 {
00072
00073 BiCGStabSolver<LevelData<T> >* newPtr = new BiCGStabSolver<LevelData<T> >;
00074
00075 m_precondBottomSolverPtr = newPtr;
00076 int bottomSolverVerbosity = 1;
00077
00078 newPtr->m_verbosity = bottomSolverVerbosity;
00079
00080 if (m_preCondSolverDepth >= 0)
00081 {
00082 m_preCondSolver.m_maxDepth = m_preCondSolverDepth;
00083 }
00084 m_preCondSolver.define(m_domains[0],
00085 *a_opFactory,
00086 m_precondBottomSolverPtr,
00087 numLevels);
00088
00089
00090
00091 Real solverEps = 1.0e-7;
00092 int maxIterations = 1;
00093 Real hang = 1.0e-7;
00094 Real normThresh = 1.0e-30;
00095
00096 int num_mg = 1;
00097 m_preCondSolver.setSolverParameters(m_num_mg_smooth,
00098 m_num_mg_smooth,
00099 m_num_mg_smooth,
00100 num_mg,
00101 maxIterations,
00102 solverEps,
00103 hang,
00104 normThresh);
00105
00106
00107 m_preCondSolver.m_verbosity = 1;
00108
00109
00110
00111 m_isPrecondSolverInitialized = false;
00112
00113 }
00114 }
00115
00116
00117 template<class T>
00118 MultilevelLinearOp<T>::~MultilevelLinearOp()
00119 {
00120 CH_TIME("MultilevelLinearOp::~MultilevelLinearOp");
00121
00122 if (m_precondBottomSolverPtr != NULL)
00123 {
00124 delete m_precondBottomSolverPtr;
00125 m_precondBottomSolverPtr = NULL;
00126 }
00127 }
00128
00129
00130
00131
00132
00133
00134
00135 template<class T>
00136 void
00137 MultilevelLinearOp<T>::residual(Vector<LevelData<T>* >& a_lhs,
00138 const Vector<LevelData<T>* >& a_phi,
00139 const Vector<LevelData<T>* >& a_rhs,
00140 bool a_homogeneous )
00141 {
00142 CH_TIME("MultilevelLinearOp::residual");
00143
00144
00145 LevelData<T> dummyT;
00146
00147 int numLevels = a_lhs.size();
00148 CH_assert (a_phi.size() == numLevels);
00149 CH_assert (a_rhs.size() == numLevels);
00150 CH_assert (m_vectOperators.size() >= numLevels);
00151
00152 int maxLevel = numLevels -1;
00153 for (int level = m_lBase; level < numLevels; level++)
00154 {
00155 const LevelData<T>& thisPhi = *(a_phi[level]);
00156 LevelData<T>& thisResid = *(a_lhs[level]);
00157 const LevelData<T>& thisRHS = *(a_rhs[level]);
00158
00159 if (level >0)
00160 {
00161 const LevelData<T>& phiCrse = *a_phi[level-1];
00162
00163 if (level < maxLevel)
00164 {
00165 const LevelData<T>& phiFine = *a_phi[level+1];
00166 m_vectOperators[level]->AMRResidual(thisResid,
00167 phiFine,
00168 thisPhi,
00169 phiCrse,
00170 thisRHS,
00171 a_homogeneous,
00172 (m_vectOperators[level+1].operator->()));
00173 }
00174 else
00175 {
00176
00177 m_vectOperators[level]->AMRResidualNF(thisResid,
00178 thisPhi,
00179 phiCrse,
00180 thisRHS,
00181 a_homogeneous);
00182 }
00183
00184 }
00185 else
00186 {
00187
00188 if (level < maxLevel)
00189 {
00190 const LevelData<T>& phiFine = *a_phi[level+1];
00191 m_vectOperators[level]->AMRResidualNC(thisResid,
00192 phiFine,
00193 thisPhi,
00194 thisRHS,
00195 a_homogeneous,
00196 (m_vectOperators[level+1].operator->()));
00197 }
00198 else
00199 {
00200
00201 m_vectOperators[level]->residual(thisResid,
00202 thisPhi,
00203 thisRHS,
00204 a_homogeneous);
00205 }
00206 }
00207 }
00208 }
00209
00210
00211
00212
00213
00214
00215 template<class T>
00216 void
00217 MultilevelLinearOp<T>::preCond(Vector<LevelData<T>* >& a_cor,
00218 const Vector<LevelData<T>* >& a_residual)
00219
00220 {
00221 CH_TIME("MultilevelLinearOp::preCond");
00222
00223
00224 int numLevels = a_cor.size();
00225 CH_assert (a_residual.size() == numLevels);
00226 CH_assert (m_vectOperators.size() <= numLevels);
00227
00228
00229 if (m_use_multigrid_preconditioner)
00230 {
00231 CH_TIME("MultilevelLinearOp::preCond::Multigrid");
00232
00233 Vector<LevelData<T>* > localCorr;
00234 Vector<LevelData<T>* > localResid;
00235 create(localCorr, a_cor);
00236 create(localResid, a_residual);
00237
00238
00239
00240 Vector<LevelData<T>* > tempCorr(numLevels, NULL);
00241 Vector<LevelData<T>* > tempResid(numLevels, NULL);
00242
00243 {
00244 CH_TIME("MultilevelLinearOp::preCond::Multigrid::NewStorage");
00245
00246
00247 for (int lev=0; lev<numLevels; lev++)
00248 {
00249 tempCorr[lev] = localCorr[lev];
00250 tempResid[lev] = const_cast<LevelData<T>*>(localResid[lev]);
00251 }
00252 }
00253
00254 int finestLevel = numLevels -1;
00255 {
00256 CH_TIME("MultilevelLinearOp::preCond::Multigrid::Initialize");
00257
00258 if (!m_isPrecondSolverInitialized)
00259 {
00260
00261 m_preCondSolver.init(tempCorr, tempResid, finestLevel,
00262 m_lBase);
00263
00264 m_preCondSolver.setBottomSolver(finestLevel, m_lBase);
00265 m_isPrecondSolverInitialized = true;
00266 }
00267 }
00268
00269 for (int iter = 0; iter<m_num_mg_iterations; iter++)
00270 {
00271 CH_TIME("MultilevelLinearOp::preCond::Multigrid::Iterate");
00272
00273
00274
00275
00276
00277
00278
00279
00280 bool homogeneous = true;
00281 residual(localResid, a_cor, a_residual, homogeneous);
00282 setToZero(localCorr);
00283
00284 m_preCondSolver.AMRVCycle(tempCorr, tempResid,
00285 finestLevel, finestLevel,
00286 m_lBase);
00287
00288
00289 incr(a_cor, localCorr, 1.0);
00290 }
00291
00292 {
00293 CH_TIME("MultilevelLinearOp::preCond::Multigrid::Clear");
00294
00295 clear(localCorr);
00296 clear(localResid);
00297 }
00298 }
00299 else
00300 {
00301 CH_TIME("MultilevelLinearOp::preCond::LevelBased");
00302
00303
00304 for (int level = m_lBase; level < numLevels; level++)
00305 {
00306 m_vectOperators[level]->preCond(*a_cor[level], *a_residual[level]);
00307 }
00308 }
00309 }
00310
00311
00312
00313
00314
00315
00316
00317 template<class T>
00318 void
00319 MultilevelLinearOp<T>::applyOp(Vector<LevelData<T>* >& a_lhs,
00320 const Vector<LevelData<T>* >& a_phi,
00321 bool a_homogeneous)
00322 {
00323 CH_TIME("MultilevelLinearOp::applyOp");
00324
00325 int numLevels = a_phi.size();
00326 int maxLevel = numLevels-1;
00327 LevelData<T> dummyLDF;
00328 for (int level = m_lBase; level<numLevels; level++)
00329 {
00330 if (level == 0)
00331 {
00332 if (level < maxLevel)
00333 {
00334 m_vectOperators[level]->AMROperatorNC(*a_lhs[level],
00335 *a_phi[level+1],
00336 *a_phi[level],
00337 a_homogeneous,
00338 m_vectOperators[level+1].operator->());
00339 }
00340 else
00341 {
00342
00343 m_vectOperators[level]->applyOp(*a_lhs[level],
00344 *a_phi[level],
00345 a_homogeneous);
00346 }
00347 }
00348 else
00349 {
00350
00351 if (level < maxLevel)
00352 {
00353 m_vectOperators[level]->AMROperator(*a_lhs[level],
00354 *a_phi[level+1],
00355 *a_phi[level],
00356 *a_phi[level-1],
00357 a_homogeneous,
00358 m_vectOperators[level+1].operator->());
00359 }
00360 else
00361 {
00362 m_vectOperators[level]->AMROperatorNF(*a_lhs[level],
00363 *a_phi[level],
00364 *a_phi[level-1],
00365 a_homogeneous);
00366 }
00367 }
00368 }
00369
00370 }
00371
00372
00373
00374
00375
00376
00377 template<class T>
00378 void
00379 MultilevelLinearOp<T>::create(Vector<LevelData<T>* >& a_lhs,
00380 const Vector<LevelData<T>* >& a_rhs)
00381 {
00382 CH_TIME("MultilevelLinearOp::create");
00383
00384 a_lhs.resize(a_rhs.size());
00385
00386 int coarsestLevel = Max(0, m_lBase-1);
00387 for (int level =coarsestLevel; level < a_rhs.size(); level++)
00388 {
00389 a_lhs[level] = new LevelData<T>;
00390 m_vectOperators[level]->create(*a_lhs[level],
00391 *a_rhs[level]);
00392 }
00393 }
00394
00395
00396
00397
00398
00399
00400 template<class T>
00401 void
00402 MultilevelLinearOp<T>::clear(Vector<LevelData<T>* >& a_lhs)
00403 {
00404 CH_TIME("MultilevelLinearOp::clear");
00405
00406 for (int level =0; level < a_lhs.size(); level++)
00407 {
00408 if (a_lhs[level] != NULL)
00409 {
00410 delete a_lhs[level];
00411 a_lhs[level] = NULL;
00412 }
00413 }
00414 }
00415
00416
00417
00418
00419
00420 template<class T>
00421 void
00422 MultilevelLinearOp<T>::assign(Vector<LevelData<T>* >& a_lhs,
00423 const Vector<LevelData<T>* >& a_rhs)
00424 {
00425 CH_TIME("MultilevelLinearOp::assign");
00426
00427 CH_assert (a_lhs.size() == a_rhs.size());
00428
00429 for (int level = m_lBase; level < a_lhs.size(); level++)
00430 {
00431 if (!(a_lhs[level] ==NULL) && !(a_rhs[level] == NULL))
00432 {
00433 m_vectOperators[level]->assign(*a_lhs[level], *a_rhs[level]);
00434 }
00435 }
00436 }
00437
00438
00439
00440
00441
00442
00443 template<class T>
00444 Real
00445 MultilevelLinearOp<T>::dotProduct(const Vector<LevelData<T>* >& a_1,
00446 const Vector<LevelData<T>* >& a_2)
00447 {
00448 CH_TIME("MultilevelLinearOp::dotProduct");
00449
00450
00451 Vector<LevelData<T>* > temp1, temp2;
00452 create(temp1, a_1);
00453 create (temp2, a_2);
00454
00455
00456
00457
00458 setToZero(temp1);
00459 setToZero(temp2);
00460
00461 assign(temp1, a_1);
00462 assign(temp2, a_2);
00463
00464
00465 for (int level =m_lBase; level<temp1.size()-1; level++)
00466 {
00467 LevelData<T>& temp1Level = *temp1[level];
00468 LevelData<T>& temp2Level = *temp2[level];
00469
00470 CH_assert(temp1[level]->getBoxes() == temp2[level]->getBoxes());
00471 CH_assert(temp1[level+1] != NULL);
00472
00473 int nRefFine = m_refRatios[level];
00474 const DisjointBoxLayout& finerGrids = temp1[level+1]->getBoxes();
00475 const DisjointBoxLayout& levelGrids = temp1[level]->getBoxes();
00476 DataIterator levelDit = levelGrids.dataIterator();
00477 LayoutIterator finerLit = finerGrids.layoutIterator();
00478
00479 for (levelDit.begin(); levelDit.ok(); ++levelDit)
00480 {
00481 const Box& thisBox = levelGrids[levelDit];
00482 for (finerLit.begin(); finerLit.ok(); ++finerLit)
00483 {
00484 Box testBox = finerGrids[finerLit];
00485 testBox.coarsen(nRefFine);
00486 testBox &= thisBox;
00487 if (!testBox.isEmpty())
00488 {
00489 temp1Level[levelDit].setVal(0.0, testBox,
00490 0, temp1Level.nComp());
00491
00492 temp2Level[levelDit].setVal(0.0, testBox,
00493 0, temp2Level.nComp());
00494 }
00495 }
00496 }
00497 }
00498
00499
00500 Real prod = 0.0;
00501 for (int level =m_lBase; level<temp1.size(); level++)
00502 {
00503 Real levelProd = m_vectOperators[level]->dotProduct(*temp1[level],
00504 *temp2[level]);
00505
00506
00507 RealVect dxLev = m_vectDx[level];
00508 Real scale = D_TERM(dxLev[0], *dxLev[1], *dxLev[2]);
00509 levelProd *= scale;
00510 prod += levelProd;
00511 }
00512
00513 clear(temp1);
00514 clear (temp2);
00515
00516 return prod;
00517 }
00518
00519
00520
00521
00522
00523 template<class T>
00524 void
00525 MultilevelLinearOp<T>::incr (Vector<LevelData<T>* >& a_lhs,
00526 const Vector<LevelData<T>* >& a_x,
00527 Real a_scale)
00528 {
00529 CH_TIME("MultilevelLinearOp::incr");
00530
00531
00532 for (int level=m_lBase; level<a_lhs.size(); level++)
00533 {
00534 m_vectOperators[level]->incr(*a_lhs[level], *a_x[level], a_scale);
00535 }
00536 }
00537
00538
00539
00540
00541
00542 template<class T>
00543 void
00544 MultilevelLinearOp<T>::axby(Vector<LevelData<T>* >& a_lhs,
00545 const Vector<LevelData<T>* >& a_x,
00546 const Vector<LevelData<T>* >& a_y,
00547 Real a_a,
00548 Real a_b)
00549 {
00550 CH_TIME("MultilevelLinearOp::axby");
00551
00552
00553 for (int level=m_lBase; level<a_lhs.size(); ++level)
00554 {
00555 m_vectOperators[level]->axby(*a_lhs[level],
00556 *a_x[level],
00557 *a_y[level],
00558 a_a, a_b);
00559 }
00560 }
00561
00562
00563
00564
00565
00566 template<class T>
00567 void
00568 MultilevelLinearOp<T>::scale(Vector<LevelData<T>* >& a_lhs,
00569 const Real& a_scale)
00570 {
00571 CH_TIME("MultilevelLinearOp::scale");
00572
00573
00574 for (int level =m_lBase; level<a_lhs.size(); level++)
00575 {
00576 m_vectOperators[level]->scale(*a_lhs[level], a_scale);
00577 }
00578
00579 }
00580
00581
00582
00583
00584
00585
00586 template<class T>
00587 Real
00588 MultilevelLinearOp<T>::norm(const Vector<LevelData<T>* >& a_rhs,
00589 int a_ord)
00590 {
00591 CH_TIME("MultilevelLinearOp::norm");
00592
00593
00594 int maxLevel = a_rhs.size()-1;
00595 Real thisNorm = 0.0;
00596 for (int level = m_lBase; level<a_rhs.size(); level++)
00597 {
00598 Real normLevel;
00599
00600 if (level < maxLevel)
00601 {
00602
00603 int refRatio = m_refRatios[level];
00604 normLevel = m_vectOperators[level]->AMRNorm(*a_rhs[level],
00605 *a_rhs[level+1],
00606 refRatio,
00607 a_ord);
00608 }
00609 else
00610 {
00611
00612 LevelData<T> tempLDF;
00613 int refRatio = -1;
00614 normLevel = m_vectOperators[level]->AMRNorm(*a_rhs[level],
00615 tempLDF,
00616 refRatio,
00617 a_ord);
00618 }
00619
00620
00621 if (a_ord > 2)
00622 {
00623 MayDay::Error("MultilevelLinearOp::norm -- undefined for order > 2");
00624 }
00625
00626 if (a_ord == 2)
00627 {
00628 normLevel = normLevel*normLevel;
00629 }
00630
00631 if (a_ord != 0)
00632 {
00633 thisNorm += normLevel;
00634 }
00635 else if (normLevel > thisNorm)
00636 {
00637 thisNorm = normLevel;
00638 }
00639 }
00640
00641 if (a_ord == 2)
00642 {
00643 thisNorm = sqrt(thisNorm);
00644 }
00645
00646 return thisNorm;
00647 }
00648
00649
00650
00651
00652
00653 template<class T>
00654 void
00655 MultilevelLinearOp<T>::setToZero(Vector<LevelData<T>* >& a_lhs)
00656 {
00657 CH_TIME("MultilevelLinearOp::setToZero");
00658
00659
00660
00661 int coarsestLevel = Max(0, m_lBase-1);
00662 for (int level = coarsestLevel; level< a_lhs.size(); level++)
00663 {
00664 m_vectOperators[level]->setToZero(*a_lhs[level]);
00665 }
00666 }
00667
00668 template<class T>
00669 void
00670 MultilevelLinearOp<T>::write(const Vector<LevelData<T>* >* a_data,
00671 const char* a_filename)
00672 {
00673 #ifdef CH_USE_HDF5
00674 writeVectorLevelName(a_data, &m_refRatios, a_filename);
00675 #else
00676 MayDay::Warning("MultilevelLinearOp<T>::write unimplemented since CH_USE_HDF5 undefined");
00677 #endif
00678 }
00679
00680 #include "NamespaceFooter.H"
00681
00682 #endif