10#ifndef MUELU_HIERARCHY_DEF_HPP
11#define MUELU_HIERARCHY_DEF_HPP
18#include <Xpetra_Matrix.hpp>
19#include <Xpetra_MultiVectorFactory.hpp>
20#include <Xpetra_Operator.hpp>
21#include <Xpetra_IO.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"
36#include "Teuchos_TimeMonitor.hpp"
40template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
51 ,
lib_(Xpetra::UseTpetra)
59template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
62 setObjectLabel(label);
63 Levels_[0]->setObjectLabel(label);
66template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
81 lib_ = A->getDomainMap()->
lib();
83 RCP<Level> Finest = rcp(
new Level);
89template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
92 setObjectLabel(label);
93 Levels_[0]->setObjectLabel(label);
96template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
99template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
103 if (level->GetLevelID() != -1 && (level->GetLevelID() != levelID))
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;
108 level->SetLevelID(levelID);
112 level->setObjectLabel(this->getObjectLabel());
115template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
118 newLevel->setlib(
lib_);
122template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
125 "MueLu::Hierarchy::GetLevel(): invalid input parameter value: LevelID = " << levelID);
129template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
134template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
136 RCP<Operator> A =
Levels_[0]->template Get<RCP<Operator>>(
"A");
137 RCP<const Teuchos::Comm<int>> comm = A->getDomainMap()->getComm();
141 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MAX, numLevels, Teuchos::ptr(&numGlobalLevels));
143 return numGlobalLevels;
146template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
148 double totalNnz = 0, lev0Nnz = 1;
151 "Operator complexity cannot be calculated because A is unavailable on level " << i);
152 RCP<Operator> A =
Levels_[i]->template Get<RCP<Operator>>(
"A");
156 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(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;
169template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
171 double node_sc = 0, global_sc = 0;
173 const size_t INVALID = Teuchos::OrdinalTraits<size_t>::invalid();
176 if (!
Levels_[0]->IsAvailable(
"A"))
return -1.0;
178 RCP<Operator> A =
Levels_[0]->template Get<RCP<Operator>>(
"A");
179 if (A.is_null())
return -1.0;
180 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
181 if (Am.is_null())
return -1.0;
182 a0_nnz = as<double>(Am->getGlobalNumEntries());
187 if (!
Levels_[i]->IsAvailable(
"PreSmoother"))
continue;
188 RCP<SmootherBase> S =
Levels_[i]->template Get<RCP<SmootherBase>>(
"PreSmoother");
189 if (S.is_null())
continue;
190 level_sc = S->getNodeSmootherComplexity();
191 if (level_sc == INVALID) {
196 node_sc += as<double>(level_sc);
200 RCP<const Teuchos::Comm<int>> comm = A->getDomainMap()->getComm();
201 Teuchos::reduceAll(*comm, Teuchos::REDUCE_SUM, node_sc, Teuchos::ptr(&global_sc));
202 Teuchos::reduceAll(*comm, Teuchos::REDUCE_MIN, node_sc, Teuchos::ptr(&min_sc));
207 return global_sc / a0_nnz;
211template <
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");
221template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
225 if (level->IsAvailable(
"A")) {
226 RCP<Operator> Aop = level->Get<RCP<Operator>>(
"A");
227 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Aop);
229 RCP<const Import> xpImporter = A->getCrsGraph()->getImporter();
230 if (!xpImporter.is_null())
231 xpImporter->setDistributorParameters(matvecParams);
232 RCP<const Export> xpExporter = A->getCrsGraph()->getExporter();
233 if (!xpExporter.is_null())
234 xpExporter->setDistributorParameters(matvecParams);
237 if (level->IsAvailable(
"P")) {
238 RCP<Matrix> P = level->Get<RCP<Matrix>>(
"P");
239 RCP<const Import> xpImporter = P->getCrsGraph()->getImporter();
240 if (!xpImporter.is_null())
241 xpImporter->setDistributorParameters(matvecParams);
242 RCP<const Export> xpExporter = P->getCrsGraph()->getExporter();
243 if (!xpExporter.is_null())
244 xpExporter->setDistributorParameters(matvecParams);
246 if (level->IsAvailable(
"R")) {
247 RCP<Matrix> R = level->Get<RCP<Matrix>>(
"R");
248 RCP<const Import> xpImporter = R->getCrsGraph()->getImporter();
249 if (!xpImporter.is_null())
250 xpImporter->setDistributorParameters(matvecParams);
251 RCP<const Export> xpExporter = R->getCrsGraph()->getExporter();
252 if (!xpExporter.is_null())
253 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.");
281 TimeMonitor m2(*
this, label + this->
ShortClassName() +
": " +
"Setup" +
" (total, level=" + Teuchos::toString(coarseLevelID) +
")");
285 "MueLu::Hierarchy::Setup(): argument coarseLevelManager cannot be null");
294 bool isFinestLevel = (fineLevelManager.is_null());
295 bool isLastLevel = (nextLevelManager.is_null());
299 RCP<Operator> A = level.
Get<RCP<Operator>>(
"A");
300 RCP<const Map> domainMap = A->getDomainMap();
301 RCP<const Teuchos::Comm<int>> comm = domainMap->getComm();
312 lib_ = domainMap->lib();
326 RCP<SetFactoryManager> SFMFine;
330 if (isFinestLevel &&
Levels_[coarseLevelID]->IsAvailable(
"Coordinates"))
339 RCP<TopSmootherFactory> coarseFact;
340 RCP<TopSmootherFactory> smootherFact = rcp(
new TopSmootherFactory(coarseLevelManager,
"Smoother"));
342 int nextLevelID = coarseLevelID + 1;
344 RCP<SetFactoryManager> SFMNext;
345 if (isLastLevel ==
false) {
382 if (coarseFact.is_null())
391 RCP<Operator> Ac = Teuchos::null;
392 TopRAPFactory coarseRAPFactory(fineLevelManager, coarseLevelManager);
395 Ac = level.
Get<RCP<Operator>>(
"A");
396 }
else if (!isFinestLevel) {
401 bool setLastLevelviaMaxCoarseSize =
false;
403 Ac = level.
Get<RCP<Operator>>(
"A");
404 RCP<Matrix> Acm = rcp_dynamic_cast<Matrix>(Ac);
408 level.
SetComm(Ac->getDomainMap()->getComm());
411 bool isOrigLastLevel = isLastLevel;
416 }
else if (Ac.is_null()) {
423 if (!Acm.is_null() && Acm->getGlobalNumRows() <=
maxCoarseSize_) {
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");
433 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(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) {
457 if (coarseFact.is_null())
465 coarseFact->Build(level);
476 smootherFact->Build(level);
481 if (isLastLevel ==
true) {
482 int actualNumLevels = nextLevelID;
483 if (isOrigLastLevel ==
false) {
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);
505 if (coarseFact.is_null())
507 Levels_[nextLevelID - 2]->Request(*coarseFact);
508 if (!(
Levels_[nextLevelID - 2]->
template Get<RCP<Matrix>>(
"A").is_null()))
509 coarseFact->Build(*(
Levels_[nextLevelID - 2]));
510 Levels_[nextLevelID - 2]->Release(*coarseFact);
512 Levels_.resize(actualNumLevels);
519 if (!isFinestLevel) {
523 level.
Release(coarseRAPFactory);
532template <
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++)
549 for (levelID = startLevel; levelID < numLevels;) {
550 bool r =
Setup(levelID,
553 (levelID + 1 != numLevels ?
levelManagers_[levelID + 1] : Teuchos::null));
574template <
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.");
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();
597 RCP<Matrix> Amat = rcp_dynamic_cast<Matrix>(A);
599 if (!Amat.is_null()) {
600 RCP<ParameterList> params = rcp(
new ParameterList());
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;
610 RCP<const FactoryManagerBase> rcpmanager = rcpFromRef(manager);
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);
637 "MueLu::Hierarchy::Setup(): number of level");
648template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
650 for (LO levelNo = 0; levelNo <
GetNumLevels(); ++levelNo) {
652 if ((level->IsAvailable(
"A") && !level->template Get<RCP<Operator>>(
"A").is_null()) && (!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;
658template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
663 for (
int iLevel = startLevel; iLevel <
GetNumLevels(); iLevel++)
667template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
674#if defined(HAVE_MUELU_EXPERIMENTAL) && defined(HAVE_MUELU_ADDITIVE_VARIANT)
675template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
677 bool InitialGuessIsZero, LO startLevel) {
678 LO nIts = conv.maxIts_;
679 MagnitudeType tol = conv.tol_;
681 std::string prefix = this->ShortClassName() +
": ";
682 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
683 std::string levelSuffix1 =
" (level=" +
toString(startLevel + 1) +
")";
686 RCP<Time> CompTime = Teuchos::TimeMonitor::getNewCounter(prefix +
"Computational Time (total)");
687 RCP<Time> Concurrent = Teuchos::TimeMonitor::getNewCounter(prefix +
"Concurrent portion");
688 RCP<Time> ApplyR = Teuchos::TimeMonitor::getNewCounter(prefix +
"R: Computational Time");
689 RCP<Time> ApplyPbar = Teuchos::TimeMonitor::getNewCounter(prefix +
"Pbar: Computational Time");
690 RCP<Time> CompFine = Teuchos::TimeMonitor::getNewCounter(prefix +
"Fine: Computational Time");
691 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix +
"Coarse: Computational Time");
692 RCP<Time> ApplySum = Teuchos::TimeMonitor::getNewCounter(prefix +
"Sum: Computational Time");
693 RCP<Time> Synchronize_beginning = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_beginning");
694 RCP<Time> Synchronize_center = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_center");
695 RCP<Time> Synchronize_end = Teuchos::TimeMonitor::getNewCounter(prefix +
"Synchronize_end");
697 RCP<Level> Fine = Levels_[0];
700 RCP<Operator> A = Fine->Get<RCP<Operator>>(
"A");
701 Teuchos::RCP<const Teuchos::Comm<int>> communicator = A->getDomainMap()->getComm();
709 SC one = STS::one(), zero = STS::zero();
711 bool zeroGuess = InitialGuessIsZero;
719 RCP<MultiVector> coarseRhs, coarseX;
721 RCP<SmootherBase> preSmoo_coarse, postSmoo_coarse;
722 bool emptyCoarseSolve =
true;
723 RCP<MultiVector> coarseX_prolonged = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
725 RCP<const Import> importer;
727 if (Levels_.size() > 1) {
729 if (Coarse->IsAvailable(
"Importer"))
730 importer = Coarse->Get<RCP<const Import>>(
"Importer");
732 R = Coarse->Get<RCP<Operator>>(
"R");
733 P = Coarse->Get<RCP<Operator>>(
"P");
736 Pbar = Coarse->Get<RCP<Operator>>(
"Pbar");
738 coarseRhs = MultiVectorFactory::Build(R->getRangeMap(), B.getNumVectors(),
true);
740 Ac = Coarse->Get<RCP<Operator>>(
"A");
743 R->apply(B, *coarseRhs, Teuchos::NO_TRANS, one, zero);
747 if (doPRrebalance_ || importer.is_null()) {
748 coarseX = MultiVectorFactory::Build(coarseRhs->getMap(), X.getNumVectors(),
true);
751 RCP<TimeMonitor> ITime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import (total)",
Timings0));
752 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import" + levelSuffix1,
Timings0));
755 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getTargetMap(), coarseRhs->getNumVectors());
756 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
757 coarseRhs.swap(coarseTmp);
759 coarseX = MultiVectorFactory::Build(importer->getTargetMap(), X.getNumVectors(),
true);
762 if (Coarse->IsAvailable(
"PreSmoother"))
763 preSmoo_coarse = Coarse->Get<RCP<SmootherBase>>(
"PreSmoother");
764 if (Coarse->IsAvailable(
"PostSmoother"))
765 postSmoo_coarse = Coarse->Get<RCP<SmootherBase>>(
"PostSmoother");
770 MagnitudeType prevNorm = STS::magnitude(STS::one()), curNorm = STS::magnitude(STS::one());
773 for (LO i = 1; i <= nIts; i++) {
774#ifdef HAVE_MUELU_DEBUG
775 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
776 std::ostringstream ss;
777 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
781 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
782 std::ostringstream ss;
783 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
789 bool emptyFineSolve =
true;
791 RCP<MultiVector> fineX;
792 fineX = MultiVectorFactory::Build(X.getMap(), X.getNumVectors(),
true);
801 if (Fine->IsAvailable(
"PreSmoother")) {
802 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>(
"PreSmoother");
804 preSmoo->Apply(*fineX, B, zeroGuess);
806 emptyFineSolve =
false;
808 if (Fine->IsAvailable(
"PostSmoother")) {
809 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>(
"PostSmoother");
811 postSmoo->Apply(*fineX, B, zeroGuess);
814 emptyFineSolve =
false;
816 if (emptyFineSolve ==
true) {
818 fineX->update(one, B, zero);
821 if (Levels_.size() > 1) {
823 if (Coarse->IsAvailable(
"PreSmoother")) {
825 preSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
827 emptyCoarseSolve =
false;
829 if (Coarse->IsAvailable(
"PostSmoother")) {
831 postSmoo_coarse->Apply(*coarseX, *coarseRhs, zeroGuess);
833 emptyCoarseSolve =
false;
835 if (emptyCoarseSolve ==
true) {
837 coarseX->update(one, *coarseRhs, zero);
844 if (!doPRrebalance_ && !importer.is_null()) {
845 RCP<TimeMonitor> ITime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export (total)",
Timings0));
846 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export" + levelSuffix1,
Timings0));
849 RCP<MultiVector> coarseTmp = MultiVectorFactory::Build(importer->getSourceMap(), coarseX->getNumVectors());
850 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
851 coarseX.swap(coarseTmp);
855 Pbar->apply(*coarseX, *coarseX_prolonged, Teuchos::NO_TRANS, one, zero);
860 X.update(1.0, *fineX, 1.0, *coarseX_prolonged, 0.0);
871template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
873 bool InitialGuessIsZero, LO startLevel) {
889 RCP<Level> Fine =
Levels_[startLevel];
893 std::string levelSuffix =
" (level=" +
toString(startLevel) +
")";
894 std::string levelSuffix1 =
" (level=" +
toString(startLevel + 1) +
")";
896 bool useStackedTimer = !Teuchos::TimeMonitor::getStackedTimer().is_null();
898 RCP<Monitor> iterateTime;
899 RCP<TimeMonitor> iterateTime1;
902 else if (!useStackedTimer)
905 std::string iterateLevelTimeLabel = prefix +
"Solve" + levelSuffix;
906 RCP<TimeMonitor> iterateLevelTime = rcp(
new TimeMonitor(*
this, iterateLevelTimeLabel,
Timings0));
908 bool zeroGuess = InitialGuessIsZero;
910 RCP<Operator> A = Fine->Get<RCP<Operator>>(
"A");
912 RCP<Time> CompCoarse = Teuchos::TimeMonitor::getNewCounter(prefix +
"Coarse: Computational Time");
925 const BlockedMultiVector* Bblocked =
dynamic_cast<const BlockedMultiVector*
>(&B);
927 ((Bblocked && !Bblocked->isSameSize(*
residual_[startLevel])) ||
928 (!Bblocked && !
residual_[startLevel]->isSameSize(B))))
933 typedef Teuchos::ScalarTraits<typename STS::magnitudeType> STM;
939 return convergenceStatus;
942 SC one = STS::one(), zero = STS::zero();
943 for (LO iteration = 1; iteration <= nIts; iteration++) {
944#ifdef HAVE_MUELU_DEBUG
946 if (A->getDomainMap()->isCompatible(*(X.getMap())) ==
false) {
947 std::ostringstream ss;
948 ss <<
"Level " << startLevel <<
": level A's domain map is not compatible with X";
952 if (A->getRangeMap()->isCompatible(*(B.getMap())) ==
false) {
953 std::ostringstream ss;
954 ss <<
"Level " << startLevel <<
": level A's range map is not compatible with B";
960 if (startLevel == as<LO>(
Levels_.size()) - 1) {
962 RCP<TimeMonitor> CLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : coarse" + levelSuffix,
Timings0));
964 bool emptySolve =
true;
967 if (Fine->IsAvailable(
"PreSmoother")) {
968 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>(
"PreSmoother");
970 preSmoo->Apply(X, B, zeroGuess);
975 if (Fine->IsAvailable(
"PostSmoother")) {
976 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>(
"PostSmoother");
978 postSmoo->Apply(X, B, zeroGuess);
983 if (emptySolve ==
true) {
985 X.update(one, B, zero);
990 RCP<Level> Coarse =
Levels_[startLevel + 1];
993 RCP<TimeMonitor> STime;
994 if (!useStackedTimer)
996 RCP<TimeMonitor> SLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : smoothing" + levelSuffix,
Timings0));
998 if (Fine->IsAvailable(
"PreSmoother")) {
999 RCP<SmootherBase> preSmoo = Fine->Get<RCP<SmootherBase>>(
"PreSmoother");
1000 preSmoo->Apply(X, B, zeroGuess);
1005 RCP<MultiVector> residual;
1007 RCP<TimeMonitor> ATime;
1008 if (!useStackedTimer)
1009 ATime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : residual calculation (total)",
Timings0));
1010 RCP<TimeMonitor> ALevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : residual calculation" + levelSuffix,
Timings0));
1021 RCP<Operator> P = Coarse->Get<RCP<Operator>>(
"P");
1022 if (Coarse->IsAvailable(
"Pbar"))
1023 P = Coarse->Get<RCP<Operator>>(
"Pbar");
1025 RCP<MultiVector> coarseRhs, coarseX;
1029 RCP<TimeMonitor> RTime;
1030 if (!useStackedTimer)
1032 RCP<TimeMonitor> RLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : restriction" + levelSuffix,
Timings0));
1036 P->apply(*residual, *coarseRhs, Teuchos::TRANS, one, zero);
1039 RCP<Operator> R = Coarse->Get<RCP<Operator>>(
"R");
1040 R->apply(*residual, *coarseRhs, Teuchos::NO_TRANS, one, zero);
1044 RCP<const Import> importer;
1045 if (Coarse->IsAvailable(
"Importer"))
1046 importer = Coarse->Get<RCP<const Import>>(
"Importer");
1050 RCP<TimeMonitor> ITime;
1051 if (!useStackedTimer)
1053 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : import" + levelSuffix1,
Timings0));
1057 coarseTmp->doImport(*coarseRhs, *importer, Xpetra::INSERT);
1058 coarseRhs.swap(coarseTmp);
1061 RCP<Operator> Ac = Coarse->Get<RCP<Operator>>(
"A");
1062 if (!Ac.is_null()) {
1063 RCP<const Map> origXMap = coarseX->getMap();
1064 RCP<const Map> origRhsMap = coarseRhs->getMap();
1067 coarseRhs->replaceMap(Ac->getRangeMap());
1068 coarseX->replaceMap(Ac->getDomainMap());
1071 iterateLevelTime = Teuchos::null;
1073 Iterate(*coarseRhs, *coarseX, 1,
true, startLevel + 1);
1076 Iterate(*coarseRhs, *coarseX, 1,
false, startLevel + 1);
1079 iterateLevelTime = rcp(
new TimeMonitor(*
this, iterateLevelTimeLabel));
1081 coarseX->replaceMap(origXMap);
1082 coarseRhs->replaceMap(origRhsMap);
1086 RCP<TimeMonitor> ITime;
1087 if (!useStackedTimer)
1089 RCP<TimeMonitor> ILevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : export" + levelSuffix1,
Timings0));
1093 coarseTmp->doExport(*coarseX, *importer, Xpetra::INSERT);
1094 coarseX.swap(coarseTmp);
1099 RCP<TimeMonitor> PTime;
1100 if (!useStackedTimer)
1102 RCP<TimeMonitor> PLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : prolongation" + levelSuffix,
Timings0));
1110 RCP<MultiVector> correction =
correction_[startLevel];
1111 P->apply(*coarseX, *correction, Teuchos::NO_TRANS, one, zero);
1118 RCP<TimeMonitor> STime;
1119 if (!useStackedTimer)
1121 RCP<TimeMonitor> SLevelTime = rcp(
new TimeMonitor(*
this, prefix +
"Solve : smoothing" + levelSuffix,
Timings0));
1123 if (Fine->IsAvailable(
"PostSmoother")) {
1124 RCP<SmootherBase> postSmoo = Fine->Get<RCP<SmootherBase>>(
"PostSmoother");
1125 postSmoo->Apply(X, B,
false);
1134 return convergenceStatus;
1141template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1143 LO startLevel = (start != -1 ? start : 0);
1144 LO endLevel = (end != -1 ? end :
Levels_.size() - 1);
1147 "MueLu::Hierarchy::Write : startLevel must be <= endLevel");
1150 "MueLu::Hierarchy::Write bad start or end level");
1152 for (LO i = startLevel; i < endLevel + 1; i++) {
1153 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(
Levels_[i]->
template Get<RCP<Operator>>(
"A")), P, R;
1155 P = rcp_dynamic_cast<Matrix>(
Levels_[i]->
template Get<RCP<Operator>>(
"P"));
1157 R = rcp_dynamic_cast<Matrix>(
Levels_[i]->
template Get<RCP<Operator>>(
"R"));
1160 if (!A.is_null()) Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"A_" +
toString(i) + suffix +
".m", *A);
1162 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"P_" +
toString(i) + suffix +
".m", *P);
1165 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(
"R_" +
toString(i) + suffix +
".m", *R);
1170template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1172 for (Array<RCP<Level>>::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1173 (*it)->Keep(ename, factory);
1176template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1178 for (Array<RCP<Level>>::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1179 (*it)->Delete(ename, factory);
1182template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1184 for (Array<RCP<Level>>::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1185 (*it)->AddKeepFlag(ename, factory, keep);
1188template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1190 for (Array<RCP<Level>>::iterator it =
Levels_.begin(); it !=
Levels_.end(); ++it)
1191 (*it)->RemoveKeepFlag(ename, factory, keep);
1194template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1197 std::ostringstream out;
1205template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1210template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1212 RCP<Operator> A0 =
Levels_[0]->template Get<RCP<Operator>>(
"A");
1213 RCP<const Teuchos::Comm<int>> comm = A0->getDomainMap()->getComm();
1216 RCP<Operator> Ac =
Levels_[numLevels - 1]->template Get<RCP<Operator>>(
"A");
1223 int root = comm->getRank();
1226 int smartData = numLevels * comm->getSize() + comm->getRank(), maxSmartData;
1227 reduceAll(*comm, Teuchos::REDUCE_MAX, smartData, Teuchos::ptr(&maxSmartData));
1228 root = maxSmartData % comm->getSize();
1232 double smoother_comp = -1.0;
1238 std::vector<Xpetra::global_size_t> nnzPerLevel;
1239 std::vector<Xpetra::global_size_t> rowsPerLevel;
1240 std::vector<int> numProcsPerLevel;
1241 bool someOpsNotMatrices =
false;
1242 const Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1243 for (
int i = 0; i < numLevels; i++) {
1245 "Operator A is not available on level " << i);
1247 RCP<Operator> A =
Levels_[i]->template Get<RCP<Operator>>(
"A");
1249 "Operator A on level " << i <<
" is null.");
1251 RCP<Matrix> Am = rcp_dynamic_cast<Matrix>(A);
1253 someOpsNotMatrices =
true;
1254 nnzPerLevel.push_back(INVALID);
1255 rowsPerLevel.push_back(A->getDomainMap()->getGlobalNumElements());
1256 numProcsPerLevel.push_back(A->getDomainMap()->getComm()->getSize());
1258 LO storageblocksize = Am->GetStorageBlockSize();
1259 Xpetra::global_size_t nnz = Am->getGlobalNumEntries() * storageblocksize * storageblocksize;
1260 nnzPerLevel.push_back(nnz);
1261 rowsPerLevel.push_back(Am->getGlobalNumRows() * storageblocksize);
1262 numProcsPerLevel.push_back(Am->getRowMap()->getComm()->getSize());
1265 if (someOpsNotMatrices)
1266 GetOStream(
Warnings0) <<
"Some level operators are not matrices, statistics calculation are incomplete" << std::endl;
1269 std::string label =
Levels_[0]->getObjectLabel();
1270 std::ostringstream oss;
1271 oss << std::setfill(
' ');
1272 oss <<
"\n--------------------------------------------------------------------------------\n";
1273 oss <<
"--- Multigrid Summary " << std::setw(28) << std::left << label <<
"---\n";
1274 oss <<
"--------------------------------------------------------------------------------" << std::endl;
1276 oss <<
"Scalar = " << Teuchos::ScalarTraits<Scalar>::name() << std::endl;
1277 oss <<
"Number of levels = " << numLevels << std::endl;
1278 oss <<
"Operator complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1279 if (!someOpsNotMatrices)
1282 oss <<
"not available (Some operators in hierarchy are not matrices.)" << std::endl;
1284 if (smoother_comp != -1.0) {
1285 oss <<
"Smoother complexity = " << std::setprecision(2) << std::setiosflags(std::ios::fixed)
1286 << smoother_comp << std::endl;
1291 oss <<
"Cycle type = V" << std::endl;
1294 oss <<
"Cycle type = W" << std::endl;
1303 Xpetra::global_size_t tt = rowsPerLevel[0];
1309 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1310 tt = nnzPerLevel[i];
1320 tt = numProcsPerLevel[0];
1326 oss <<
"level " << std::setw(rowspacer) <<
" rows " << std::setw(nnzspacer) <<
" nnz "
1327 <<
" nnz/row" << std::setw(npspacer) <<
" c ratio"
1328 <<
" procs" << std::endl;
1329 for (
size_t i = 0; i < nnzPerLevel.size(); ++i) {
1330 oss <<
" " << i <<
" ";
1331 oss << std::setw(rowspacer) << rowsPerLevel[i];
1332 if (nnzPerLevel[i] != INVALID) {
1333 oss << std::setw(nnzspacer) << nnzPerLevel[i];
1334 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1335 oss << std::setw(9) << as<double>(nnzPerLevel[i]) / rowsPerLevel[i];
1337 oss << std::setw(nnzspacer) <<
"Operator";
1338 oss << std::setprecision(2) << std::setiosflags(std::ios::fixed);
1339 oss << std::setw(9) <<
" ";
1342 oss << std::setw(9) << as<double>(rowsPerLevel[i - 1]) / rowsPerLevel[i];
1344 oss << std::setw(9) <<
" ";
1345 oss <<
" " << std::setw(npspacer) << numProcsPerLevel[i] << std::endl;
1349 RCP<SmootherBase> preSmoo, postSmoo;
1350 if (
Levels_[i]->IsAvailable(
"PreSmoother"))
1351 preSmoo =
Levels_[i]->template Get<RCP<SmootherBase>>(
"PreSmoother");
1352 if (
Levels_[i]->IsAvailable(
"PostSmoother"))
1353 postSmoo =
Levels_[i]->template Get<RCP<SmootherBase>>(
"PostSmoother");
1355 if (preSmoo != null && preSmoo == postSmoo)
1356 oss <<
"Smoother (level " << i <<
") both : " << preSmoo->description() << std::endl;
1358 oss <<
"Smoother (level " << i <<
") pre : "
1359 << (preSmoo != null ? preSmoo->description() :
"no smoother") << std::endl;
1360 oss <<
"Smoother (level " << i <<
") post : "
1361 << (postSmoo != null ? postSmoo->description() :
"no smoother") << std::endl;
1372 RCP<const Teuchos::MpiComm<int>> mpiComm = rcp_dynamic_cast<const Teuchos::MpiComm<int>>(comm);
1373 MPI_Comm rawComm = (*mpiComm->getRawMpiComm())();
1375 int strLength = outstr.size();
1376 MPI_Bcast(&strLength, 1, MPI_INT, root, rawComm);
1377 if (comm->getRank() != root)
1378 outstr.resize(strLength);
1379 MPI_Bcast(&outstr[0], strLength, MPI_CHAR, root, rawComm);
1386template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1388 Teuchos::OSTab tab2(out);
1393template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1398template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1402#if defined(HAVE_MUELU_BOOST) && defined(HAVE_MUELU_BOOST_FOR_REAL) && defined(BOOST_VERSION) && (BOOST_VERSION >= 104400)
1407 dp.property(
"label", boost::get(boost::vertex_name, graph));
1408 dp.property(
"id", boost::get(boost::vertex_index, graph));
1409 dp.property(
"label", boost::get(boost::edge_name, graph));
1410 dp.property(
"color", boost::get(boost::edge_color, graph));
1413 std::map<const FactoryBase*, BoostVertex> vindices;
1414 typedef std::map<std::pair<BoostVertex, BoostVertex>, std::string> emap;
1417 static int call_id = 0;
1419 RCP<Operator> A =
Levels_[0]->template Get<RCP<Operator>>(
"A");
1420 int rank = A->getDomainMap()->getComm()->getRank();
1423 for (
int i = currLevel; i <= currLevel + 1 && i <
GetNumLevels(); i++) {
1425 Levels_[i]->UpdateGraph(vindices, edges, dp, graph);
1427 for (emap::const_iterator eit = edges.begin(); eit != edges.end(); eit++) {
1428 std::pair<BoostEdge, bool> boost_edge = boost::add_edge(eit->first.first, eit->first.second, graph);
1431 if (eit->second == std::string(
"Graph"))
1432 boost::put(
"label", dp, boost_edge.first, std::string(
"Graph_"));
1434 boost::put(
"label", dp, boost_edge.first, eit->second);
1436 boost::put(
"color", dp, boost_edge.first, std::string(
"red"));
1438 boost::put(
"color", dp, boost_edge.first, std::string(
"blue"));
1442 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"));
1443 boost::write_graphviz_dp(out, graph, dp, std::string(
"id"));
1447 GetOStream(
Errors) <<
"Dependency graph output requires boost and MueLu_ENABLE_Boost_for_real" << std::endl;
1452template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1454 RCP<Operator> Ao = level.
Get<RCP<Operator>>(
"A");
1455 RCP<Matrix> A = rcp_dynamic_cast<Matrix>(Ao);
1457 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is not a matrix, skipping..." << std::endl;
1460 if (Teuchos::rcp_dynamic_cast<BlockedCrsMatrix>(A) != Teuchos::null) {
1461 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: operator is a BlockedCrsMatrix, skipping..." << std::endl;
1465 typedef Xpetra::MultiVector<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LO, GO, NO> xdMV;
1467 RCP<xdMV> coords = level.
Get<RCP<xdMV>>(
"Coordinates");
1469 if (A->getRowMap()->isSameAs(*(coords->getMap()))) {
1470 GetOStream(
Runtime1) <<
"Hierarchy::ReplaceCoordinateMap: matrix and coordinates maps are same, skipping..." << std::endl;
1474 if (A->IsView(
"stridedMaps") && rcp_dynamic_cast<const StridedMap>(A->getRowMap(
"stridedMaps")) != Teuchos::null) {
1475 RCP<const StridedMap> stridedRowMap = rcp_dynamic_cast<const StridedMap>(A->getRowMap(
"stridedMaps"));
1478 TEUCHOS_TEST_FOR_EXCEPTION(stridedRowMap->getStridedBlockId() != -1 || stridedRowMap->getOffset() != 0,
1479 Exceptions::RuntimeError,
"Hierarchy::ReplaceCoordinateMap: nontrivial maps (block id = " << stridedRowMap->getStridedBlockId() <<
", offset = " << stridedRowMap->getOffset() <<
")");
1483 TEUCHOS_TEST_FOR_EXCEPTION(A->GetFixedBlockSize() % A->GetStorageBlockSize() != 0,
Exceptions::RuntimeError,
"Hierarchy::ReplaceCoordinateMap: Storage block size does not evenly divide fixed block size");
1485 size_t blkSize = A->GetFixedBlockSize() / A->GetStorageBlockSize();
1487 RCP<const Map> nodeMap = A->getRowMap();
1490 RCP<const Map> dofMap = A->getRowMap();
1491 GO indexBase = dofMap->getIndexBase();
1492 size_t numLocalDOFs = dofMap->getLocalNumElements();
1494 "Hierarchy::ReplaceCoordinateMap: block size (" << blkSize <<
") is incompatible with the number of local dofs in a row map (" << numLocalDOFs);
1495 ArrayView<const GO> GIDs = dofMap->getLocalElementList();
1497 Array<GO> nodeGIDs(numLocalDOFs / blkSize);
1498 for (
size_t i = 0; i < numLocalDOFs; i += blkSize)
1499 nodeGIDs[i / blkSize] = (GIDs[i] - indexBase) / blkSize + indexBase;
1501 Xpetra::global_size_t INVALID = Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid();
1502 nodeMap = MapFactory::Build(dofMap->lib(), INVALID, nodeGIDs(), indexBase, dofMap->getComm());
1508 if (coords->getLocalLength() != A->getRowMap()->getLocalNumElements()) {
1509 GetOStream(
Warnings) <<
"Coordinate vector does not match row map of matrix A!" << std::endl;
1514 Array<ArrayView<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>> coordDataView;
1515 std::vector<ArrayRCP<const typename Teuchos::ScalarTraits<Scalar>::coordinateType>> coordData;
1516 for (
size_t i = 0; i < coords->getNumVectors(); i++) {
1517 coordData.push_back(coords->getData(i));
1518 coordDataView.push_back(coordData[i]());
1521 RCP<xdMV> newCoords = Xpetra::MultiVectorFactory<typename Teuchos::ScalarTraits<Scalar>::coordinateType, LO, GO, NO>::Build(nodeMap, coordDataView(), coords->getNumVectors());
1522 level.
Set(
"Coordinates", newCoords);
1525template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1542 for (
int i = 0; i < N; i++) {
1543 RCP<Operator> A =
Levels_[i]->template Get<RCP<Operator>>(
"A");
1546 RCP<const BlockedCrsMatrix> A_as_blocked = Teuchos::rcp_dynamic_cast<const BlockedCrsMatrix>(A);
1547 RCP<const Map> Arm = A->getRangeMap();
1548 RCP<const Map> Adm = A->getDomainMap();
1549 if (!A_as_blocked.is_null()) {
1550 Adm = A_as_blocked->getFullDomainMap();
1555 residual_[i] = MultiVectorFactory::Build(Arm, numvecs,
true);
1557 correction_[i] = MultiVectorFactory::Build(Adm, numvecs,
false);
1563 RCP<Operator> P =
Levels_[i + 1]->template Get<RCP<Operator>>(
"P");
1565 RCP<const Map> map = P->getDomainMap();
1567 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1570 RCP<Operator> R =
Levels_[i + 1]->template Get<RCP<Operator>>(
"R");
1572 RCP<const Map> map = R->getRangeMap();
1574 coarseRhs_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1578 RCP<const Import> importer;
1579 if (
Levels_[i + 1]->IsAvailable(
"Importer"))
1580 importer =
Levels_[i + 1]->template Get<RCP<const Import>>(
"Importer");
1582 RCP<const Map> map =
coarseRhs_[i]->getMap();
1584 coarseX_[i] = MultiVectorFactory::Build(map, numvecs,
true);
1587 map = importer->getTargetMap();
1589 coarseImport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1590 coarseX_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1592 map = importer->getSourceMap();
1594 coarseExport_[i] = MultiVectorFactory::Build(map, numvecs,
false);
1601template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1613template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1615 const LO startLevel,
const ConvData& conv)
const {
1619template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1621 const Teuchos::Array<MagnitudeType>& residualNorm,
const MagnitudeType convergenceTolerance)
const {
1624 if (convergenceTolerance > Teuchos::ScalarTraits<MagnitudeType>::zero()) {
1626 for (LO k = 0; k < residualNorm.size(); k++)
1627 if (residualNorm[k] >= convergenceTolerance)
1636 return convergenceStatus;
1639template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1641 const LO iteration,
const Teuchos::Array<MagnitudeType>& residualNorm)
const {
1643 << std::setiosflags(std::ios::left)
1644 << std::setprecision(3) << std::setw(4) << iteration
1646 << std::setprecision(10) << residualNorm
1650template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
1652 const Operator& A,
const MultiVector& X,
const MultiVector& B,
const LO iteration,
1654 Teuchos::Array<MagnitudeType> residualNorm;
1658 rate_ = currentResidualNorm / previousResidualNorm;
1659 previousResidualNorm = currentResidualNorm;
virtual std::string ShortClassName() const
Return the class name of the object, without template parameters and without namespace.
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
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.
Array< RCP< MultiVector > > residual_
virtual ~Hierarchy()
Destructor.
void CheckLevel(Level &level, int levelID)
Helper function.
int WCycleStartLevel_
Level at which to start W-cycle.
std::string description() const
Return a simple one-line description of this object.
void CheckForEmptySmoothersAndCoarseSolve()
std::string description_
cache description to avoid recreating in each call to description() - use ResetDescription() to force...
void IsPreconditioner(const bool flag)
static CycleType GetDefaultCycle()
Array< RCP< Level > > Levels_
Container for Level objects.
static bool GetDefaultPRrebalance()
Array< RCP< MultiVector > > coarseRhs_
Array< RCP< MultiVector > > coarseExport_
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.
void ResetDescription()
force recreation of cached description_ next time description() is called:
STS::magnitudeType MagnitudeType
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
bool isPreconditioner_
Hierarchy may be used in a standalone mode, or as a preconditioner.
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
int sizeOfAllocatedLevelMultiVectors_
Caching (Multi)Vectors used in Hierarchy::Iterate().
static bool GetDefaultImplicitTranspose()
void SetMatvecParams(RCP< ParameterList > matvecParams)
static Xpetra::global_size_t GetDefaultMaxCoarseSize()
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.
CycleType Cycle_
V- or W-cycle.
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
MagnitudeType rate_
Convergece rate.
void PrintResidualHistory(const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const
Print residualNorm for this iteration to the screen.
Array< RCP< MultiVector > > coarseX_
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.
double scalingFactor_
Scaling factor to be applied to coarse grid correction.
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.
Xpetra::UnderlyingLib lib()
static bool GetDefaultFuseProlongationAndUpdate()
bool doPRViaCopyrebalance_
Array< RCP< MultiVector > > coarseImport_
bool isDumpingEnabled_
Graph dumping.
Xpetra::global_size_t maxCoarseSize_
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
Call Level::AddKeepFlag for each level of the Hierarchy.
Array< RCP< MultiVector > > correction_
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.
bool fuseProlongationAndUpdate_
Array< RCP< const FactoryManagerBase > > levelManagers_
Level managers used during the Setup.
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, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
static void SetRandomSeed(const Teuchos::Comm< int > &comm)
static RCP< MultiVector > Residual(const Xpetra::Operator< Scalar, DefaultLocalOrdinal, DefaultGlobalOrdinal, DefaultNode > &Op, const MultiVector &X, const MultiVector &RHS)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
VerbLevel GetVerbLevel() const
Get the verbosity level.
int SetProcRankVerbose(int procRank) const
Set proc rank used for printing.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
int GetProcRankVerbose() const
Get proc rank used for printing. Do not use this information for any other purpose....
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.