10 #ifndef MUELU_HIERARCHY_DEF_HPP
11 #define MUELU_HIERARCHY_DEF_HPP
19 #include <Xpetra_MultiVectorFactory.hpp>
25 #include "MueLu_FactoryManager.hpp"
26 #include "MueLu_HierarchyUtils.hpp"
27 #include "MueLu_TopRAPFactory.hpp"
28 #include "MueLu_TopSmootherFactory.hpp"
31 #include "MueLu_PerfUtils.hpp"
32 #include "MueLu_PFactory.hpp"
33 #include "MueLu_SmootherFactory.hpp"
40 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
42 : maxCoarseSize_(GetDefaultMaxCoarseSize())
43 , implicitTranspose_(GetDefaultImplicitTranspose())
44 , fuseProlongationAndUpdate_(GetDefaultFuseProlongationAndUpdate())
45 , doPRrebalance_(GetDefaultPRrebalance())
46 , doPRViaCopyrebalance_(false)
47 , isPreconditioner_(true)
48 , Cycle_(GetDefaultCycle())
49 , WCycleStartLevel_(0)
50 , scalingFactor_(Teuchos::ScalarTraits<double>::one())
52 , isDumpingEnabled_(false)
55 , sizeOfAllocatedLevelMultiVectors_(0) {
59 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
63 Levels_[0]->setObjectLabel(label);
66 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
68 : maxCoarseSize_(GetDefaultMaxCoarseSize())
69 , implicitTranspose_(GetDefaultImplicitTranspose())
70 , fuseProlongationAndUpdate_(GetDefaultFuseProlongationAndUpdate())
71 , doPRrebalance_(GetDefaultPRrebalance())
72 , doPRViaCopyrebalance_(false)
73 , isPreconditioner_(true)
74 , Cycle_(GetDefaultCycle())
75 , WCycleStartLevel_(0)
76 , scalingFactor_(Teuchos::ScalarTraits<double>::one())
77 , isDumpingEnabled_(false)
80 , sizeOfAllocatedLevelMultiVectors_(0) {
81 lib_ = A->getDomainMap()->lib();
89 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
93 Levels_[0]->setObjectLabel(label);
96 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
99 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
101 int levelID = LastLevelID() + 1;
104 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"
105 <<
" because last level ID of the hierarchy was " << LastLevelID() <<
"." << std::endl;
107 Levels_.push_back(level);
111 level->
SetPreviousLevel((levelID == 0) ? Teuchos::null : Levels_[LastLevelID() - 1]);
115 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
117 RCP<Level> newLevel = Levels_[LastLevelID()]->Build();
119 this->AddLevel(newLevel);
122 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 "MueLu::Hierarchy::GetLevel(): invalid input parameter value: LevelID = " << levelID);
126 return Levels_[levelID];
129 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
131 return Levels_.size();
134 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
136 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator> >(
"A");
139 int numLevels = GetNumLevels();
143 return numGlobalLevels;
146 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
148 double totalNnz = 0, lev0Nnz = 1;
149 for (
int i = 0; i < GetNumLevels(); ++i) {
151 "Operator complexity cannot be calculated because A is unavailable on level " << i);
152 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator> >(
"A");
158 GetOStream(
Warnings0) <<
"Some level operators are not matrices, operator complexity calculation aborted" << std::endl;
162 totalNnz += as<double>(Am->getGlobalNumEntries());
166 return totalNnz / lev0Nnz;
169 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
171 double node_sc = 0, global_sc = 0;
175 if (GetNumLevels() <= 0)
return -1.0;
176 if (!Levels_[0]->IsAvailable(
"A"))
return -1.0;
178 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator> >(
"A");
182 a0_nnz = as<double>(Am->getGlobalNumEntries());
185 for (
int i = 0; i < GetNumLevels(); ++i) {
187 if (!Levels_[i]->IsAvailable(
"PreSmoother"))
continue;
188 RCP<SmootherBase> S = Levels_[i]->template Get<RCP<SmootherBase> >(
"PreSmoother");
191 if (level_sc == INVALID) {
196 node_sc += as<double>(level_sc);
207 return global_sc / a0_nnz;
211 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
214 "MueLu::Hierarchy::CheckLevel(): wrong underlying linear algebra library.");
216 "MueLu::Hierarchy::CheckLevel(): wrong level ID");
218 "MueLu::Hierarchy::Setup(): wrong level parent");
221 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
223 for (
int i = 0; i < GetNumLevels(); ++i) {
231 xpImporter->setDistributorParameters(matvecParams);
234 xpExporter->setDistributorParameters(matvecParams);
241 xpImporter->setDistributorParameters(matvecParams);
244 xpExporter->setDistributorParameters(matvecParams);
250 xpImporter->setDistributorParameters(matvecParams);
253 xpExporter->setDistributorParameters(matvecParams);
258 xpImporter->setDistributorParameters(matvecParams);
265 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
274 "MueLu::Hierarchy:Setup(): level " << coarseLevelID <<
" (specified by coarseLevelID argument) "
275 "must be built before calling this function.");
277 Level& level = *Levels_[coarseLevelID];
280 TimeMonitor m1(*
this, label + this->ShortClassName() +
": " +
"Setup (total)");
285 "MueLu::Hierarchy::Setup(): argument coarseLevelManager cannot be null");
290 if (levelManagers_.size() < coarseLevelID + 1)
291 levelManagers_.resize(coarseLevelID + 1);
292 levelManagers_[coarseLevelID] = coarseLevelManager;
294 bool isFinestLevel = (fineLevelManager.
is_null());
295 bool isLastLevel = (nextLevelManager.
is_null());
308 oldRank = SetProcRankVerbose(comm->getRank());
312 lib_ = domainMap->lib();
319 Level& prevLevel = *Levels_[coarseLevelID - 1];
320 oldRank = SetProcRankVerbose(prevLevel.
GetComm()->getRank());
323 CheckLevel(level, coarseLevelID);
330 if (isFinestLevel && Levels_[coarseLevelID]->IsAvailable(
"Coordinates"))
331 ReplaceCoordinateMap(*Levels_[coarseLevelID]);
336 if (isDumpingEnabled_ && (dumpLevel_ == 0 || dumpLevel_ == -1) && coarseLevelID == 1)
342 int nextLevelID = coarseLevelID + 1;
345 if (isLastLevel ==
false) {
347 if (nextLevelID > LastLevelID())
349 CheckLevel(*Levels_[nextLevelID], nextLevelID);
353 Levels_[nextLevelID]->Request(TopRAPFactory(coarseLevelManager, nextLevelManager));
383 coarseFact =
rcp(
new TopSmootherFactory(coarseLevelManager,
"CoarseSolver"));
392 TopRAPFactory coarseRAPFactory(fineLevelManager, coarseLevelManager);
396 }
else if (!isFinestLevel) {
401 bool setLastLevelviaMaxCoarseSize =
false;
408 level.
SetComm(Ac->getDomainMap()->getComm());
411 bool isOrigLastLevel = isLastLevel;
423 if (!Acm.
is_null() && Acm->getGlobalNumRows() <= maxCoarseSize_) {
425 GetOStream(
Runtime0) <<
"Max coarse size (<= " << maxCoarseSize_ <<
") achieved" << std::endl;
427 if (Acm->getGlobalNumRows() != 0) setLastLevelviaMaxCoarseSize =
true;
431 if (!Ac.
is_null() && !isFinestLevel) {
432 RCP<Operator> A = Levels_[coarseLevelID - 1]->template Get<RCP<Operator> >(
"A");
435 const double maxCoarse2FineRatio = 0.8;
436 if (!Acm.
is_null() && !Am.is_null() && Acm->getGlobalNumRows() > maxCoarse2FineRatio * Am->getGlobalNumRows()) {
444 GetOStream(
Warnings0) <<
"Aggregation stagnated. Please check your matrix and/or adjust your configuration file."
445 <<
"Possible fixes:\n"
446 <<
" - reduce the maximum number of levels\n"
447 <<
" - enable repartitioning\n"
448 <<
" - increase the minimum coarse size." << std::endl;
453 if (!isOrigLastLevel) {
458 coarseFact =
rcp(
new TopSmootherFactory(coarseLevelManager,
"CoarseSolver"));
465 coarseFact->
Build(level);
476 smootherFact->
Build(level);
481 if (isLastLevel ==
true) {
482 int actualNumLevels = nextLevelID;
483 if (isOrigLastLevel ==
false) {
486 Levels_[nextLevelID]->Release(TopRAPFactory(coarseLevelManager, nextLevelManager));
492 if (!setLastLevelviaMaxCoarseSize) {
493 if (Levels_[nextLevelID - 1]->IsAvailable(
"P")) {
494 if (Levels_[nextLevelID - 1]->
template Get<
RCP<Matrix> >(
"P") == Teuchos::null) actualNumLevels = nextLevelID - 1;
496 actualNumLevels = nextLevelID - 1;
499 if (actualNumLevels == nextLevelID - 1) {
501 Levels_[nextLevelID - 2]->Release(*smootherFact);
503 if (Levels_[nextLevelID - 2]->IsAvailable(
"PreSmoother")) Levels_[nextLevelID - 2]->RemoveKeepFlag(
"PreSmoother",
NoFactory::get());
504 if (Levels_[nextLevelID - 2]->IsAvailable(
"PostSmoother")) Levels_[nextLevelID - 2]->RemoveKeepFlag(
"PostSmoother",
NoFactory::get());
506 coarseFact =
rcp(
new TopSmootherFactory(coarseLevelManager,
"CoarseSolver"));
507 Levels_[nextLevelID - 2]->Request(*coarseFact);
509 coarseFact->
Build(*(Levels_[nextLevelID - 2]));
510 Levels_[nextLevelID - 2]->Release(*coarseFact);
512 Levels_.resize(actualNumLevels);
516 if (isDumpingEnabled_ && ((dumpLevel_ > 0 && coarseLevelID == dumpLevel_) || dumpLevel_ == -1))
517 DumpCurrentGraph(coarseLevelID);
519 if (!isFinestLevel) {
523 level.
Release(coarseRAPFactory);
527 SetProcRankVerbose(oldRank);
532 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
534 int numLevels = Levels_.size();
536 "Hierarchy::SetupRe: " << Levels_.size() <<
" levels, but " << levelManagers_.size() <<
" level factory managers");
538 const int startLevel = 0;
541 #ifdef HAVE_MUELU_DEBUG
543 for (
int i = 0; i < numLevels; i++)
544 levelManagers_[i]->ResetDebugData();
549 for (levelID = startLevel; levelID < numLevels;) {
550 bool r = Setup(levelID,
551 (levelID != 0 ? levelManagers_[levelID - 1] : Teuchos::null),
552 levelManagers_[levelID],
553 (levelID + 1 != numLevels ? levelManagers_[levelID + 1] : Teuchos::null));
559 Levels_.resize(levelID);
560 levelManagers_.resize(levelID);
562 int sizeOfVecs = sizeOfAllocatedLevelMultiVectors_;
564 AllocateLevelMultiVectors(sizeOfVecs,
true);
571 CheckForEmptySmoothersAndCoarseSolve();
574 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
583 "MueLu::Hierarchy::Setup(): fine level (" << startLevel <<
") does not exist");
586 "Constructing non-positive (" << numDesiredLevels <<
") number of levels does not make sense.");
589 TEUCHOS_TEST_FOR_EXCEPTION(!Levels_[startLevel]->IsAvailable(
"A"), Exceptions::RuntimeError,
590 "MueLu::Hierarchy::Setup(): fine level (" << startLevel <<
") has no matrix A! "
591 "Set fine level matrix A using Level.Set()");
593 RCP<Operator> A = Levels_[startLevel]->template Get<RCP<Operator> >(
"A");
594 lib_ = A->getDomainMap()->lib();
601 params->
set(
"printLoadBalancingInfo",
true);
602 params->
set(
"printCommInfo",
true);
606 GetOStream(
Warnings1) <<
"Fine level operator is not a matrix, statistics are not available" << std::endl;
612 const int lastLevel = startLevel + numDesiredLevels - 1;
613 GetOStream(
Runtime0) <<
"Setup loop: startLevel = " << startLevel <<
", lastLevel = " << lastLevel
614 <<
" (stop if numLevels = " << numDesiredLevels <<
" or Ac.size() < " << maxCoarseSize_ <<
")" << std::endl;
618 if (numDesiredLevels == 1) {
620 Setup(startLevel, Teuchos::null, rcpmanager, Teuchos::null);
623 bool bIsLastLevel = Setup(startLevel, Teuchos::null, rcpmanager, rcpmanager);
624 if (bIsLastLevel ==
false) {
625 for (iLevel = startLevel + 1; iLevel < lastLevel; iLevel++) {
626 bIsLastLevel = Setup(iLevel, rcpmanager, rcpmanager, rcpmanager);
627 if (bIsLastLevel ==
true)
630 if (bIsLastLevel ==
false)
631 Setup(lastLevel, rcpmanager, rcpmanager, Teuchos::null);
636 TEUCHOS_TEST_FOR_EXCEPTION(iLevel != Levels_.size() - 1, Exceptions::RuntimeError,
637 "MueLu::Hierarchy::Setup(): number of level");
645 CheckForEmptySmoothersAndCoarseSolve();
648 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
650 for (
LO levelNo = 0; levelNo < as<LO>(Levels_.size()); ++levelNo) {
651 auto level = Levels_[levelNo];
652 if ((!level->IsAvailable(
"PreSmoother")) && (!level->IsAvailable(
"PostSmoother")))
653 GetOStream(
Warnings1) <<
"No " << (levelNo == as<LO>(Levels_.size()) - 1 ?
"coarse grid solver" :
"smoother") <<
" on level " << level->GetLevelID() << std::endl;
657 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
659 if (startLevel < GetNumLevels())
660 GetOStream(
Runtime0) <<
"Clearing old data (if any)" << std::endl;
662 for (
int iLevel = startLevel; iLevel < GetNumLevels(); iLevel++)
663 Levels_[iLevel]->Clear();
666 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
668 GetOStream(
Runtime0) <<
"Clearing old data (expert)" << std::endl;
669 for (
int iLevel = 0; iLevel < GetNumLevels(); iLevel++)
670 Levels_[iLevel]->ExpertClear();
673 #if defined(HAVE_MUELU_EXPERIMENTAL) && defined(HAVE_MUELU_ADDITIVE_VARIANT)
674 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
676 bool InitialGuessIsZero,
LO startLevel) {
677 LO nIts = conv.maxIts_;
678 MagnitudeType tol = conv.tol_;
680 std::string prefix = this->ShortClassName() +
": ";
681 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
682 std::string levelSuffix1 =
" (level=" +
toString(startLevel + 1) +
")";
684 using namespace Teuchos;
708 SC one = STS::one(), zero = STS::zero();
710 bool zeroGuess = InitialGuessIsZero;
721 bool emptyCoarseSolve =
true;
722 RCP<MultiVector> coarseX_prolonged = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
726 if (Levels_.size() > 1) {
728 if (Coarse->IsAvailable(
"Importer"))
737 coarseRhs = MultiVectorFactory::Build(R->getRangeMap(), B.getNumVectors(),
true);
746 if (doPRrebalance_ || importer.is_null()) {
747 coarseX = MultiVectorFactory::Build(coarseRhs->getMap(), X.getNumVectors(),
true);
754 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getTargetMap(), coarseRhs->getNumVectors());
756 coarseRhs.
swap(coarseTmp);
758 coarseX = MultiVectorFactory::Build(importer->getTargetMap(), X.getNumVectors(),
true);
761 if (Coarse->IsAvailable(
"PreSmoother"))
763 if (Coarse->IsAvailable(
"PostSmoother"))
769 MagnitudeType prevNorm = STS::magnitude(STS::one()), curNorm = STS::magnitude(STS::one());
772 for (
LO i = 1; i <= nIts; i++) {
773 #ifdef HAVE_MUELU_DEBUG
774 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
775 std::ostringstream ss;
776 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
777 throw Exceptions::Incompatible(ss.str());
780 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
781 std::ostringstream ss;
782 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
783 throw Exceptions::Incompatible(ss.str());
788 bool emptyFineSolve =
true;
791 fineX = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
800 if (Fine->IsAvailable(
"PreSmoother")) {
803 preSmoo->Apply(*fineX, B, zeroGuess);
805 emptyFineSolve =
false;
807 if (Fine->IsAvailable(
"PostSmoother")) {
810 postSmoo->Apply(*fineX, B, zeroGuess);
813 emptyFineSolve =
false;
815 if (emptyFineSolve ==
true) {
817 fineX->update(one, B, zero);
820 if (Levels_.size() > 1) {
822 if (Coarse->IsAvailable(
"PreSmoother")) {
824 preSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
826 emptyCoarseSolve =
false;
828 if (Coarse->IsAvailable(
"PostSmoother")) {
830 postSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
832 emptyCoarseSolve =
false;
834 if (emptyCoarseSolve ==
true) {
836 coarseX->update(one, *coarseRhs, zero);
843 if (!doPRrebalance_ && !importer.is_null()) {
848 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getSourceMap(), coarseX->getNumVectors());
850 coarseX.
swap(coarseTmp);
859 X.update(1.0, *fineX, 1.0, *coarseX_prolonged, 0.0);
870 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
872 bool InitialGuessIsZero,
LO startLevel) {
891 std::string prefix = label + this->ShortClassName() +
": ";
892 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
893 std::string levelSuffix1 =
" (level=" +
toString(startLevel + 1) +
")";
901 else if (!useStackedTimer)
904 std::string iterateLevelTimeLabel = prefix +
"Solve" + levelSuffix;
907 bool zeroGuess = InitialGuessIsZero;
910 using namespace Teuchos;
924 const BlockedMultiVector* Bblocked =
dynamic_cast<const BlockedMultiVector*
>(&B);
925 if (residual_.size() > startLevel &&
926 ((Bblocked && !Bblocked->isSameSize(*residual_[startLevel])) ||
927 (!Bblocked && !residual_[startLevel]->isSameSize(B))))
928 DeleteLevelMultiVectors();
929 AllocateLevelMultiVectors(X.getNumVectors());
935 if (IsCalculationOfResidualRequired(startLevel, conv)) {
938 return convergenceStatus;
941 SC one = STS::one(), zero = STS::zero();
942 for (
LO iteration = 1; iteration <= nIts; iteration++) {
943 #ifdef HAVE_MUELU_DEBUG
945 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
946 std::ostringstream ss;
947 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
951 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
952 std::ostringstream ss;
953 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
959 if (startLevel == as<LO>(Levels_.size()) - 1) {
963 bool emptySolve =
true;
966 if (Fine->IsAvailable(
"PreSmoother")) {
969 preSmoo->Apply(X, B, zeroGuess);
974 if (Fine->IsAvailable(
"PostSmoother")) {
977 postSmoo->Apply(X, B, zeroGuess);
982 if (emptySolve ==
true) {
984 X.update(one, B, zero);
993 if (!useStackedTimer)
997 if (Fine->IsAvailable(
"PreSmoother")) {
999 preSmoo->Apply(X, B, zeroGuess);
1007 if (!useStackedTimer)
1017 residual = residual_[startLevel];
1021 if (Coarse->IsAvailable(
"Pbar"))
1029 if (!useStackedTimer)
1032 coarseRhs = coarseRhs_[startLevel];
1034 if (implicitTranspose_) {
1044 if (Coarse->IsAvailable(
"Importer"))
1047 coarseX = coarseX_[startLevel];
1048 if (!doPRrebalance_ && !importer.
is_null()) {
1050 if (!useStackedTimer)
1057 coarseRhs.
swap(coarseTmp);
1066 coarseRhs->replaceMap(Ac->getRangeMap());
1067 coarseX->replaceMap(Ac->getDomainMap());
1070 iterateLevelTime = Teuchos::null;
1072 Iterate(*coarseRhs, *coarseX, 1,
true, startLevel + 1);
1074 if (Cycle_ ==
WCYCLE && WCycleStartLevel_ >= startLevel)
1075 Iterate(*coarseRhs, *coarseX, 1,
false, startLevel + 1);
1078 iterateLevelTime =
rcp(
new TimeMonitor(*
this, iterateLevelTimeLabel));
1080 coarseX->replaceMap(origXMap);
1081 coarseRhs->replaceMap(origRhsMap);
1084 if (!doPRrebalance_ && !importer.
is_null()) {
1086 if (!useStackedTimer)
1093 coarseX.
swap(coarseTmp);
1099 if (!useStackedTimer)
1106 if (fuseProlongationAndUpdate_) {
1111 X.update(scalingFactor_, *correction, one);
1118 if (!useStackedTimer)
1122 if (Fine->IsAvailable(
"PostSmoother")) {
1124 postSmoo->Apply(X, B,
false);
1130 if (IsCalculationOfResidualRequired(startLevel, conv)) {
1131 ConvergenceStatus convergenceStatus = ComputeResidualAndPrintHistory(*A, X, B, iteration, startLevel, conv, prevNorm);
1133 return convergenceStatus;
1140 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1142 LO startLevel = (start != -1 ? start : 0);
1143 LO endLevel = (end != -1 ? end : Levels_.size() - 1);
1146 "MueLu::Hierarchy::Write : startLevel must be <= endLevel");
1149 "MueLu::Hierarchy::Write bad start or end level");
1151 for (
LO i = startLevel; i < endLevel + 1; i++) {
1152 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Levels_[i]->template Get<RCP<Operator> >(
"A")), P, R;
1154 P = rcp_dynamic_cast<Matrix>(Levels_[i]->template Get<RCP<Operator> >(
"P"));
1155 if (!implicitTranspose_)
1156 R = rcp_dynamic_cast<Matrix>(Levels_[i]->template Get<RCP<Operator> >(
"R"));
1169 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1171 for (
Array<
RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1172 (*it)->Keep(ename, factory);
1175 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1177 for (
Array<
RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1178 (*it)->Delete(ename, factory);
1181 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1183 for (
Array<
RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1184 (*it)->AddKeepFlag(ename, factory, keep);
1187 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1189 for (
Array<
RCP<Level> >::iterator it = Levels_.begin(); it != Levels_.end(); ++it)
1190 (*it)->RemoveKeepFlag(ename, factory, keep);
1193 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1195 if (description_.empty()) {
1196 std::ostringstream out;
1198 out <<
"{#levels = " << GetGlobalNumLevels() <<
", complexity = " << GetOperatorComplexity() <<
"}";
1199 description_ = out.str();
1201 return description_;
1204 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1209 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1211 RCP<Operator> A0 = Levels_[0]->template Get<RCP<Operator> >(
"A");
1214 int numLevels = GetNumLevels();
1215 RCP<Operator> Ac = Levels_[numLevels - 1]->template Get<RCP<Operator> >(
"A");
1222 int root = comm->getRank();
1225 int smartData = numLevels * comm->getSize() + comm->getRank(), maxSmartData;
1227 root = maxSmartData % comm->getSize();
1231 double smoother_comp = -1.0;
1233 smoother_comp = GetSmootherComplexity();
1237 std::vector<Xpetra::global_size_t> nnzPerLevel;
1238 std::vector<Xpetra::global_size_t> rowsPerLevel;
1239 std::vector<int> numProcsPerLevel;
1240 bool someOpsNotMatrices =
false;
1242 for (
int i = 0; i < numLevels; i++) {
1244 "Operator A is not available on level " << i);
1246 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator> >(
"A");
1248 "Operator A on level " << i <<
" is null.");
1252 someOpsNotMatrices =
true;
1253 nnzPerLevel.push_back(INVALID);
1254 rowsPerLevel.push_back(A->getDomainMap()->getGlobalNumElements());
1255 numProcsPerLevel.push_back(A->getDomainMap()->getComm()->getSize());
1257 LO storageblocksize = Am->GetStorageBlockSize();
1259 nnzPerLevel.push_back(nnz);
1260 rowsPerLevel.push_back(Am->getGlobalNumRows() * storageblocksize);
1261 numProcsPerLevel.push_back(Am->getRowMap()->getComm()->getSize());
1264 if (someOpsNotMatrices)
1265 GetOStream(
Warnings0) <<
"Some level operators are not matrices, statistics calculation are incomplete" << std::endl;
1268 std::string label = Levels_[0]->getObjectLabel();
1269 std::ostringstream oss;
1270 oss << std::setfill(
' ');
1271 oss <<
"\n--------------------------------------------------------------------------------\n";
1272 oss <<
"--- Multigrid Summary " << std::setw(28) << std::left << label <<
"---\n";
1273 oss <<
"--------------------------------------------------------------------------------" << std::endl;
1276 oss <<
"Number of levels = " << numLevels << std::endl;
1277 oss <<
"Operator complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1278 if (!someOpsNotMatrices)
1279 oss << GetOperatorComplexity() << std::endl;
1281 oss <<
"not available (Some operators in hierarchy are not matrices.)" << std::endl;
1283 if (smoother_comp != -1.0) {
1284 oss <<
"Smoother complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed)
1285 << smoother_comp << std::endl;
1290 oss <<
"Cycle type = V" << std::endl;
1293 oss <<
"Cycle type = W" << std::endl;
1294 if (WCycleStartLevel_ > 0)
1295 oss <<
"Cycle start level = " << WCycleStartLevel_ << std::endl;
1308 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1309 tt = nnzPerLevel[i];
1319 tt = numProcsPerLevel[0];
1325 oss <<
"level " << std::setw(rowspacer) <<
" rows " << std::setw(nnzspacer) <<
" nnz "
1326 <<
" nnz/row" << std::setw(npspacer) <<
" c ratio"
1327 <<
" procs" << std::endl;
1328 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1329 oss <<
" " << i <<
" ";
1330 oss << std::setw(rowspacer) << rowsPerLevel[i];
1331 if (nnzPerLevel[i] != INVALID) {
1332 oss << std::setw(nnzspacer) << nnzPerLevel[i];
1333 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1334 oss << std::setw(9) << as<double>(nnzPerLevel[i]) / rowsPerLevel[i];
1336 oss << std::setw(nnzspacer) <<
"Operator";
1337 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1338 oss << std::setw(9) <<
" ";
1341 oss << std::setw(9) << as<double>(rowsPerLevel[i - 1]) / rowsPerLevel[i];
1343 oss << std::setw(9) <<
" ";
1344 oss <<
" " << std::setw(npspacer) << numProcsPerLevel[i] << std::endl;
1347 for (
int i = 0; i < GetNumLevels(); ++i) {
1349 if (Levels_[i]->IsAvailable(
"PreSmoother"))
1351 if (Levels_[i]->IsAvailable(
"PostSmoother"))
1354 if (preSmoo != null && preSmoo == postSmoo)
1355 oss <<
"Smoother (level " << i <<
") both : " << preSmoo->
description() << std::endl;
1357 oss <<
"Smoother (level " << i <<
") pre : "
1358 << (preSmoo != null ? preSmoo->
description() :
"no smoother") << std::endl;
1359 oss <<
"Smoother (level " << i <<
") post : "
1360 << (postSmoo != null ? postSmoo->
description() :
"no smoother") << std::endl;
1372 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
1374 int strLength = outstr.size();
1375 MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
1376 if (comm->getRank() != root)
1377 outstr.resize(strLength);
1378 MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
1385 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1388 for (
int i = 0; i < GetNumLevels(); ++i)
1389 Levels_[i]->print(out, verbLevel);
1392 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1394 isPreconditioner_ = flag;
1397 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1399 if (GetProcRankVerbose() != 0)
1401 #if defined(HAVE_MUELU_BOOST) && defined(HAVE_MUELU_BOOST_FOR_REAL) && defined(BOOST_VERSION) && (BOOST_VERSION >= 104400)
1406 dp.property(
"label", boost::get(boost::vertex_name, graph));
1407 dp.property(
"id", boost::get(boost::vertex_index, graph));
1408 dp.property(
"label", boost::get(boost::edge_name, graph));
1409 dp.property(
"color", boost::get(boost::edge_color, graph));
1412 std::map<const FactoryBase*, BoostVertex> vindices;
1413 typedef std::map<std::pair<BoostVertex, BoostVertex>, std::string> emap;
1416 static int call_id = 0;
1418 RCP<Operator> A = Levels_[0]->template Get<RCP<Operator> >(
"A");
1419 int rank = A->getDomainMap()->getComm()->getRank();
1422 for (
int i = currLevel; i <= currLevel + 1 && i < GetNumLevels(); i++) {
1424 Levels_[i]->UpdateGraph(vindices, edges, dp, graph);
1426 for (emap::const_iterator eit = edges.begin(); eit != edges.end(); eit++) {
1427 std::pair<BoostEdge, bool> boost_edge = boost::add_edge(eit->first.first, eit->first.second, graph);
1430 if (eit->second == std::string(
"Graph"))
1431 boost::put(
"label", dp, boost_edge.first, std::string(
"Graph_"));
1433 boost::put(
"label", dp, boost_edge.first, eit->second);
1435 boost::put(
"color", dp, boost_edge.first, std::string(
"red"));
1437 boost::put(
"color", dp, boost_edge.first, std::string(
"blue"));
1441 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"));
1442 boost::write_graphviz_dp(out, graph, dp, std::string(
"id"));
1446 GetOStream(
Errors) <<
"Dependency graph output requires boost and MueLu_ENABLE_Boost_for_real" << std::endl;
1451 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1456 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is not a matrix, skipping..." << std::endl;
1459 if (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A) != Teuchos::null) {
1460 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is a BlockedCrsMatrix, skipping..." << std::endl;
1468 if (A->getRowMap()->isSameAs(*(coords->getMap()))) {
1469 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: matrix and coordinates maps are same, skipping..." << std::endl;
1473 if (A->IsView(
"stridedMaps") && rcp_dynamic_cast<
const StridedMap>(A->getRowMap(
"stridedMaps")) != Teuchos::null) {
1474 RCP<const StridedMap> stridedRowMap = rcp_dynamic_cast<
const StridedMap>(A->getRowMap(
"stridedMaps"));
1478 Exceptions::RuntimeError,
"Hierarchy::ReplaceCoordinateMap: nontrivial maps (block id = " << stridedRowMap->getStridedBlockId() <<
", offset = " << stridedRowMap->getOffset() <<
")");
1481 GetOStream(
Runtime1) <<
"Replacing coordinate map" << std::endl;
1484 size_t blkSize = A->GetFixedBlockSize() / A->GetStorageBlockSize();
1490 GO indexBase = dofMap->getIndexBase();
1491 size_t numLocalDOFs = dofMap->getLocalNumElements();
1493 "Hierarchy::ReplaceCoordinateMap: block size (" << blkSize <<
") is incompatible with the number of local dofs in a row map (" << numLocalDOFs);
1496 Array<GO> nodeGIDs(numLocalDOFs / blkSize);
1497 for (
size_t i = 0; i < numLocalDOFs; i += blkSize)
1498 nodeGIDs[i / blkSize] = (GIDs[i] - indexBase) / blkSize + indexBase;
1501 nodeMap = MapFactory::Build(dofMap->lib(), INVALID, nodeGIDs(), indexBase, dofMap->getComm());
1507 if (coords->getLocalLength() != A->getRowMap()->getLocalNumElements()) {
1508 GetOStream(
Warnings) <<
"Coordinate vector does not match row map of matrix A!" << std::endl;
1514 std::vector<ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType> > coordData;
1515 for (
size_t i = 0; i < coords->getNumVectors(); i++) {
1516 coordData.
push_back(coords->getData(i));
1517 coordDataView.
push_back(coordData[i]());
1521 level.
Set(
"Coordinates", newCoords);
1524 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1526 int N = Levels_.size();
1527 if (((sizeOfAllocatedLevelMultiVectors_ == numvecs && residual_.size() == N) || numvecs <= 0) && !forceMapCheck)
return;
1530 if (residual_.size() != N) {
1531 DeleteLevelMultiVectors();
1533 residual_.resize(N);
1534 coarseRhs_.resize(N);
1536 coarseImport_.resize(N);
1537 coarseExport_.resize(N);
1538 correction_.resize(N);
1541 for (
int i = 0; i < N; i++) {
1542 RCP<Operator> A = Levels_[i]->template Get<RCP<Operator> >(
"A");
1548 if (!A_as_blocked.
is_null()) {
1549 Adm = A_as_blocked->getFullDomainMap();
1552 if (residual_[i].
is_null() || !residual_[i]->getMap()->isSameAs(*Arm))
1554 residual_[i] = MultiVectorFactory::Build(Arm, numvecs,
true);
1555 if (correction_[i].
is_null() || !correction_[i]->getMap()->isSameAs(*Adm))
1556 correction_[i] = MultiVectorFactory::Build(Adm, numvecs,
false);
1561 if (implicitTranspose_) {
1562 RCP<Operator> P = Levels_[i + 1]->template Get<RCP<Operator> >(
"P");
1565 if (coarseRhs_[i].
is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1566 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1569 RCP<Operator> R = Levels_[i + 1]->template Get<RCP<Operator> >(
"R");
1572 if (coarseRhs_[i].
is_null() || !coarseRhs_[i]->getMap()->isSameAs(*map))
1573 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1578 if (Levels_[i + 1]->IsAvailable(
"Importer"))
1579 importer = Levels_[i + 1]->template Get<RCP<const Import> >(
"Importer");
1580 if (doPRrebalance_ || importer.
is_null()) {
1582 if (coarseX_[i].
is_null() || !coarseX_[i]->getMap()->isSameAs(*map))
1583 coarseX_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1586 map = importer->getTargetMap();
1587 if (coarseImport_[i].
is_null() || !coarseImport_[i]->getMap()->isSameAs(*map)) {
1588 coarseImport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1589 coarseX_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1591 map = importer->getSourceMap();
1592 if (coarseExport_[i].
is_null() || !coarseExport_[i]->getMap()->isSameAs(*map))
1593 coarseExport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1597 sizeOfAllocatedLevelMultiVectors_ = numvecs;
1600 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1602 if (sizeOfAllocatedLevelMultiVectors_ == 0)
return;
1603 residual_.resize(0);
1604 coarseRhs_.resize(0);
1606 coarseImport_.resize(0);
1607 coarseExport_.resize(0);
1608 correction_.resize(0);
1609 sizeOfAllocatedLevelMultiVectors_ = 0;
1612 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1614 const LO startLevel,
const ConvData& conv)
const {
1615 return (startLevel == 0 && !isPreconditioner_ && (IsPrint(
Statistics1) || conv.
tol_ > 0));
1618 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1625 for (
LO k = 0; k < residualNorm.
size(); k++)
1626 if (residualNorm[k] >= convergenceTolerance)
1635 return convergenceStatus;
1638 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1642 << std::setiosflags(std::ios::left)
1643 << std::setprecision(3) << std::setw(4) << iteration
1645 << std::setprecision(10) << residualNorm
1649 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1651 const Operator& A,
const MultiVector& X,
const MultiVector& B,
const LO iteration,
1657 rate_ = currentResidualNorm / previousResidualNorm;
1658 previousResidualNorm = currentResidualNorm;
1661 PrintResidualHistory(iteration, residualNorm);
1663 return IsConverged(residualNorm, conv.
tol_);
1668 #endif // MUELU_HIERARCHY_DEF_HPP
Important warning messages (one line)
void Build(Level &level) const
Build an object with this factory.
void IsPreconditioner(const bool flag)
static Teuchos::RCP< Teuchos::StackedTimer > getStackedTimer()
High level timing information (use Teuchos::TimeMonitor::summarize() to print)
virtual std::string getObjectLabel() const
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level->Get< RCP<Matrix> >("A", factory) if factory == NULL => use default factory.
bool is_null(const boost::shared_ptr< T > &p)
RCP< Level > & GetLevel(const int levelID=0)
Retrieve a certain level from hierarchy.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
Hierarchy()
Default constructor.
virtual size_t getNodeSmootherComplexity() const =0
Compute a rough estimate of the cost to apply this smoother on this MPI rank. Return Teuchos::Ordinal...
void DeleteLevelMultiVectors()
void Release(const FactoryBase &factory)
Decrement the storage counter for all the inputs of a factory.
void AddLevel(const RCP< Level > &level)
Add a level at the end of the hierarchy.
void CheckForEmptySmoothersAndCoarseSolve()
void PrintResidualHistory(const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const
Print residualNorm for this iteration to the screen.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static void Write(const std::string &fileName, const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > &M)
void DumpCurrentGraph(int level) const
void AddNewLevel()
Add a new level at the end of the hierarchy.
VerbLevel toMueLuVerbLevel(const Teuchos::EVerbosityLevel verbLevel)
Translate Teuchos verbosity level to MueLu verbosity level.
void swap(RCP< T > &r_ptr)
One-liner description of what is happening.
virtual void Clean() const
void Clear(int startLevel=0)
Clear impermanent data from previous setup.
void SetPreviousLevel(const RCP< Level > &previousLevel)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Exception throws to report incompatible objects (like maps).
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
Set seed for random number generator.
Integrates Teuchos::TimeMonitor with MueLu verbosity system.
static std::string name()
static const NoFactory * get()
ConvergenceStatus Iterate(const MultiVector &B, MultiVector &X, ConvData conv=ConvData(), bool InitialGuessIsZero=false, LO startLevel=0)
Apply the multigrid preconditioner.
void ReplaceCoordinateMap(Level &level)
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
Call Level::AddKeepFlag for each level of the Hierarchy.
Print even more statistics.
void setlib(Xpetra::UnderlyingLib lib2)
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...
double GetSmootherComplexity() const
Base class for factories (e.g., R, P, and A_coarse).
Print statistics that do not involve significant additional computation.
Xpetra::UnderlyingLib lib_
Epetra/Tpetra mode.
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void Write(const LO &start=-1, const LO &end=-1, const std::string &suffix="")
Print matrices in the multigrid hierarchy to file.
Class that holds all level-specific information.
std::string description() const
Return a simple one-line description of this object.
Class that provides default factories within Needs class.
void RemoveKeepFlag(const std::string &ename, const FactoryBase *factory, KeepType keep=MueLu::All)
Call Level::RemoveKeepFlag for each level of the Hierarchy.
void CheckLevel(Level &level, int levelID)
Helper function.
Xpetra::UnderlyingLib lib()
virtual void setObjectLabel(const std::string &objectLabel)
void AllocateLevelMultiVectors(int numvecs, bool forceMapCheck=false)
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void Set(const std::string &ename, const T &entry, const FactoryBase *factory=NoFactory::get())
void SetMatvecParams(RCP< ParameterList > matvecParams)
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
int GetGlobalNumLevels() const
void push_back(const value_type &x)
RCP< Level > & GetPreviousLevel()
Previous level.
static Teuchos::Array< Magnitude > ResidualNorm(const Xpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, const MultiVector &X, const MultiVector &RHS)
Data struct for defining stopping criteria of multigrid iteration.
void SetLevelID(int levelID)
Set level number.
STS::magnitudeType MagnitudeType
void Delete(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Delete(ename, factory) for each level of the Hierarchy.
Timer to be used in non-factories.
Array< RCP< Level > > Levels_
Container for Level objects.
bool IsCalculationOfResidualRequired(const LO startLevel, const ConvData &conv) const
Decide if the residual needs to be computed.
int GetLevelID() const
Return level number.
Print class parameters (more parameters, more verbose)
virtual ~Hierarchy()
Destructor.
void residual(const Operator< SC, LO, GO, NO > &Aop, const MultiVector< SC, LO, GO, NO > &X_in, const MultiVector< SC, LO, GO, NO > &B_in, MultiVector< SC, LO, GO, NO > &R_in)
Exception throws to report errors in the internal logical of the program.
An exception safe way to call the method 'Level::SetFactoryManager()'.
void SetComm(RCP< const Teuchos::Comm< int > > const &comm)
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
Print all warning messages.
Description of what is happening (more verbose)
void Keep(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Keep(ename, factory) for each level of the Hierarchy.
ConvergenceStatus IsConverged(const Teuchos::Array< MagnitudeType > &residualNorm, const MagnitudeType convergenceTolerance) const
Decide if the multigrid iteration is converged.
double GetOperatorComplexity() const
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.
RCP< const Teuchos::Comm< int > > GetComm() const
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.
virtual std::string description() const
Return a simple one-line description of this object.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
std::string toString(const T &t)
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
void Request(const FactoryBase &factory)
Increment the storage counter for all the inputs of a factory.