10#ifndef MUELU_HIERARCHYMANAGER_DEF_HPP
11#define MUELU_HIERARCHYMANAGER_DEF_HPP
16#include <Teuchos_Array.hpp>
18#include <Xpetra_Operator.hpp>
19#include <Xpetra_IO.hpp>
25#include "MueLu_Aggregates.hpp"
26#include "MueLu_Hierarchy.hpp"
30#include "MueLu_PerfUtils.hpp"
32#ifdef HAVE_MUELU_INTREPID2
33#include "Kokkos_DynRankView.hpp"
38template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
40 : numDesiredLevel_(numDesiredLevel)
41 , maxCoarseSize_(
MasterList::getDefault<int>(
"coarse: max size"))
43 , doPRrebalance_(
MasterList::getDefault<bool>(
"repartition: rebalance P and R"))
44 , doPRViaCopyrebalance_(
MasterList::getDefault<bool>(
"repartition: explicit via new copy rebalance P and R"))
45 , implicitTranspose_(
MasterList::getDefault<bool>(
"transpose: use implicit"))
46 , fuseProlongationAndUpdate_(
MasterList::getDefault<bool>(
"fuse prolongation and update"))
47 , suppressNullspaceDimensionCheck_(
MasterList::getDefault<bool>(
"nullspace: suppress dimension check"))
48 , sizeOfMultiVectors_(
MasterList::getDefault<int>(
"number of vectors"))
49 , graphOutputLevel_(-2) {}
51template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
53 const int lastLevel = startLevel + numDesiredLevel - 1;
54 if (levelManagers_.size() < lastLevel + 1)
55 levelManagers_.resize(lastLevel + 1);
57 for (
int iLevel = startLevel; iLevel <= lastLevel; iLevel++)
58 levelManagers_[iLevel] = manager;
61template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
64 return (levelID >= levelManagers_.size() ? levelManagers_[levelManagers_.size() - 1] : levelManagers_[levelID]);
67template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
69 return levelManagers_.size();
72template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
74 for (
int i = 0; i < levelManagers_.size(); i++)
75 TEUCHOS_TEST_FOR_EXCEPTION(levelManagers_[i] == Teuchos::null, Exceptions::RuntimeError,
"MueLu:HierarchyConfig::CheckConfig(): Undefined configuration for level:");
78template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
80 return rcp(
new Hierarchy());
83template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
85 return rcp(
new Hierarchy(label));
88template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
90 TEUCHOS_TEST_FOR_EXCEPTION(!H.GetLevel(0)->IsAvailable(
"A"), Exceptions::RuntimeError,
"No fine level operator");
92 RCP<Level> l0 = H.GetLevel(0);
93 RCP<Operator> Op = l0->Get<RCP<Operator>>(
"A");
96 if (l0->IsAvailable(
"Nullspace")) {
97 RCP<Matrix> A = Teuchos::rcp_dynamic_cast<Matrix>(Op);
98 if (A != Teuchos::null) {
99 RCP<MultiVector> nullspace = l0->Get<RCP<MultiVector>>(
"Nullspace");
101 if (
static_cast<size_t>(A->GetFixedBlockSize()) > nullspace->getNumVectors()) {
102 std::stringstream msg;
103 msg <<
"User-provided nullspace has fewer vectors ("
104 << nullspace->getNumVectors() <<
") than number of PDE equations ("
105 << A->GetFixedBlockSize() <<
"). ";
107 if (suppressNullspaceDimensionCheck_) {
108 msg <<
"It depends on the PDE, if this is a problem or not.";
109 this->GetOStream(Warnings0) << msg.str() << std::endl;
111 msg <<
"Add the missing nullspace vectors! (You can suppress this check. See the MueLu user guide for details.)";
112 TEUCHOS_TEST_FOR_EXCEPTION(
static_cast<size_t>(A->GetFixedBlockSize()) > nullspace->getNumVectors(), Exceptions::RuntimeError, msg.str());
116 this->GetOStream(Warnings0) <<
"Skipping dimension check of user-supplied nullspace because user-supplied operator is not a matrix" << std::endl;
120#ifdef HAVE_MUELU_DEBUG
122 for (
int i = 0; i < levelManagers_.size(); i++)
123 levelManagers_[i]->ResetDebugData();
130 Xpetra::UnderlyingLib lib = Op->getDomainMap()->lib();
137 H.SetMaxCoarseSize(maxCoarseSize_);
138 VerboseObject::SetDefaultVerbLevel(verbosity_);
139 if (graphOutputLevel_ >= 0 || graphOutputLevel_ == -1)
140 H.EnableGraphDumping(
"dep_graph", graphOutputLevel_);
142 if (VerboseObject::IsPrint(Statistics2)) {
143 RCP<Matrix> Amat = rcp_dynamic_cast<Matrix>(Op);
145 if (!Amat.is_null()) {
146 RCP<ParameterList> params = rcp(
new ParameterList());
147 params->set(
"printLoadBalancingInfo",
true);
148 params->set(
"printCommInfo",
true);
150 VerboseObject::GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Amat,
"A0", params);
152 VerboseObject::GetOStream(Warnings1) <<
"Fine level operator is not a matrix, statistics are not available" << std::endl;
156 H.SetPRrebalance(doPRrebalance_);
157 H.SetPRViaCopyrebalance(doPRViaCopyrebalance_);
158 H.SetImplicitTranspose(implicitTranspose_);
159 H.SetFuseProlongationAndUpdate(fuseProlongationAndUpdate_);
180 for (
int i = 0; i < numDesiredLevel_; i++) {
181 std::map<int, std::vector<keep_pair>>::const_iterator it = keep_.find(i);
182 if (it != keep_.end()) {
183 RCP<Level> l = H.GetLevel(i);
184 const std::vector<keep_pair>& keeps = it->second;
185 for (
size_t j = 0; j < keeps.size(); j++)
186 l->Keep(keeps[j].first, keeps[j].second);
188 if (i < numDesiredLevel_ - 1) {
189 RCP<Level> newLevel = rcp(
new Level());
190 H.AddLevel(newLevel);
195 for (
auto iter = matricesToPrint_.begin(); iter != matricesToPrint_.end(); iter++)
196 ExportDataSetKeepFlags(H, iter->second, iter->first);
199 ExportDataSetKeepFlags(H, nullspaceToPrint_,
"Nullspace");
200 ExportDataSetKeepFlags(H, coordinatesToPrint_,
"Coordinates");
201 ExportDataSetKeepFlags(H, materialToPrint_,
"Material");
203 ExportDataSetKeepFlagsNextLevel(H, aggregatesToPrint_,
"Aggregates");
204#ifdef HAVE_MUELU_INTREPID2
205 ExportDataSetKeepFlags(H, elementToNodeMapsToPrint_,
"pcoarsen: element to node map");
209 for (
int i = 0; i < dataToSave_.size(); i++)
210 ExportDataSetKeepFlagsAll(H, dataToSave_[i]);
213 int lastLevelID = numDesiredLevel_ - 1;
214 bool isLastLevel =
false;
216 while (!isLastLevel) {
217 bool r = H.Setup(levelID,
218 LvlMngr(levelID - 1, lastLevelID),
219 LvlMngr(levelID, lastLevelID),
220 LvlMngr(levelID + 1, lastLevelID));
221 if (levelID < H.GetNumLevels())
222 H.GetLevel(levelID)->print(H.GetOStream(Developer), verbosity_);
224 isLastLevel = r || (levelID == lastLevelID);
227 if (!matvecParams_.is_null())
228 H.SetMatvecParams(matvecParams_);
229 H.AllocateLevelMultiVectors(sizeOfMultiVectors_);
233 H.describe(H.GetOStream(Runtime0), verbosity_);
234 H.CheckForEmptySmoothersAndCoarseSolve();
245 numDesiredLevel_ = levelID;
248 for (
auto iter = matricesToPrint_.begin(); iter != matricesToPrint_.end(); iter++) {
249 WriteData<Matrix>(H, iter->second, iter->first);
253 WriteData<MultiVector>(H, nullspaceToPrint_,
"Nullspace");
254 WriteData<MultiVector>(H, coordinatesToPrint_,
"Coordinates");
255 WriteData<MultiVector>(H, materialToPrint_,
"Material");
256 WriteDataAggregates(H, aggregatesToPrint_,
"Aggregates");
258#ifdef HAVE_MUELU_INTREPID2
259 typedef Kokkos::DynRankView<LocalOrdinal, typename Node::device_type> FCi;
260 WriteDataFC<FCi>(H, elementToNodeMapsToPrint_,
"pcoarsen: element to node map",
"el2node");
265template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
269 return Teuchos::null;
271 if (levelID == lastLevelID + 1)
272 return Teuchos::null;
274 if (levelManagers_.size() == 0) {
276 static RCP<FactoryManagerBase> defaultMngr = rcp(
new FactoryManager());
280 return GetFactoryManager(levelID);
283template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
285 for (
int i = 0; i < data.size(); ++i) {
286 if (data[i] < H.GetNumLevels()) {
287 RCP<Level> L = H.GetLevel(data[i]);
288 if (!L.is_null() && data[i] < levelManagers_.size())
289 L->AddKeepFlag(name, &*levelManagers_[data[i]]->GetFactory(name));
294template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
296 for (
int i = 0; i < data.size(); ++i) {
297 if (data[i] < H.GetNumLevels()) {
298 RCP<Level> L = H.GetLevel(data[i]);
299 if (!L.is_null() && data[i] + 1 < levelManagers_.size())
300 L->AddKeepFlag(name, &*levelManagers_[data[i] + 1]->GetFactory(name));
305template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
307 for (
int i = 0; i < H.GetNumLevels(); i++) {
308 RCP<Level> L = H.GetLevel(i);
309 if (!L.is_null() && i < levelManagers_.size())
310 L->AddKeepFlag(name, &*levelManagers_[i]->GetFactory(name));
314template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
317 for (
int i = 0; i < data.size(); ++i) {
318 std::string fileName;
319 if (H.getObjectLabel() !=
"")
320 fileName = H.getObjectLabel() +
"_" + name +
"_" + Teuchos::toString(data[i]) +
".m";
322 fileName = name +
"_" + Teuchos::toString(data[i]) +
".m";
324 if (data[i] < H.GetNumLevels()) {
325 RCP<Level> L = H.GetLevel(data[i]);
326 if (data[i] < levelManagers_.size() && L->IsAvailable(name, &*levelManagers_[data[i]]->GetFactory(name))) {
328 RCP<T> M = L->template Get<RCP<T>>(name, &*levelManagers_[data[i]]->GetFactory(name));
330 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(fileName, *M);
332 }
else if (L->IsAvailable(name)) {
334 RCP<T> M = L->template Get<RCP<T>>(name);
336 Xpetra::IO<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Write(fileName, *M);
343template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
345 for (
int i = 0; i < data.size(); ++i) {
346 const std::string fileName = name +
"_" + Teuchos::toString(data[i]) +
".m";
348 if (data[i] < H.GetNumLevels()) {
349 RCP<Level> L = H.GetLevel(data[i]);
353 if (data[i] + 1 < H.GetNumLevels() && L->IsAvailable(name, &*levelManagers_[data[i] + 1]->GetFactory(name))) {
355 agg = L->template Get<RCP<Aggregates>>(name, &*levelManagers_[data[i] + 1]->GetFactory(name));
356 }
else if (L->IsAvailable(name)) {
357 agg = L->template Get<RCP<Aggregates>>(
"Aggregates");
359 if (!agg.is_null()) {
360 std::ofstream ofs(fileName);
361 Teuchos::FancyOStream fofs(rcp(&ofs,
false));
362 agg->print(fofs, Teuchos::VERB_EXTREME);
368template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
371 for (
int i = 0; i < data.size(); ++i) {
372 const std::string fileName = ofname +
"_" + Teuchos::toString(data[i]) +
".m";
374 if (data[i] < H.GetNumLevels()) {
375 RCP<Level> L = H.GetLevel(data[i]);
377 if (L->IsAvailable(name)) {
378 RCP<T> M = L->template Get<RCP<T>>(name);
380 RCP<Matrix> A = L->template Get<RCP<Matrix>>(
"A");
381 RCP<const CrsGraph> AG = A->getCrsGraph();
382 WriteFieldContainer<T>(fileName, *M, *AG->getColMap());
389template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
392 size_t num_els = (size_t)fcont.extent(0);
393 size_t num_vecs = (size_t)fcont.extent(1);
396 Teuchos::RCP<const Map> rowMap = Xpetra::MapFactory<LO, GO, NO>::Build(colMap.lib(), Teuchos::OrdinalTraits<Xpetra::global_size_t>::invalid(), fcont.extent(0), colMap.getIndexBase(), colMap.getComm());
399 RCP<GOMultiVector> vec = Xpetra::MultiVectorFactory<GO, LO, GO, NO>::Build(rowMap, num_vecs);
401 for (
size_t j = 0; j < num_vecs; j++) {
402 Teuchos::ArrayRCP<GO> v = vec->getDataNonConst(j);
403 for (
size_t i = 0; i < num_els; i++)
404 v[i] = colMap.getGlobalElement(fcont(i, j));
407 Xpetra::IO<SC, LO, GO, NO>::WriteGOMV(fileName, *vec);
void WriteFieldContainer(const std::string &fileName, T &fcont, const Map &colMap) const
void ExportDataSetKeepFlagsNextLevel(Hierarchy &H, const Teuchos::Array< int > &data, const std::string &name) const
void ExportDataSetKeepFlags(Hierarchy &H, const Teuchos::Array< int > &data, const std::string &name) const
void WriteDataFC(Hierarchy &H, const Teuchos::Array< int > &data, const std::string &name, const std::string &ofname) const
void WriteData(Hierarchy &H, const Teuchos::Array< int > &data, const std::string &name) const
RCP< FactoryManagerBase > GetFactoryManager(int levelID) const
void WriteDataAggregates(Hierarchy &H, const Teuchos::Array< int > &data, const std::string &name) const
void AddFactoryManager(int startLevel, int numDesiredLevel, RCP< FactoryManagerBase > manager)
virtual void SetupHierarchy(Hierarchy &H) const
Setup Hierarchy object.
virtual RCP< Hierarchy > CreateHierarchy() const
Create an empty Hierarchy object.
void ExportDataSetKeepFlagsAll(Hierarchy &H, const std::string &name) const
HierarchyManager(int numDesiredLevel=MasterList::getDefault< int >("max levels"))
Constructor.
Teuchos::RCP< FactoryManagerBase > LvlMngr(int levelID, int lastLevelID) const
size_t getNumFactoryManagers() const
returns number of factory managers stored in levelManagers_ vector.
Static class that holds the complete list of valid MueLu parameters.
Namespace for MueLu classes and methods.