MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_RAPFactory_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef MUELU_RAPFACTORY_DEF_HPP
11 #define MUELU_RAPFACTORY_DEF_HPP
12 
13 #include <sstream>
14 
15 #include <Xpetra_Matrix.hpp>
16 #include <Xpetra_MatrixFactory.hpp>
17 #include <Xpetra_MatrixMatrix.hpp>
18 #include <Xpetra_MatrixUtils.hpp>
20 
22 
23 #include "MueLu_MasterList.hpp"
24 #include "MueLu_Monitor.hpp"
25 #include "MueLu_PerfUtils.hpp"
26 
27 namespace MueLu {
28 
29 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31  : hasDeclaredInput_(false) {}
32 
33 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
35 
36 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
38  RCP<ParameterList> validParamList = rcp(new ParameterList());
39 
40 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
41  SET_VALID_ENTRY("transpose: use implicit");
42  SET_VALID_ENTRY("rap: triple product");
43  SET_VALID_ENTRY("rap: fix zero diagonals");
44  SET_VALID_ENTRY("rap: fix zero diagonals threshold");
45  SET_VALID_ENTRY("rap: fix zero diagonals replacement");
46  SET_VALID_ENTRY("rap: relative diagonal floor");
47 #undef SET_VALID_ENTRY
48  validParamList->set<RCP<const FactoryBase> >("A", null, "Generating factory of the matrix A used during the prolongator smoothing process");
49  validParamList->set<RCP<const FactoryBase> >("P", null, "Prolongator factory");
50  validParamList->set<RCP<const FactoryBase> >("R", null, "Restrictor factory");
51 
52  validParamList->set<bool>("CheckMainDiagonal", false, "Check main diagonal for zeros");
53  validParamList->set<bool>("RepairMainDiagonal", false, "Repair zeros on main diagonal");
54 
55  // Make sure we don't recursively validate options for the matrixmatrix kernels
56  ParameterList norecurse;
57  norecurse.disableRecursiveValidation();
58  validParamList->set<ParameterList>("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
59 
60  return validParamList;
61 }
62 
63 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
65  const Teuchos::ParameterList& pL = GetParameterList();
66  if (pL.get<bool>("transpose: use implicit") == false)
67  Input(coarseLevel, "R");
68 
69  Input(fineLevel, "A");
70  Input(coarseLevel, "P");
71 
72  // call DeclareInput of all user-given transfer factories
73  for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it)
74  (*it)->CallDeclareInput(coarseLevel);
75 
76  hasDeclaredInput_ = true;
77 }
78 
79 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
81  const bool doTranspose = true;
82  const bool doFillComplete = true;
83  const bool doOptimizeStorage = true;
84  RCP<Matrix> Ac;
85  {
86  FactoryMonitor m(*this, "Computing Ac", coarseLevel);
87  std::ostringstream levelstr;
88  levelstr << coarseLevel.GetLevelID();
89  std::string labelstr = FormattingHelper::getColonLabel(coarseLevel.getObjectLabel());
90 
91  TEUCHOS_TEST_FOR_EXCEPTION(hasDeclaredInput_ == false, Exceptions::RuntimeError,
92  "MueLu::RAPFactory::Build(): CallDeclareInput has not been called before Build!");
93 
94  const Teuchos::ParameterList& pL = GetParameterList();
95  RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel, "A");
96  RCP<Matrix> P = Get<RCP<Matrix> >(coarseLevel, "P"), AP;
97  // We don't have a valid P (e.g., # global aggregates = 0) so we bail.
98  // This level will ultimately be removed in MueLu_Hierarchy_defs.h via a resize()
99  if (P == Teuchos::null) {
100  Ac = Teuchos::null;
101  Set(coarseLevel, "A", Ac);
102  return;
103  }
104 
105  bool isEpetra = A->getRowMap()->lib() == Xpetra::UseEpetra;
106  bool isGPU =
107 #ifdef KOKKOS_ENABLE_CUDA
108  (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosCudaWrapperNode).name()) ||
109 #endif
110 #ifdef KOKKOS_ENABLE_HIP
111  (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosHIPWrapperNode).name()) ||
112 #endif
113 #ifdef KOKKOS_ENABLE_SYCL
114  (typeid(Node).name() == typeid(Tpetra::KokkosCompat::KokkosSYCLWrapperNode).name()) ||
115 #endif
116  false;
117 
118  if (pL.get<bool>("rap: triple product") == false || isEpetra || isGPU) {
119  if (pL.get<bool>("rap: triple product") && isEpetra)
120  GetOStream(Warnings1) << "Switching from triple product to R x (A x P) since triple product has not been implemented for Epetra.\n";
121 #if defined(KOKKOS_ENABLE_CUDA) || defined(KOKKOS_ENABLE_HIP) || defined(KOKKOS_ENABLE_SYCL)
122  if (pL.get<bool>("rap: triple product") && isGPU)
123  GetOStream(Warnings1) << "Switching from triple product to R x (A x P) since triple product has not been implemented for "
124  << Node::execution_space::name() << std::endl;
125 #endif
126 
127  // Reuse pattern if available (multiple solve)
128  RCP<ParameterList> APparams = rcp(new ParameterList);
129  if (pL.isSublist("matrixmatrix: kernel params"))
130  APparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
131 
132  // By default, we don't need global constants for A*P
133  APparams->set("compute global constants: temporaries", APparams->get("compute global constants: temporaries", false));
134  APparams->set("compute global constants", APparams->get("compute global constants", false));
135 
136  if (coarseLevel.IsAvailable("AP reuse data", this)) {
137  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous AP data" << std::endl;
138 
139  APparams = coarseLevel.Get<RCP<ParameterList> >("AP reuse data", this);
140 
141  if (APparams->isParameter("graph"))
142  AP = APparams->get<RCP<Matrix> >("graph");
143  }
144 
145  {
146  SubFactoryMonitor subM(*this, "MxM: A x P", coarseLevel);
147 
148  AP = MatrixMatrix::Multiply(*A, !doTranspose, *P, !doTranspose, AP, GetOStream(Statistics2),
149  doFillComplete, doOptimizeStorage, labelstr + std::string("MueLu::A*P-") + levelstr.str(), APparams);
150  }
151 
152  // Reuse coarse matrix memory if available (multiple solve)
153  RCP<ParameterList> RAPparams = rcp(new ParameterList);
154  if (pL.isSublist("matrixmatrix: kernel params"))
155  RAPparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
156 
157  if (coarseLevel.IsAvailable("RAP reuse data", this)) {
158  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous RAP data" << std::endl;
159 
160  RAPparams = coarseLevel.Get<RCP<ParameterList> >("RAP reuse data", this);
161 
162  if (RAPparams->isParameter("graph"))
163  Ac = RAPparams->get<RCP<Matrix> >("graph");
164 
165  // Some eigenvalue may have been cached with the matrix in the previous run.
166  // As the matrix values will be updated, we need to reset the eigenvalue.
167  Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
168  }
169 
170  // We *always* need global constants for the RAP, but not for the temps
171  RAPparams->set("compute global constants: temporaries", RAPparams->get("compute global constants: temporaries", false));
172  RAPparams->set("compute global constants", true);
173 
174  // Allow optimization of storage.
175  // This is necessary for new faster Epetra MM kernels.
176  // Seems to work with matrix modifications to repair diagonal entries.
177 
178  if (pL.get<bool>("transpose: use implicit") == true) {
179  SubFactoryMonitor m2(*this, "MxM: P' x (AP) (implicit)", coarseLevel);
180 
181  Ac = MatrixMatrix::Multiply(*P, doTranspose, *AP, !doTranspose, Ac, GetOStream(Statistics2),
182  doFillComplete, doOptimizeStorage, labelstr + std::string("MueLu::R*(AP)-implicit-") + levelstr.str(), RAPparams);
183 
184  } else {
185  RCP<Matrix> R = Get<RCP<Matrix> >(coarseLevel, "R");
186 
187  SubFactoryMonitor m2(*this, "MxM: R x (AP) (explicit)", coarseLevel);
188 
189  Ac = MatrixMatrix::Multiply(*R, !doTranspose, *AP, !doTranspose, Ac, GetOStream(Statistics2),
190  doFillComplete, doOptimizeStorage, labelstr + std::string("MueLu::R*(AP)-explicit-") + levelstr.str(), RAPparams);
191  }
192 
193  Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >("rap: relative diagonal floor")();
194  if (relativeFloor.size() > 0) {
196  }
197 
198  bool repairZeroDiagonals = pL.get<bool>("RepairMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
199  bool checkAc = pL.get<bool>("CheckMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
200  ;
201  if (checkAc || repairZeroDiagonals) {
202  using magnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
203  magnitudeType threshold;
204  if (pL.isType<magnitudeType>("rap: fix zero diagonals threshold"))
205  threshold = pL.get<magnitudeType>("rap: fix zero diagonals threshold");
206  else
207  threshold = Teuchos::as<magnitudeType>(pL.get<double>("rap: fix zero diagonals threshold"));
208  Scalar replacement = Teuchos::as<Scalar>(pL.get<double>("rap: fix zero diagonals replacement"));
209  Xpetra::MatrixUtils<SC, LO, GO, NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(Warnings1), threshold, replacement);
210  }
211 
212  if (IsPrint(Statistics2)) {
213  RCP<ParameterList> params = rcp(new ParameterList());
214  ;
215  params->set("printLoadBalancingInfo", true);
216  params->set("printCommInfo", true);
217  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
218  }
219 
220  if (!Ac.is_null()) {
221  std::ostringstream oss;
222  oss << "A_" << coarseLevel.GetLevelID();
223  Ac->setObjectLabel(oss.str());
224  }
225  Set(coarseLevel, "A", Ac);
226 
227  if (!isGPU) {
228  APparams->set("graph", AP);
229  Set(coarseLevel, "AP reuse data", APparams);
230  }
231  if (!isGPU) {
232  RAPparams->set("graph", Ac);
233  Set(coarseLevel, "RAP reuse data", RAPparams);
234  }
235  } else {
236  RCP<ParameterList> RAPparams = rcp(new ParameterList);
237  if (pL.isSublist("matrixmatrix: kernel params"))
238  RAPparams->sublist("matrixmatrix: kernel params") = pL.sublist("matrixmatrix: kernel params");
239 
240  if (coarseLevel.IsAvailable("RAP reuse data", this)) {
241  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous RAP data" << std::endl;
242 
243  RAPparams = coarseLevel.Get<RCP<ParameterList> >("RAP reuse data", this);
244 
245  if (RAPparams->isParameter("graph"))
246  Ac = RAPparams->get<RCP<Matrix> >("graph");
247 
248  // Some eigenvalue may have been cached with the matrix in the previous run.
249  // As the matrix values will be updated, we need to reset the eigenvalue.
250  Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
251  }
252 
253  // We *always* need global constants for the RAP, but not for the temps
254  RAPparams->set("compute global constants: temporaries", RAPparams->get("compute global constants: temporaries", false));
255  RAPparams->set("compute global constants", true);
256 
257  if (pL.get<bool>("transpose: use implicit") == true) {
258  Ac = MatrixFactory::Build(P->getDomainMap(), Teuchos::as<LO>(0));
259 
260  SubFactoryMonitor m2(*this, "MxMxM: R x A x P (implicit)", coarseLevel);
261 
263  MultiplyRAP(*P, doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
264  doOptimizeStorage, labelstr + std::string("MueLu::R*A*P-implicit-") + levelstr.str(),
265  RAPparams);
266  } else {
267  RCP<Matrix> R = Get<RCP<Matrix> >(coarseLevel, "R");
268  Ac = MatrixFactory::Build(R->getRowMap(), Teuchos::as<LO>(0));
269 
270  SubFactoryMonitor m2(*this, "MxMxM: R x A x P (explicit)", coarseLevel);
271 
273  MultiplyRAP(*R, !doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
274  doOptimizeStorage, labelstr + std::string("MueLu::R*A*P-explicit-") + levelstr.str(),
275  RAPparams);
276  }
277 
278  Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >("rap: relative diagonal floor")();
279  if (relativeFloor.size() > 0) {
281  }
282 
283  bool repairZeroDiagonals = pL.get<bool>("RepairMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
284  bool checkAc = pL.get<bool>("CheckMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
285  ;
286  if (checkAc || repairZeroDiagonals) {
287  using magnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
288  magnitudeType threshold;
289  if (pL.isType<magnitudeType>("rap: fix zero diagonals threshold"))
290  threshold = pL.get<magnitudeType>("rap: fix zero diagonals threshold");
291  else
292  threshold = Teuchos::as<magnitudeType>(pL.get<double>("rap: fix zero diagonals threshold"));
293  Scalar replacement = Teuchos::as<Scalar>(pL.get<double>("rap: fix zero diagonals replacement"));
294  Xpetra::MatrixUtils<SC, LO, GO, NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(Warnings1), threshold, replacement);
295  }
296 
297  if (IsPrint(Statistics2)) {
298  RCP<ParameterList> params = rcp(new ParameterList());
299  ;
300  params->set("printLoadBalancingInfo", true);
301  params->set("printCommInfo", true);
302  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
303  }
304 
305  if (!Ac.is_null()) {
306  std::ostringstream oss;
307  oss << "A_" << coarseLevel.GetLevelID();
308  Ac->setObjectLabel(oss.str());
309  }
310  Set(coarseLevel, "A", Ac);
311 
312  if (!isGPU) {
313  RAPparams->set("graph", Ac);
314  Set(coarseLevel, "RAP reuse data", RAPparams);
315  }
316  }
317  }
318 
319 #ifdef HAVE_MUELU_DEBUG
320  MatrixUtils::checkLocalRowMapMatchesColMap(*Ac);
321 #endif // HAVE_MUELU_DEBUG
322 
323  if (transferFacts_.begin() != transferFacts_.end()) {
324  SubFactoryMonitor m(*this, "Projections", coarseLevel);
325 
326  // call Build of all user-given transfer factories
327  for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it) {
328  RCP<const FactoryBase> fac = *it;
329  GetOStream(Runtime0) << "RAPFactory: call transfer factory: " << fac->description() << std::endl;
330  fac->CallBuild(coarseLevel);
331  // Coordinates transfer is marginally different from all other operations
332  // because it is *optional*, and not required. For instance, we may need
333  // coordinates only on level 4 if we start repartitioning from that level,
334  // but we don't need them on level 1,2,3. As our current Hierarchy setup
335  // assumes propagation of dependencies only through three levels, this
336  // means that we need to rely on other methods to propagate optional data.
337  //
338  // The method currently used is through RAP transfer factories, which are
339  // simply factories which are called at the end of RAP with a single goal:
340  // transfer some fine data to coarser level. Because these factories are
341  // kind of outside of the mainline factories, they behave different. In
342  // particular, we call their Build method explicitly, rather than through
343  // Get calls. This difference is significant, as the Get call is smart
344  // enough to know when to release all factory dependencies, and Build is
345  // dumb. This led to the following CoordinatesTransferFactory sequence:
346  // 1. Request level 0
347  // 2. Request level 1
348  // 3. Request level 0
349  // 4. Release level 0
350  // 5. Release level 1
351  //
352  // The problem is missing "6. Release level 0". Because it was missing,
353  // we had outstanding request on "Coordinates", "Aggregates" and
354  // "CoarseMap" on level 0.
355  //
356  // This was fixed by explicitly calling Release on transfer factories in
357  // RAPFactory. I am still unsure how exactly it works, but now we have
358  // clear data requests for all levels.
359  coarseLevel.Release(*fac);
360  }
361  }
362 }
363 
364 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
366  // check if it's a TwoLevelFactoryBase based transfer factory
367  TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::rcp_dynamic_cast<const TwoLevelFactoryBase>(factory) == Teuchos::null, Exceptions::BadCast,
368  "MueLu::RAPFactory::AddTransferFactory: Transfer factory is not derived from TwoLevelFactoryBase. "
369  "This is very strange. (Note: you can remove this exception if there's a good reason for)");
370  TEUCHOS_TEST_FOR_EXCEPTION(hasDeclaredInput_, Exceptions::RuntimeError, "MueLu::RAPFactory::AddTransferFactory: Factory is being added after we have already declared input");
371  transferFacts_.push_back(factory);
372 }
373 
374 } // namespace MueLu
375 
376 #define MUELU_RAPFACTORY_SHORT
377 #endif // MUELU_RAPFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
Exception indicating invalid cast attempted.
virtual std::string getObjectLabel() const
static void RelativeDiagonalBoost(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &A, const Teuchos::ArrayView< const double > &relativeThreshold, Teuchos::FancyOStream &fos)
virtual void CallBuild(Level &requestedLevel) const =0
ParameterList & disableRecursiveValidation()
void Release(const FactoryBase &factory)
Decrement the storage counter for all the inputs of a factory.
T & get(const std::string &name, T def_value)
static void MultiplyRAP(const Matrix &R, bool transposeR, const Matrix &A, bool transposeA, const Matrix &P, bool transposeP, Matrix &Ac, bool call_FillComplete_on_result=true, bool doOptimizeStorage=true, const std::string &label=std::string(), const RCP< ParameterList > &params=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
size_type size() const
One-liner description of what is happening.
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
MueLu::DefaultNode Node
static void CheckRepairMainDiagonal(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &Ac, bool const &repairZeroDiagonals, Teuchos::FancyOStream &fos, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType threshold=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< Scalar >::magnitudeType >::zero())
static std::string getColonLabel(const std::string &label)
Helper function for object label.
Print even more statistics.
Additional warnings.
bool isParameter(const std::string &name) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
bool isSublist(const std::string &name) const
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
virtual ~RAPFactory()
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
bool isType(const std::string &name) const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:51
Exception throws to report errors in the internal logical of the program.
virtual std::string description() const
Return a simple one-line description of this object.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
bool is_null() const