10#ifndef MUELU_HIERARCHY_DEF_HPP
11#define MUELU_HIERARCHY_DEF_HPP
19#include <Xpetra_Matrix.hpp>
20#include <Xpetra_MultiVectorFactory.hpp>
21#include <Xpetra_Operator.hpp>
22#include <Xpetra_IO.hpp>
26#include "MueLu_FactoryManager.hpp"
27#include "MueLu_HierarchyUtils.hpp"
28#include "MueLu_TopRAPFactory.hpp"
29#include "MueLu_TopSmootherFactory.hpp"
32#include "MueLu_PerfUtils.hpp"
33#include "MueLu_PFactory.hpp"
34#include "MueLu_SmootherFactory.hpp"
38#include "Teuchos_TimeMonitor.hpp"
42template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
44 : maxCoarseSize_(GetDefaultMaxCoarseSize())
45 , implicitTranspose_(GetDefaultImplicitTranspose())
46 , fuseProlongationAndUpdate_(GetDefaultFuseProlongationAndUpdate())
47 , doPRrebalance_(GetDefaultPRrebalance())
48 , doPRViaCopyrebalance_(false)
49 , isPreconditioner_(true)
50 , Cycle_(GetDefaultCycle())
51 , WCycleStartLevel_(0)
52 , scalingFactor_(
Teuchos::ScalarTraits<double>::one())
54 , isDumpingEnabled_(false)
57 , sizeOfAllocatedLevelMultiVectors_(0) {
61template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
65 setObjectLabel(label);
66 Levels_[0]->setObjectLabel(label);
69template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
71 : maxCoarseSize_(GetDefaultMaxCoarseSize())
72 , implicitTranspose_(GetDefaultImplicitTranspose())
73 , fuseProlongationAndUpdate_(GetDefaultFuseProlongationAndUpdate())
74 , doPRrebalance_(GetDefaultPRrebalance())
75 , doPRViaCopyrebalance_(false)
76 , isPreconditioner_(true)
77 , Cycle_(GetDefaultCycle())
78 , WCycleStartLevel_(0)
79 , scalingFactor_(
Teuchos::ScalarTraits<double>::one())
80 , isDumpingEnabled_(false)
83 , sizeOfAllocatedLevelMultiVectors_(0) {
84 lib_ = A->getDomainMap()->lib();
86 RCP<Level> Finest = rcp(
new Level);
92template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
95 setObjectLabel(label);
96 Levels_[0]->setObjectLabel(label);
99template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
102template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
104 int levelID = LastLevelID() + 1;
106 if (level->GetLevelID() != -1 && (level->GetLevelID() != levelID))
107 GetOStream(
Warnings1) <<
"Hierarchy::AddLevel(): Level with ID=" << level->GetLevelID() <<
" have been added at the end of the hierarchy\n but its ID have been redefined"
108 <<
" because last level ID of the hierarchy was " << LastLevelID() <<
"." << std::endl;
110 Levels_.push_back(level);
111 level->SetLevelID(levelID);
114 level->SetPreviousLevel((levelID == 0) ? Teuchos::null : Levels_[LastLevelID() - 1]);
115 level->setObjectLabel(this->getObjectLabel());
118template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
120 RCP<Level> newLevel = Levels_[LastLevelID()]->Build();
121 newLevel->setlib(lib_);
122 this->AddLevel(newLevel);
125template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
128 "MueLu::Hierarchy::GetLevel(): invalid input parameter value: LevelID = " << levelID);
129 return Levels_[levelID];
132template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134 return Levels_.size();
137template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
139 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator>>(
"A");
140 RCP<const Teuchos::Comm<int>> comm = A->getDomainMap()->getComm();
142 int numLevels = GetNumLevels();
144 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, numLevels, Teuchos::ptr(&numGlobalLevels));
146 return numGlobalLevels;
149template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
151 double totalNnz = 0, lev0Nnz = 1;
152 for (
int i = 0; i < GetNumLevels(); ++i) {
154 "Operator complexity cannot be calculated because A is unavailable on level " << i);
155 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator>>(
"A");
159 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
161 GetOStream(
Warnings0) <<
"Some level operators are not matrices, operator complexity calculation aborted" << std::endl;
165 totalNnz += as<double>(Am->getGlobalNumEntries());
169 return totalNnz / lev0Nnz;
172template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
174 double node_sc = 0, global_sc = 0;
176 const size_t INVALID = Teuchos::OrdinalTraits<size_t>::invalid();
178 if (GetNumLevels() <= 0)
return -1.0;
179 if (!Levels_[0]->IsAvailable(
"A"))
return -1.0;
181 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator>>(
"A");
182 if (A.is_null())
return -1.0;
183 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
184 if (Am.is_null())
return -1.0;
185 a0_nnz = as<double>(Am->getGlobalNumEntries());
188 for (
int i = 0; i < GetNumLevels(); ++i) {
190 if (!Levels_[i]->IsAvailable(
"PreSmoother"))
continue;
191 RCP<SmootherBase> S = Levels_[i]->template Get<RCP<SmootherBase>>(
"PreSmoother");
192 if (S.is_null())
continue;
193 level_sc = S->getNodeSmootherComplexity();
194 if (level_sc == INVALID) {
199 node_sc += as<double>(level_sc);
203 RCP<const Teuchos::Comm<int>> comm = A->getDomainMap()->getComm();
204 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, node_sc, Teuchos::ptr(&global_sc));
205 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MIN, node_sc, Teuchos::ptr(&min_sc));
210 return global_sc / a0_nnz;
214template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
217 "MueLu::Hierarchy::CheckLevel(): wrong underlying linear algebra library.");
219 "MueLu::Hierarchy::CheckLevel(): wrong level ID");
221 "MueLu::Hierarchy::Setup(): wrong level parent");
224template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
226 for (
int i = 0; i < GetNumLevels(); ++i) {
227 RCP<Level> level = Levels_[i];
228 if (level->IsAvailable(
"A")) {
229 RCP<Operator> Aop = level->Get<RCP<Operator>>(
"A");
230 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Aop);
232 RCP<const Import> xpImporter = A->getCrsGraph()->getImporter();
233 if (!xpImporter.is_null())
234 xpImporter->setDistributorParameters(matvecParams);
235 RCP<const Export> xpExporter = A->getCrsGraph()->getExporter();
236 if (!xpExporter.is_null())
237 xpExporter->setDistributorParameters(matvecParams);
240 const std::list<std::string> matrices = {
"P",
"R",
"D0",
"NodeMatrix"};
241 for (
auto it = matrices.begin(); it != matrices.end(); ++it) {
242 if (level->IsAvailable(*it)) {
243 RCP<Matrix> mat = level->Get<RCP<Matrix>>(*it);
244 if (!mat.is_null()) {
245 RCP<const Import> xpImporter = mat->getCrsGraph()->getImporter();
246 if (!xpImporter.is_null()) {
247 xpImporter->setDistributorParameters(matvecParams);
249 RCP<const Export> xpExporter = mat->getCrsGraph()->getExporter();
250 if (!xpExporter.is_null())
251 xpExporter->setDistributorParameters(matvecParams);
255 if (level->IsAvailable(
"Importer")) {
256 RCP<const Import> xpImporter = level->Get<RCP<const Import>>(
"Importer");
257 if (!xpImporter.is_null())
258 xpImporter->setDistributorParameters(matvecParams);
265template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
267 const RCP<const FactoryManagerBase> fineLevelManager,
268 const RCP<const FactoryManagerBase> coarseLevelManager,
269 const RCP<const FactoryManagerBase> nextLevelManager) {
274 "MueLu::Hierarchy:Setup(): level " << coarseLevelID <<
" (specified by coarseLevelID argument) "
275 "must be built before calling this function.");
277 Level& level = *Levels_[coarseLevelID];
279 bool useStackedTimer = !Teuchos::TimeMonitor::stackedTimerNameIsDefault();
283 if (!useStackedTimer)
284 m1 = rcp(
new TimeMonitor(*
this, label + this->ShortClassName() +
": " +
"Setup (total)"));
285 TimeMonitor m2(*
this, label + this->ShortClassName() +
": " +
"Setup" +
" (total, level=" + Teuchos::toString(coarseLevelID) +
")");
289 "MueLu::Hierarchy::Setup(): argument coarseLevelManager cannot be null");
294 if (levelManagers_.size() < coarseLevelID + 1)
295 levelManagers_.resize(coarseLevelID + 1);
296 levelManagers_[coarseLevelID] = coarseLevelManager;
298 bool isFinestLevel = (fineLevelManager.is_null());
299 bool isLastLevel = (nextLevelManager.is_null());
303 RCP<Operator> A = level.
Get<RCP<Operator>>(
"A");
304 RCP<const Map> domainMap = A->getDomainMap();
305 RCP<const Teuchos::Comm<int>> comm = domainMap->getComm();
312 oldRank = SetProcRankVerbose(comm->getRank());
316 lib_ = domainMap->lib();
323 Level& prevLevel = *Levels_[coarseLevelID - 1];
324 oldRank = SetProcRankVerbose(prevLevel.
GetComm()->getRank());
327 CheckLevel(level, coarseLevelID);
330 RCP<SetFactoryManager> SFMFine;
332 SFMFine = rcp(
new SetFactoryManager(Levels_[coarseLevelID - 1], fineLevelManager));
334 if (isFinestLevel && Levels_[coarseLevelID]->IsAvailable(
"Coordinates"))
335 ReplaceCoordinateMap(*Levels_[coarseLevelID]);
340 if (isDumpingEnabled_ && (dumpLevel_ == 0 || dumpLevel_ == -1) && coarseLevelID == 1)
343 RCP<TopSmootherFactory> coarseFact;
344 RCP<TopSmootherFactory> smootherFact = rcp(
new TopSmootherFactory(coarseLevelManager,
"Smoother"));
346 int nextLevelID = coarseLevelID + 1;
348 RCP<SetFactoryManager> SFMNext;
349 if (isLastLevel ==
false) {
351 if (nextLevelID > LastLevelID())
353 CheckLevel(*Levels_[nextLevelID], nextLevelID);
357 Levels_[nextLevelID]->Request(
TopRAPFactory(coarseLevelManager, nextLevelManager));
386 if (coarseFact.is_null())
395 RCP<Operator> Ac = Teuchos::null;
396 TopRAPFactory coarseRAPFactory(fineLevelManager, coarseLevelManager);
399 Ac = level.
Get<RCP<Operator>>(
"A");
400 }
else if (!isFinestLevel) {
405 bool setLastLevelviaMaxCoarseSize =
false;
407 Ac = level.
Get<RCP<Operator>>(
"A");
408 RCP<Matrix> Acm = rcp_dynamic_cast<Matrix>(Ac);
412 level.
SetComm(Ac->getDomainMap()->getComm());
415 bool isOrigLastLevel = isLastLevel;
420 }
else if (Ac.is_null()) {
427 if (!Acm.is_null() && Acm->getGlobalNumRows() <= maxCoarseSize_) {
429 GetOStream(
Runtime0) <<
"Max coarse size (<= " << maxCoarseSize_ <<
") achieved" << std::endl;
431 if (Acm->getGlobalNumRows() != 0) setLastLevelviaMaxCoarseSize =
true;
435 if (!Ac.is_null() && !isFinestLevel) {
436 RCP<Operator> A = Levels_[coarseLevelID - 1]->template Get<RCP<Operator>>(
"A");
437 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
439 const double maxCoarse2FineRatio = 0.8;
440 if (!Acm.is_null() && !Am.is_null() && Acm->getGlobalNumRows() > maxCoarse2FineRatio * Am->getGlobalNumRows()) {
448 GetOStream(
Warnings0) <<
"Aggregation stagnated. Please check your matrix and/or adjust your configuration file."
449 <<
"Possible fixes:\n"
450 <<
" - reduce the maximum number of levels\n"
451 <<
" - enable repartitioning\n"
452 <<
" - increase the minimum coarse size." << std::endl;
457 if (!isOrigLastLevel) {
461 if (coarseFact.is_null())
469 coarseFact->Build(level);
480 smootherFact->Build(level);
485 if (isLastLevel ==
true) {
486 int actualNumLevels = nextLevelID;
487 if (isOrigLastLevel ==
false) {
490 Levels_[nextLevelID]->Release(
TopRAPFactory(coarseLevelManager, nextLevelManager));
496 if (!setLastLevelviaMaxCoarseSize) {
497 if (Levels_[nextLevelID - 1]->IsAvailable(
"P")) {
498 if (Levels_[nextLevelID - 1]->
template Get<RCP<Matrix>>(
"P") == Teuchos::null) actualNumLevels = nextLevelID - 1;
500 actualNumLevels = nextLevelID - 1;
503 if (actualNumLevels == nextLevelID - 1) {
505 Levels_[nextLevelID - 2]->Release(*smootherFact);
507 if (Levels_[nextLevelID - 2]->IsAvailable(
"PreSmoother")) Levels_[nextLevelID - 2]->RemoveKeepFlag(
"PreSmoother",
NoFactory::get());
508 if (Levels_[nextLevelID - 2]->IsAvailable(
"PostSmoother")) Levels_[nextLevelID - 2]->RemoveKeepFlag(
"PostSmoother",
NoFactory::get());
509 if (coarseFact.is_null())
511 Levels_[nextLevelID - 2]->Request(*coarseFact);
512 if (!(Levels_[nextLevelID - 2]->
template Get<RCP<Matrix>>(
"A").is_null()))
513 coarseFact->Build(*(Levels_[nextLevelID - 2]));
514 Levels_[nextLevelID - 2]->Release(*coarseFact);
516 Levels_.resize(actualNumLevels);
517 levelManagers_.resize(actualNumLevels);
521 if (isDumpingEnabled_ && ((dumpLevel_ > 0 && coarseLevelID == dumpLevel_) || dumpLevel_ == -1))
522 DumpCurrentGraph(coarseLevelID);
524 if (!isFinestLevel) {
528 level.
Release(coarseRAPFactory);
532 SetProcRankVerbose(oldRank);
537template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
539 int numLevels = Levels_.size();
541 "Hierarchy::SetupRe: " << Levels_.size() <<
" levels, but " << levelManagers_.size() <<
" level factory managers");
543 const int startLevel = 0;
546#ifdef HAVE_MUELU_DEBUG
548 for (
int i = 0; i < numLevels; i++)
549 levelManagers_[i]->ResetDebugData();
554 for (levelID = startLevel; levelID < numLevels;) {
555 bool r = Setup(levelID,
556 (levelID != 0 ? levelManagers_[levelID - 1] : Teuchos::null),
557 levelManagers_[levelID],
558 (levelID + 1 != numLevels ? levelManagers_[levelID + 1] : Teuchos::null));
564 Levels_.resize(levelID);
565 levelManagers_.resize(levelID);
567 int sizeOfVecs = sizeOfAllocatedLevelMultiVectors_;
569 AllocateLevelMultiVectors(sizeOfVecs,
true);
576 CheckForEmptySmoothersAndCoarseSolve();
579template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
588 "MueLu::Hierarchy::Setup(): fine level (" << startLevel <<
") does not exist");
591 "Constructing non-positive (" << numDesiredLevels <<
") number of levels does not make sense.");
595 "MueLu::Hierarchy::Setup(): fine level (" << startLevel <<
") has no matrix A! "
596 "Set fine level matrix A using Level.Set()");
598 RCP<Operator> A = Levels_[startLevel]->template Get<RCP<Operator>>(
"A");
599 lib_ = A->getDomainMap()->lib();
602 RCP<Matrix> Amat = rcp_dynamic_cast<Matrix>(A);
604 if (!Amat.is_null()) {
605 RCP<ParameterList> params = rcp(
new ParameterList());
606 params->set(
"printLoadBalancingInfo",
true);
607 params->set(
"printCommInfo",
true);
611 GetOStream(
Warnings1) <<
"Fine level operator is not a matrix, statistics are not available" << std::endl;
615 RCP<const FactoryManagerBase> rcpmanager = rcpFromRef(manager);
617 const int lastLevel = startLevel + numDesiredLevels - 1;
618 GetOStream(
Runtime0) <<
"Setup loop: startLevel = " << startLevel <<
", lastLevel = " << lastLevel
619 <<
" (stop if numLevels = " << numDesiredLevels <<
" or Ac.size() < " << maxCoarseSize_ <<
")" << std::endl;
623 if (numDesiredLevels == 1) {
625 Setup(startLevel, Teuchos::null, rcpmanager, Teuchos::null);
628 bool bIsLastLevel = Setup(startLevel, Teuchos::null, rcpmanager, rcpmanager);
629 if (bIsLastLevel ==
false) {
630 for (iLevel = startLevel + 1; iLevel < lastLevel; iLevel++) {
631 bIsLastLevel = Setup(iLevel, rcpmanager, rcpmanager, rcpmanager);
632 if (bIsLastLevel ==
true)
635 if (bIsLastLevel ==
false)
636 Setup(lastLevel, rcpmanager, rcpmanager, Teuchos::null);
642 "MueLu::Hierarchy::Setup(): number of level");
650 CheckForEmptySmoothersAndCoarseSolve();
653template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
655 for (LO levelNo = 0; levelNo < GetNumLevels(); ++levelNo) {
656 auto level = Levels_[levelNo];
657 if ((level->IsAvailable(
"A") && !level->template Get<RCP<Operator>>(
"A").is_null()) && (!level->IsAvailable(
"PreSmoother")) && (!level->IsAvailable(
"PostSmoother"))) {
658 GetOStream(
Warnings1) <<
"No " << (levelNo == as<LO>(Levels_.size()) - 1 ?
"coarse grid solver" :
"smoother") <<
" on level " << level->GetLevelID() << std::endl;
663template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
665 if (startLevel < GetNumLevels())
666 GetOStream(
Runtime0) <<
"Clearing old data (if any)" << std::endl;
668 for (
int iLevel = startLevel; iLevel < GetNumLevels(); iLevel++)
669 Levels_[iLevel]->Clear();
672template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
674 GetOStream(
Runtime0) <<
"Clearing old data (expert)" << std::endl;
675 for (
int iLevel = 0; iLevel < GetNumLevels(); iLevel++)
676 Levels_[iLevel]->ExpertClear();
679#if defined(HAVE_MUELU_EXPERIMENTAL) && defined(HAVE_MUELU_ADDITIVE_VARIANT)
680template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
682 bool InitialGuessIsZero, LO startLevel) {
683 LO nIts = conv.maxIts_;
684 MagnitudeType tol = conv.tol_;
686 std::string prefix = this->ShortClassName() +
": ";
687 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
688 std::string levelSuffix1 =
" (level=" +
toString(startLevel + 1) +
")";
691 RCP<Time> CompTime = Teuchos::TimeMonitor::getNewCounter(prefix +
"Computational Time (total)");
692 RCP<Time> Concurrent = Teuchos::TimeMonitor::getNewCounter(prefix +
"Concurrent portion");
693 RCP<Time> ApplyR = Teuchos::TimeMonitor::getNewCounter(prefix +
"R: Computational Time");
694 RCP<Time> ApplyPbar = Teuchos::TimeMonitor::getNewCounter(prefix +
"Pbar: Computational Time");
695 RCP<Time> CompFine = Teuchos::TimeMonitor::getNewCounter(prefix +
"Fine: Computational Time");
696 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix +
"Coarse: Computational Time");
697 RCP<Time> ApplySum = Teuchos::TimeMonitor::getNewCounter(prefix +
"Sum: Computational Time");
698 RCP<Time> Synchronize_beginning = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_beginning");
699 RCP<Time> Synchronize_center = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_center");
700 RCP<Time> Synchronize_end = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_end");
702 RCP<Level> Fine = Levels_[0];
705 RCP<Operator> A = Fine->Get<RCP<Operator>>(
"A");
706 Teuchos::RCP<const Teuchos::Comm<int>> communicator = A->getDomainMap()->getComm();
714 SC one = STS::one(), zero = STS::zero();
716 bool zeroGuess = InitialGuessIsZero;
724 RCP<MultiVector> coarseRhs, coarseX;
726 RCP<SmootherBase> preSmoo_coarse, postSmoo_coarse;
727 bool emptyCoarseSolve =
true;
728 RCP<MultiVector> coarseX_prolonged = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
730 RCP<const Import> importer;
732 if (Levels_.size() > 1) {
734 if (Coarse->IsAvailable(
"Importer"))
735 importer = Coarse->Get<RCP<const Import>>(
"Importer");
737 R = Coarse->Get<RCP<Operator>>(
"R");
738 P = Coarse->Get<RCP<Operator>>(
"P");
741 Pbar = Coarse->Get<RCP<Operator>>(
"Pbar");
743 coarseRhs = MultiVectorFactory::Build(R->getRangeMap(), B.getNumVectors(),
true);
745 Ac = Coarse->Get<RCP<Operator>>(
"A");
748 R->apply(B, *coarseRhs, Teuchos::NO_TRANS, one, zero);
752 if (doPRrebalance_ || importer.is_null()) {
753 coarseX = MultiVectorFactory::Build(coarseRhs->getMap(), X.getNumVectors(),
true);
756 RCP<TimeMonitor> ITime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import (total)",
Timings0));
757 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import" + levelSuffix1,
Timings0));
760 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getTargetMap(), coarseRhs->getNumVectors());
761 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
762 coarseRhs.swap(coarseTmp);
764 coarseX = MultiVectorFactory::Build(importer->getTargetMap(), X.getNumVectors(),
true);
767 if (Coarse->IsAvailable(
"PreSmoother"))
768 preSmoo_coarse = Coarse->Get<RCP<SmootherBase>>(
"PreSmoother");
769 if (Coarse->IsAvailable(
"PostSmoother"))
770 postSmoo_coarse = Coarse->Get<RCP<SmootherBase>>(
"PostSmoother");
775 MagnitudeType prevNorm = STS::magnitude(STS::one()), curNorm = STS::magnitude(STS::one());
779 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
780 std::ostringstream ss;
781 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
782 throw Exceptions::Incompatible(ss.str());
785 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
786 std::ostringstream ss;
787 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
788 throw Exceptions::Incompatible(ss.str());
792 bool emptyFineSolve =
true;
794 RCP<MultiVector> fineX;
795 fineX = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
804 if (Fine->IsAvailable(
"PreSmoother")) {
805 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>(
"PreSmoother");
807 preSmoo->Apply(*fineX, B, zeroGuess);
809 emptyFineSolve =
false;
811 if (Fine->IsAvailable(
"PostSmoother")) {
812 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>(
"PostSmoother");
814 postSmoo->Apply(*fineX, B, zeroGuess);
817 emptyFineSolve =
false;
819 if (emptyFineSolve ==
true) {
821 fineX->update(one, B, zero);
824 if (Levels_.size() > 1) {
826 if (Coarse->IsAvailable(
"PreSmoother")) {
828 preSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
830 emptyCoarseSolve =
false;
832 if (Coarse->IsAvailable(
"PostSmoother")) {
834 postSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
836 emptyCoarseSolve =
false;
838 if (emptyCoarseSolve ==
true) {
840 coarseX->update(one, *coarseRhs, zero);
847 if (!doPRrebalance_ && !importer.is_null()) {
848 RCP<TimeMonitor> ITime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export (total)",
Timings0));
849 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export" + levelSuffix1,
Timings0));
852 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getSourceMap(), coarseX->getNumVectors());
853 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
854 coarseX.swap(coarseTmp);
858 Pbar->apply(*coarseX, *coarseX_prolonged, Teuchos::NO_TRANS, one, zero);
863 X.update(1.0, *fineX, 1.0, *coarseX_prolonged, 0.0);
874template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
876 bool InitialGuessIsZero, LO startLevel) {
892 RCP<Level> Fine = Levels_[startLevel];
895 std::string prefix = label + this->ShortClassName() +
": ";
896 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
897 std::string levelSuffix1 =
" (level=" +
toString(startLevel + 1) +
")";
899 bool useStackedTimer = !Teuchos::TimeMonitor::stackedTimerNameIsDefault();
901 RCP<Monitor> iterateTime;
902 RCP<TimeMonitor> iterateTime1;
905 else if (!useStackedTimer)
908 std::string iterateLevelTimeLabel = prefix +
"Solve" + levelSuffix;
909 RCP<TimeMonitor> iterateLevelTime = rcp(
new TimeMonitor(*
this, iterateLevelTimeLabel,
Timings0));
911 bool zeroGuess = InitialGuessIsZero;
913 RCP<Operator> A = Fine->Get<RCP<Operator>>(
"A");
915 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix +
"Coarse: Computational Time");
928 const BlockedMultiVector* Bblocked =
dynamic_cast<const BlockedMultiVector*
>(&B);
929 if (residual_.size() > startLevel &&
930 ((Bblocked && !Bblocked->isSameSize(*residual_[startLevel])) ||
931 (!Bblocked && !residual_[startLevel]->isSameSize(B))))
932 DeleteLevelMultiVectors();
933 AllocateLevelMultiVectors(X.getNumVectors());
936 typedef Teuchos::ScalarTraits<typename STS::magnitudeType> STM;
939 if (IsCalculationOfResidualRequired(startLevel, conv)) {
940 ConvergenceStatus convergenceStatus = ComputeResidualAndPrintHistory(*A, X, B, Teuchos::ScalarTraits<LO>::zero(), startLevel, conv, prevNorm);
942 return convergenceStatus;
945 SC one = STS::one(), zero = STS::zero();
946 for (LO iteration = 1; iteration <= nIts; iteration++) {
947#ifdef HAVE_MUELU_DEBUG
949 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
950 std::ostringstream ss;
951 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
955 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
956 std::ostringstream ss;
957 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
963 if (startLevel == as<LO>(Levels_.size()) - 1) {
965 RCP<TimeMonitor> CLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : coarse" + levelSuffix,
Timings0));
967 bool emptySolve =
true;
970 if (Fine->IsAvailable(
"PreSmoother")) {
971 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>(
"PreSmoother");
973 preSmoo->Apply(X, B, zeroGuess);
978 if (Fine->IsAvailable(
"PostSmoother")) {
979 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>(
"PostSmoother");
981 postSmoo->Apply(X, B, zeroGuess);
986 if (emptySolve ==
true) {
988 X.update(one, B, zero);
993 RCP<Level> Coarse = Levels_[startLevel + 1];
996 RCP<TimeMonitor> STime;
997 if (!useStackedTimer)
999 RCP<TimeMonitor> SLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : smoothing" + levelSuffix,
Timings0));
1001 if (Fine->IsAvailable(
"PreSmoother")) {
1002 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>(
"PreSmoother");
1003 preSmoo->Apply(X, B, zeroGuess);
1008 RCP<MultiVector> residual;
1010 RCP<TimeMonitor> ATime;
1011 if (!useStackedTimer)
1012 ATime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : residual calculation (total)",
Timings0));
1013 RCP<TimeMonitor> ALevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : residual calculation" + levelSuffix,
Timings0));
1021 residual = residual_[startLevel];
1024 RCP<Operator> P = Coarse->Get<RCP<Operator>>(
"P");
1025 if (Coarse->IsAvailable(
"Pbar"))
1026 P = Coarse->Get<RCP<Operator>>(
"Pbar");
1028 RCP<MultiVector> coarseRhs, coarseX;
1032 RCP<TimeMonitor> RTime;
1033 if (!useStackedTimer)
1035 RCP<TimeMonitor> RLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : restriction" + levelSuffix,
Timings0));
1036 coarseRhs = coarseRhs_[startLevel];
1038 if (implicitTranspose_) {
1039 P->apply(*residual, *coarseRhs, Teuchos::TRANS, one, zero);
1042 RCP<Operator> R = Coarse->Get<RCP<Operator>>(
"R");
1043 R->apply(*residual, *coarseRhs, Teuchos::NO_TRANS, one, zero);
1047 RCP<const Import> importer;
1048 if (Coarse->IsAvailable(
"Importer"))
1049 importer = Coarse->Get<RCP<const Import>>(
"Importer");
1051 coarseX = coarseX_[startLevel];
1052 if (!doPRrebalance_ && !importer.is_null()) {
1053 RCP<TimeMonitor> ITime;
1054 if (!useStackedTimer)
1056 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import" + levelSuffix1,
Timings0));
1059 RCP<MultiVector> coarseTmp = coarseImport_[startLevel];
1060 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
1061 coarseRhs.swap(coarseTmp);
1064 RCP<Operator> Ac = Coarse->Get<RCP<Operator>>(
"A");
1065 if (!Ac.is_null()) {
1066 RCP<const Map> origXMap = coarseX->getMap();
1067 RCP<const Map> origRhsMap = coarseRhs->getMap();
1070 coarseRhs->replaceMap(Ac->getRangeMap());
1071 coarseX->replaceMap(Ac->getDomainMap());
1074 iterateLevelTime = Teuchos::null;
1076 Iterate(*coarseRhs, *coarseX, 1,
true, startLevel + 1);
1078 if (Cycle_ ==
WCYCLE && WCycleStartLevel_ <= startLevel)
1079 Iterate(*coarseRhs, *coarseX, 1,
false, startLevel + 1);
1082 iterateLevelTime = rcp(
new TimeMonitor(*
this, iterateLevelTimeLabel));
1084 coarseX->replaceMap(origXMap);
1085 coarseRhs->replaceMap(origRhsMap);
1088 if (!doPRrebalance_ && !importer.is_null()) {
1089 RCP<TimeMonitor> ITime;
1090 if (!useStackedTimer)
1092 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export" + levelSuffix1,
Timings0));
1095 RCP<MultiVector> coarseTmp = coarseExport_[startLevel];
1096 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
1097 coarseX.swap(coarseTmp);
1102 RCP<TimeMonitor> PTime;
1103 if (!useStackedTimer)
1105 RCP<TimeMonitor> PLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : prolongation" + levelSuffix,
Timings0));
1110 if (fuseProlongationAndUpdate_) {
1111 P->apply(*coarseX, X, Teuchos::NO_TRANS, scalingFactor_, one);
1113 RCP<MultiVector> correction = correction_[startLevel];
1114 P->apply(*coarseX, *correction, Teuchos::NO_TRANS, one, zero);
1115 X.update(scalingFactor_, *correction, one);
1121 RCP<TimeMonitor> STime;
1122 if (!useStackedTimer)
1124 RCP<TimeMonitor> SLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : smoothing" + levelSuffix,
Timings0));
1126 if (Fine->IsAvailable(
"PostSmoother")) {
1127 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>(
"PostSmoother");
1128 postSmoo->Apply(X, B,
false);
1134 if (IsCalculationOfResidualRequired(startLevel, conv)) {
1135 ConvergenceStatus convergenceStatus = ComputeResidualAndPrintHistory(*A, X, B, iteration, startLevel, conv, prevNorm);
1137 return convergenceStatus;
1144template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1146 LO startLevel = (start != -1 ? start : 0);
1147 LO endLevel = (end != -1 ? end : Levels_.size() - 1);
1150 "MueLu::Hierarchy::Write : startLevel must be <= endLevel");
1153 "MueLu::Hierarchy::Write bad start or end level");
1155 for (LO i = startLevel; i < endLevel + 1; i++) {
1156 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Levels_[i]->
template Get<RCP<Operator>>(
"A")), P, R;
1158 P = rcp_dynamic_cast<Matrix>(Levels_[i]->
template Get<RCP<Operator>>(
"P"));
1159 if (!implicitTranspose_)
1160 R = rcp_dynamic_cast<Matrix>(Levels_[i]->
template Get<RCP<Operator>>(
"R"));
1163 if (!A.is_null()) Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"A_" +
toString(i) + suffix +
".m", *A);
1165 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"P_" +
toString(i) + suffix +
".m", *P);
1168 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"R_" +
toString(i) + suffix +
".m", *R);
1173template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1175 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1176 (*it)->Keep(ename, factory);
1179template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1181 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1182 (*it)->Delete(ename, factory);
1185template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1187 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1188 (*it)->AddKeepFlag(ename, factory, keep);
1191template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1193 for (Array<RCP<Level>>::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1194 (*it)->RemoveKeepFlag(ename, factory, keep);
1197template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1199 if (description_.empty()) {
1200 std::ostringstream out;
1202 out <<
"{#levels = " << GetGlobalNumLevels() <<
", complexity = " << GetOperatorComplexity() <<
"}";
1203 description_ = out.str();
1205 return description_;
1208template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1213template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1215 RCP<Operator> A0 = Levels_[0]->template Get<RCP<Operator>>(
"A");
1216 RCP<const Teuchos::Comm<int>> comm = A0->getDomainMap()->getComm();
1218 int numLevels = GetNumLevels();
1219 RCP<Operator> Ac = Levels_[numLevels - 1]->template Get<RCP<Operator>>(
"A");
1226 int root = comm->getRank();
1229 int smartData = numLevels * comm->getSize() + comm->getRank(), maxSmartData;
1230 reduceAll(*comm, Teuchos::REDUCE_MAX, smartData, Teuchos::ptr(&maxSmartData));
1231 root = maxSmartData % comm->getSize();
1235 double smoother_comp = -1.0;
1237 smoother_comp = GetSmootherComplexity();
1241 std::vector<Xpetra::global_size_t> nnzPerLevel;
1242 std::vector<Xpetra::global_size_t> rowsPerLevel;
1243 std::vector<int> numProcsPerLevel;
1244 bool someOpsNotMatrices =
false;
1245 const Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1246 for (
int i = 0; i < numLevels; i++) {
1248 "Operator A is not available on level " << i);
1250 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator>>(
"A");
1252 "Operator A on level " << i <<
" is null.");
1254 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
1256 someOpsNotMatrices =
true;
1257 nnzPerLevel.push_back(INVALID);
1258 rowsPerLevel.push_back(A->getDomainMap()->getGlobalNumElements());
1259 numProcsPerLevel.push_back(A->getDomainMap()->getComm()->getSize());
1261 LO storageblocksize = Am->GetStorageBlockSize();
1262 Xpetra::global_size_t nnz = Am->getGlobalNumEntries() * storageblocksize * storageblocksize;
1263 nnzPerLevel.push_back(nnz);
1264 rowsPerLevel.push_back(Am->getGlobalNumRows() * storageblocksize);
1265 numProcsPerLevel.push_back(Am->getRowMap()->getComm()->getSize());
1268 if (someOpsNotMatrices)
1269 GetOStream(
Warnings0) <<
"Some level operators are not matrices, statistics calculation are incomplete" << std::endl;
1272 std::string label = Levels_[0]->getObjectLabel();
1273 std::ostringstream oss;
1274 oss << std::setfill(
' ');
1275 oss <<
"\n--------------------------------------------------------------------------------\n";
1276 oss <<
"--- Multigrid Summary " << std::setw(32) <<
"---\n";
1277 oss <<
"--------------------------------------------------------------------------------" << std::endl;
1278 if (hierarchyLabel_ !=
"") oss <<
"Label = " << hierarchyLabel_ << std::endl;
1280 oss <<
"Scalar = " << Teuchos::ScalarTraits<Scalar>::name() << std::endl;
1281 oss <<
"Number of levels = " << numLevels << std::endl;
1282 oss <<
"Operator complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1283 if (!someOpsNotMatrices)
1284 oss << GetOperatorComplexity() << std::endl;
1286 oss <<
"not available (Some operators in hierarchy are not matrices.)" << std::endl;
1288 if (smoother_comp != -1.0) {
1289 oss <<
"Smoother complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed)
1290 << smoother_comp << std::endl;
1295 oss <<
"Cycle type = V" << std::endl;
1298 oss <<
"Cycle type = W" << std::endl;
1299 if (WCycleStartLevel_ > 0)
1300 oss <<
"Cycle start level = " << WCycleStartLevel_ << std::endl;
1307 Xpetra::global_size_t tt = rowsPerLevel[0];
1313 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1314 tt = nnzPerLevel[i];
1324 tt = numProcsPerLevel[0];
1330 oss <<
"level " << std::setw(rowspacer) <<
" rows " << std::setw(nnzspacer) <<
" nnz "
1331 <<
" nnz/row" << std::setw(npspacer) <<
" c ratio"
1332 <<
" procs" << std::endl;
1333 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1334 oss <<
" " << i <<
" ";
1335 oss << std::setw(rowspacer) << rowsPerLevel[i];
1336 if (nnzPerLevel[i] != INVALID) {
1337 oss << std::setw(nnzspacer) << nnzPerLevel[i];
1338 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1339 oss << std::setw(9) << as<double>(nnzPerLevel[i]) / rowsPerLevel[i];
1341 oss << std::setw(nnzspacer) <<
"Operator";
1342 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1343 oss << std::setw(9) <<
" ";
1346 oss << std::setw(9) << as<double>(rowsPerLevel[i - 1]) / rowsPerLevel[i];
1348 oss << std::setw(9) <<
" ";
1349 oss <<
" " << std::setw(npspacer) << numProcsPerLevel[i] << std::endl;
1352 for (
int i = 0; i < GetNumLevels(); ++i) {
1353 RCP<SmootherBase> preSmoo, postSmoo;
1354 if (Levels_[i]->IsAvailable(
"PreSmoother"))
1355 preSmoo = Levels_[i]->
template Get<RCP<SmootherBase>>(
"PreSmoother");
1356 if (Levels_[i]->IsAvailable(
"PostSmoother"))
1357 postSmoo = Levels_[i]->
template Get<RCP<SmootherBase>>(
"PostSmoother");
1359 if (preSmoo != null && preSmoo == postSmoo)
1360 oss <<
"Smoother (level " << i <<
") both : " << preSmoo->description() << std::endl;
1362 oss <<
"Smoother (level " << i <<
") pre : "
1363 << (preSmoo != null ? preSmoo->description() :
"no smoother") << std::endl;
1364 oss <<
"Smoother (level " << i <<
") post : "
1365 << (postSmoo != null ? postSmoo->description() :
"no smoother") << std::endl;
1376 RCP<const Teuchos::MpiComm<int>> mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int>>(comm);
1377 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
1379 int strLength = outstr.size();
1380 MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
1381 if (comm->getRank() != root)
1382 outstr.resize(strLength);
1383 MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
1390template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1392 Teuchos::OSTab tab2(out);
1393 for (
int i = 0; i < GetNumLevels(); ++i)
1394 Levels_[i]->print(out, verbLevel);
1397template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1399 isPreconditioner_ = flag;
1402template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1404 if (GetProcRankVerbose() != 0)
1406#if defined(HAVE_MUELU_BOOST) && defined(HAVE_MUELU_BOOST_FOR_REAL) && defined(BOOST_VERSION) && (BOOST_VERSION >= 104400)
1411 dp.property(
"label", boost::get(boost::vertex_name, graph));
1412 dp.property(
"id", boost::get(boost::vertex_index, graph));
1413 dp.property(
"label", boost::get(boost::edge_name, graph));
1414 dp.property(
"color", boost::get(boost::edge_color, graph));
1417 std::map<const FactoryBase*, BoostVertex> vindices;
1418 typedef std::map<std::pair<BoostVertex, BoostVertex>, std::string> emap;
1421 static int call_id = 0;
1423 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator>>(
"A");
1424 int rank = A->getDomainMap()->getComm()->getRank();
1427 for (
int i = currLevel; i <= currLevel + 1 && i < GetNumLevels(); i++) {
1429 Levels_[i]->UpdateGraph(vindices, edges, dp, graph);
1431 for (emap::const_iterator eit = edges.begin(); eit != edges.end(); eit++) {
1432 std::pair<BoostEdge, bool> boost_edge = boost::add_edge(eit->first.first, eit->first.second, graph);
1435 if (eit->second == std::string(
"Graph"))
1436 boost::put(
"label", dp, boost_edge.first, std::string(
"Graph_"));
1438 boost::put(
"label", dp, boost_edge.first, eit->second);
1440 boost::put(
"color", dp, boost_edge.first, std::string(
"red"));
1442 boost::put(
"color", dp, boost_edge.first, std::string(
"blue"));
1446 std::ofstream out(dumpFile_.c_str() + std::string(
"_") + std::to_string(currLevel) + std::string(
"_") + std::to_string(call_id) + std::string(
"_") + std::to_string(rank) + std::string(
".dot"));
1447 boost::write_graphviz_dp(out, graph, dp, std::string(
"id"));
1451 GetOStream(
Errors) <<
"Dependency graph output requires boost and MueLu_ENABLE_Boost_for_real" << std::endl;
1456template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1458 RCP<Operator> Ao = level.
Get<RCP<Operator>>(
"A");
1459 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Ao);
1461 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is not a matrix, skipping..." << std::endl;
1464 if (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A) != Teuchos::null) {
1465 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is a BlockedCrsMatrix, skipping..." << std::endl;
1469 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LO, GO, NO> xdMV;
1471 RCP<xdMV> coords = level.
Get<RCP<xdMV>>(
"Coordinates");
1473 if (A->getRowMap()->isSameAs(*(coords->getMap()))) {
1474 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: matrix and coordinates maps are same, skipping..." << std::endl;
1478 if (A->IsView(
"stridedMaps") && rcp_dynamic_cast<const StridedMap>(A->getRowMap(
"stridedMaps")) != Teuchos::null) {
1479 RCP<const StridedMap> stridedRowMap = rcp_dynamic_cast<const StridedMap>(A->getRowMap(
"stridedMaps"));
1482 TEUCHOS_TEST_FOR_EXCEPTION(stridedRowMap->getStridedBlockId() != -1 || stridedRowMap->getOffset() != 0,
1483 Exceptions::RuntimeError,
"Hierarchy::ReplaceCoordinateMap: nontrivial maps (block id = " << stridedRowMap->getStridedBlockId() <<
", offset = " << stridedRowMap->getOffset() <<
")");
1486 GetOStream(
Runtime1) <<
"Replacing coordinate map" << std::endl;
1487 TEUCHOS_TEST_FOR_EXCEPTION(A->GetFixedBlockSize() % A->GetStorageBlockSize() != 0,
Exceptions::RuntimeError,
"Hierarchy::ReplaceCoordinateMap: Storage block size does not evenly divide fixed block size");
1489 size_t blkSize = A->GetFixedBlockSize() / A->GetStorageBlockSize();
1491 RCP<const Map> nodeMap = A->getRowMap();
1494 RCP<const Map> dofMap = A->getRowMap();
1495 GO indexBase = dofMap->getIndexBase();
1496 size_t numLocalDOFs = dofMap->getLocalNumElements();
1498 "Hierarchy::ReplaceCoordinateMap: block size (" << blkSize <<
") is incompatible with the number of local dofs in a row map (" << numLocalDOFs);
1499 ArrayView<const GO> GIDs = dofMap->getLocalElementList();
1501 Array<GO> nodeGIDs(numLocalDOFs / blkSize);
1502 for (
size_t i = 0; i < numLocalDOFs; i += blkSize)
1503 nodeGIDs[i / blkSize] = (GIDs[i] - indexBase) / blkSize + indexBase;
1505 Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1506 nodeMap = MapFactory::Build(dofMap->lib(), INVALID, nodeGIDs(), indexBase, dofMap->getComm());
1512 if (coords->getLocalLength() != A->getRowMap()->getLocalNumElements()) {
1513 GetOStream(
Warnings) <<
"Coordinate vector does not match row map of matrix A!" << std::endl;
1518 Array<ArrayView<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>> coordDataView;
1519 std::vector<ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>> coordData;
1520 for (
size_t i = 0; i < coords->getNumVectors(); i++) {
1521 coordData.push_back(coords->getData(i));
1522 coordDataView.push_back(coordData[i]());
1525 RCP<xdMV> newCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LO, GO, NO>::Build(nodeMap, coordDataView(), coords->getNumVectors());
1526 level.
Set(
"Coordinates", newCoords);
1529template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1531 int N = Levels_.size();
1532 if (((sizeOfAllocatedLevelMultiVectors_ == numvecs && residual_.size() == N) || numvecs <= 0) && !forceMapCheck)
return;
1535 if (residual_.size() != N) {
1536 DeleteLevelMultiVectors();
1538 residual_.resize(N);
1539 coarseRhs_.resize(N);
1541 coarseImport_.resize(N);
1542 coarseExport_.resize(N);
1543 correction_.resize(N);
1546 for (
int i = 0; i < N; i++) {
1547 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator>>(
"A");
1550 RCP<const BlockedCrsMatrix> A_as_blocked = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(A);
1551 RCP<const Map> Arm = A->getRangeMap();
1552 RCP<const Map> Adm = A->getDomainMap();
1553 if (!A_as_blocked.is_null()) {
1554 Adm = A_as_blocked->getFullDomainMap();
1557 if (residual_[i].is_null() || !residual_[i]->getMap()->isSameAs(*Arm))
1559 residual_[i] = MultiVectorFactory::Build(Arm, numvecs,
true);
1560 if (correction_[i].is_null() || !correction_[i]->getMap()->isSameAs(*Adm))
1561 correction_[i] = MultiVectorFactory::Build(Adm, numvecs,
false);
1566 if (implicitTranspose_) {
1567 RCP<Operator> P = Levels_[i + 1]->template Get<RCP<Operator>>(
"P");
1569 RCP<const Map> map = P->getDomainMap();
1570 if (coarseRhs_[i].is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1571 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1574 RCP<Operator> R = Levels_[i + 1]->template Get<RCP<Operator>>(
"R");
1576 RCP<const Map> map = R->getRangeMap();
1577 if (coarseRhs_[i].is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1578 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1582 RCP<const Import> importer;
1583 if (Levels_[i + 1]->IsAvailable(
"Importer"))
1584 importer = Levels_[i + 1]->
template Get<RCP<const Import>>(
"Importer");
1585 if (doPRrebalance_ || importer.is_null()) {
1586 RCP<const Map> map = coarseRhs_[i]->getMap();
1587 if (coarseX_[i].is_null() || !coarseX_[i]->getMap()->isSameAs(*map))
1588 coarseX_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1591 map = importer->getTargetMap();
1592 if (coarseImport_[i].is_null() || !coarseImport_[i]->getMap()->isSameAs(*map)) {
1593 coarseImport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1594 coarseX_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1596 map = importer->getSourceMap();
1597 if (coarseExport_[i].is_null() || !coarseExport_[i]->getMap()->isSameAs(*map))
1598 coarseExport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1602 sizeOfAllocatedLevelMultiVectors_ = numvecs;
1605template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1607 if (sizeOfAllocatedLevelMultiVectors_ == 0)
return;
1608 residual_.resize(0);
1609 coarseRhs_.resize(0);
1611 coarseImport_.resize(0);
1612 coarseExport_.resize(0);
1613 correction_.resize(0);
1614 sizeOfAllocatedLevelMultiVectors_ = 0;
1617template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1619 const LO startLevel,
const ConvData& conv)
const {
1620 return (startLevel == 0 && !isPreconditioner_ && (IsPrint(
Statistics1) || conv.
tol_ > 0));
1623template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1625 const Teuchos::Array<MagnitudeType>& residualNorm,
const MagnitudeType convergenceTolerance)
const {
1628 if (convergenceTolerance > Teuchos::ScalarTraits<MagnitudeType>::zero()) {
1630 for (LO k = 0; k < residualNorm.size(); k++)
1631 if (residualNorm[k] >= convergenceTolerance)
1640 return convergenceStatus;
1643template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1645 const LO iteration,
const Teuchos::Array<MagnitudeType>& residualNorm)
const {
1647 << std::setiosflags(std::ios::left)
1648 << std::setprecision(3) << std::setw(4) << iteration
1650 << std::setprecision(10) << residualNorm
1654template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1656 const Operator& A,
const MultiVector& X,
const MultiVector& B,
const LO iteration,
1658 Teuchos::Array<MagnitudeType> residualNorm;
1662 rate_ = currentResidualNorm / previousResidualNorm;
1663 previousResidualNorm = currentResidualNorm;
1666 PrintResidualHistory(iteration, residualNorm);
1668 return IsConverged(residualNorm, conv.
tol_);
static bool debug()
Whether MueLu is in debug mode.
virtual std::string description() const
Return a simple one-line description of this object.
Exception throws to report incompatible objects (like maps).
Exception throws to report errors in the internal logical of the program.
Base class for factories (e.g., R, P, and A_coarse).
Class that provides default factories within Needs class.
virtual void Clean() const
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
void AddLevel(const RCP< Level > &level)
Add a level at the end of the hierarchy.
double GetSmootherComplexity() const
void Write(const LO &start=-1, const LO &end=-1, const std::string &suffix="")
Print matrices in the multigrid hierarchy to file.
RCP< Level > & GetLevel(const int levelID=0)
Retrieve a certain level from hierarchy.
virtual ~Hierarchy()
Destructor.
void CheckLevel(Level &level, int levelID)
Helper function.
std::string description() const
Return a simple one-line description of this object.
void CheckForEmptySmoothersAndCoarseSolve()
void IsPreconditioner(const bool flag)
Array< RCP< Level > > Levels_
Container for Level objects.
bool Setup(int coarseLevelID, const RCP< const FactoryManagerBase > fineLevelManager, const RCP< const FactoryManagerBase > coarseLevelManager, const RCP< const FactoryManagerBase > nextLevelManager=Teuchos::null)
Multi-level setup phase: build a new level of the hierarchy.
STS::magnitudeType MagnitudeType
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
ConvergenceStatus IsConverged(const Teuchos::Array< MagnitudeType > &residualNorm, const MagnitudeType convergenceTolerance) const
Decide if the multigrid iteration is converged.
void DeleteLevelMultiVectors()
ConvergenceStatus Iterate(const MultiVector &B, MultiVector &X, ConvData conv=ConvData(), bool InitialGuessIsZero=false, LO startLevel=0)
Apply the multigrid preconditioner.
void DumpCurrentGraph(int level) const
void SetMatvecParams(RCP< ParameterList > matvecParams)
Xpetra::UnderlyingLib lib_
Epetra/Tpetra mode.
void Clear(int startLevel=0)
Clear impermanent data from previous setup.
bool IsCalculationOfResidualRequired(const LO startLevel, const ConvData &conv) const
Decide if the residual needs to be computed.
ConvergenceStatus ComputeResidualAndPrintHistory(const Operator &A, const MultiVector &X, const MultiVector &B, const LO iteration, const LO startLevel, const ConvData &conv, MagnitudeType &previousResidualNorm)
Compute the residual norm and print it depending on the verbosity level.
double GetOperatorComplexity() const
void PrintResidualHistory(const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const
Print residualNorm for this iteration to the screen.
void AllocateLevelMultiVectors(int numvecs, bool forceMapCheck=false)
void print(std::ostream &out=std::cout, const VerbLevel verbLevel=(MueLu::Parameters|MueLu::Statistics0)) const
Hierarchy::print is local hierarchy function, thus the statistics can be different from global ones.
void Delete(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Delete(ename, factory) for each level of the Hierarchy.
int GetGlobalNumLevels() const
void Keep(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Keep(ename, factory) for each level of the Hierarchy.
void SetLabel(const std::string &hierarchyLabel)
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
Call Level::AddKeepFlag for each level of the Hierarchy.
void AddNewLevel()
Add a new level at the end of the hierarchy.
void RemoveKeepFlag(const std::string &ename, const FactoryBase *factory, KeepType keep=MueLu::All)
Call Level::RemoveKeepFlag for each level of the Hierarchy.
void ReplaceCoordinateMap(Level &level)
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void SetComm(RCP< const Teuchos::Comm< int > > const &comm)
void setlib(Xpetra::UnderlyingLib lib2)
RCP< Level > & GetPreviousLevel()
Previous level.
void Release(const FactoryBase &factory)
Decrement the storage counter for all the inputs of a factory.
RCP< const Teuchos::Comm< int > > GetComm() const
int GetLevelID() const
Return level number.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
void Request(const FactoryBase &factory)
Increment the storage counter for all the inputs of a factory.
Xpetra::UnderlyingLib lib()
Timer to be used in non-factories.
static const NoFactory * get()
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
An exception safe way to call the method 'Level::SetFactoryManager()'.
Integrates Teuchos::TimeMonitor with MueLu verbosity system.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
Namespace for MueLu classes and methods.
@ Warnings0
Important warning messages (one line)
@ Statistics2
Print even more statistics.
@ Warnings
Print all warning messages.
@ Statistics1
Print more statistics.
@ Timings0
High level timing information (use Teuchos::TimeMonitor::summarize() to print)
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose)
@ Warnings1
Additional warnings.
@ Statistics0
Print statistics that do not involve significant additional computation.
@ Parameters1
Print class parameters (more parameters, more verbose)
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
VerbLevel toMueLuVerbLevel(const Teuchos::EVerbosityLevel verbLevel)
Translate Teuchos verbosity level to MueLu verbosity level.
Data struct for defining stopping criteria of multigrid iteration.