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 = Node::is_gpu;
107 
108  if (pL.get<bool>("rap: triple product") == false || isEpetra || isGPU) {
109  if (pL.get<bool>("rap: triple product") && isEpetra)
110  GetOStream(Warnings1) << "Switching from triple product to R x (A x P) since triple product has not been implemented for Epetra.\n";
111  if (pL.get<bool>("rap: triple product") && isGPU)
112  GetOStream(Warnings1) << "Switching from triple product to R x (A x P) since triple product has not been implemented for "
113  << Node::execution_space::name() << std::endl;
114 
115  // Reuse pattern if available (multiple solve)
116  RCP<ParameterList> APparams = rcp(new ParameterList);
117  if (pL.isSublist("matrixmatrix: kernel params"))
118  APparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
119 
120  // By default, we don't need global constants for A*P
121  APparams->set("compute global constants: temporaries", APparams->get("compute global constants: temporaries", false));
122  APparams->set("compute global constants", APparams->get("compute global constants", false));
123 
124  if (coarseLevel.IsAvailable("AP reuse data", this)) {
125  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous AP data" << std::endl;
126 
127  APparams = coarseLevel.Get<RCP<ParameterList> >("AP reuse data", this);
128 
129  if (APparams->isParameter("graph"))
130  AP = APparams->get<RCP<Matrix> >("graph");
131  }
132 
133  {
134  SubFactoryMonitor subM(*this, "MxM: A x P", coarseLevel);
135 
136  AP = MatrixMatrix::Multiply(*A, !doTranspose, *P, !doTranspose, AP, GetOStream(Statistics2),
137  doFillComplete, doOptimizeStorage, labelstr + std::string("MueLu::A*P-") + levelstr.str(), APparams);
138  }
139 
140  // Reuse coarse matrix memory if available (multiple solve)
141  RCP<ParameterList> RAPparams = rcp(new ParameterList);
142  if (pL.isSublist("matrixmatrix: kernel params"))
143  RAPparams = rcp(new ParameterList(pL.sublist("matrixmatrix: kernel params")));
144 
145  if (coarseLevel.IsAvailable("RAP reuse data", this)) {
146  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous RAP data" << std::endl;
147 
148  RAPparams = coarseLevel.Get<RCP<ParameterList> >("RAP reuse data", this);
149 
150  if (RAPparams->isParameter("graph"))
151  Ac = RAPparams->get<RCP<Matrix> >("graph");
152 
153  // Some eigenvalue may have been cached with the matrix in the previous run.
154  // As the matrix values will be updated, we need to reset the eigenvalue.
155  Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
156  }
157 
158  // We *always* need global constants for the RAP, but not for the temps
159  RAPparams->set("compute global constants: temporaries", RAPparams->get("compute global constants: temporaries", false));
160  RAPparams->set("compute global constants", true);
161 
162  // Allow optimization of storage.
163  // This is necessary for new faster Epetra MM kernels.
164  // Seems to work with matrix modifications to repair diagonal entries.
165 
166  if (pL.get<bool>("transpose: use implicit") == true) {
167  SubFactoryMonitor m2(*this, "MxM: P' x (AP) (implicit)", coarseLevel);
168 
169  Ac = MatrixMatrix::Multiply(*P, doTranspose, *AP, !doTranspose, Ac, GetOStream(Statistics2),
170  doFillComplete, doOptimizeStorage, labelstr + std::string("MueLu::R*(AP)-implicit-") + levelstr.str(), RAPparams);
171 
172  } else {
173  RCP<Matrix> R = Get<RCP<Matrix> >(coarseLevel, "R");
174 
175  SubFactoryMonitor m2(*this, "MxM: R x (AP) (explicit)", coarseLevel);
176 
177  Ac = MatrixMatrix::Multiply(*R, !doTranspose, *AP, !doTranspose, Ac, GetOStream(Statistics2),
178  doFillComplete, doOptimizeStorage, labelstr + std::string("MueLu::R*(AP)-explicit-") + levelstr.str(), RAPparams);
179  }
180 
181  Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >("rap: relative diagonal floor")();
182  if (relativeFloor.size() > 0) {
184  }
185 
186  bool repairZeroDiagonals = pL.get<bool>("RepairMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
187  bool checkAc = pL.get<bool>("CheckMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
188  ;
189  if (checkAc || repairZeroDiagonals) {
190  using magnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
191  magnitudeType threshold;
192  if (pL.isType<magnitudeType>("rap: fix zero diagonals threshold"))
193  threshold = pL.get<magnitudeType>("rap: fix zero diagonals threshold");
194  else
195  threshold = Teuchos::as<magnitudeType>(pL.get<double>("rap: fix zero diagonals threshold"));
196  Scalar replacement = Teuchos::as<Scalar>(pL.get<double>("rap: fix zero diagonals replacement"));
197  Xpetra::MatrixUtils<SC, LO, GO, NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(Warnings1), threshold, replacement);
198  }
199 
200  if (IsPrint(Statistics2)) {
201  RCP<ParameterList> params = rcp(new ParameterList());
202  ;
203  params->set("printLoadBalancingInfo", true);
204  params->set("printCommInfo", true);
205  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
206  }
207 
208  if (!Ac.is_null()) {
209  std::ostringstream oss;
210  oss << "A_" << coarseLevel.GetLevelID();
211  Ac->setObjectLabel(oss.str());
212  }
213  Set(coarseLevel, "A", Ac);
214 
215  if (!isGPU) {
216  APparams->set("graph", AP);
217  Set(coarseLevel, "AP reuse data", APparams);
218  }
219  if (!isGPU) {
220  RAPparams->set("graph", Ac);
221  Set(coarseLevel, "RAP reuse data", RAPparams);
222  }
223  } else {
224  RCP<ParameterList> RAPparams = rcp(new ParameterList);
225  if (pL.isSublist("matrixmatrix: kernel params"))
226  RAPparams->sublist("matrixmatrix: kernel params") = pL.sublist("matrixmatrix: kernel params");
227 
228  if (coarseLevel.IsAvailable("RAP reuse data", this)) {
229  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous RAP data" << std::endl;
230 
231  RAPparams = coarseLevel.Get<RCP<ParameterList> >("RAP reuse data", this);
232 
233  if (RAPparams->isParameter("graph"))
234  Ac = RAPparams->get<RCP<Matrix> >("graph");
235 
236  // Some eigenvalue may have been cached with the matrix in the previous run.
237  // As the matrix values will be updated, we need to reset the eigenvalue.
238  Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
239  }
240 
241  // We *always* need global constants for the RAP, but not for the temps
242  RAPparams->set("compute global constants: temporaries", RAPparams->get("compute global constants: temporaries", false));
243  RAPparams->set("compute global constants", true);
244 
245  if (pL.get<bool>("transpose: use implicit") == true) {
246  Ac = MatrixFactory::Build(P->getDomainMap(), Teuchos::as<LO>(0));
247 
248  SubFactoryMonitor m2(*this, "MxMxM: R x A x P (implicit)", coarseLevel);
249 
251  MultiplyRAP(*P, doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
252  doOptimizeStorage, labelstr + std::string("MueLu::R*A*P-implicit-") + levelstr.str(),
253  RAPparams);
254  } else {
255  RCP<Matrix> R = Get<RCP<Matrix> >(coarseLevel, "R");
256  Ac = MatrixFactory::Build(R->getRowMap(), Teuchos::as<LO>(0));
257 
258  SubFactoryMonitor m2(*this, "MxMxM: R x A x P (explicit)", coarseLevel);
259 
261  MultiplyRAP(*R, !doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
262  doOptimizeStorage, labelstr + std::string("MueLu::R*A*P-explicit-") + levelstr.str(),
263  RAPparams);
264  }
265 
266  Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >("rap: relative diagonal floor")();
267  if (relativeFloor.size() > 0) {
269  }
270 
271  bool repairZeroDiagonals = pL.get<bool>("RepairMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
272  bool checkAc = pL.get<bool>("CheckMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
273  ;
274  if (checkAc || repairZeroDiagonals) {
275  using magnitudeType = typename Teuchos::ScalarTraits<Scalar>::magnitudeType;
276  magnitudeType threshold;
277  if (pL.isType<magnitudeType>("rap: fix zero diagonals threshold"))
278  threshold = pL.get<magnitudeType>("rap: fix zero diagonals threshold");
279  else
280  threshold = Teuchos::as<magnitudeType>(pL.get<double>("rap: fix zero diagonals threshold"));
281  Scalar replacement = Teuchos::as<Scalar>(pL.get<double>("rap: fix zero diagonals replacement"));
282  Xpetra::MatrixUtils<SC, LO, GO, NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(Warnings1), threshold, replacement);
283  }
284 
285  if (IsPrint(Statistics2)) {
286  RCP<ParameterList> params = rcp(new ParameterList());
287  ;
288  params->set("printLoadBalancingInfo", true);
289  params->set("printCommInfo", true);
290  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
291  }
292 
293  if (!Ac.is_null()) {
294  std::ostringstream oss;
295  oss << "A_" << coarseLevel.GetLevelID();
296  Ac->setObjectLabel(oss.str());
297  }
298  Set(coarseLevel, "A", Ac);
299 
300  if (!isGPU) {
301  RAPparams->set("graph", Ac);
302  Set(coarseLevel, "RAP reuse data", RAPparams);
303  }
304  }
305  }
306 
307 #ifdef HAVE_MUELU_DEBUG
308  MatrixUtils::checkLocalRowMapMatchesColMap(*Ac);
309 #endif // HAVE_MUELU_DEBUG
310 
311  if (transferFacts_.begin() != transferFacts_.end()) {
312  SubFactoryMonitor m(*this, "Projections", coarseLevel);
313 
314  // call Build of all user-given transfer factories
315  for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it) {
316  RCP<const FactoryBase> fac = *it;
317  GetOStream(Runtime0) << "RAPFactory: call transfer factory: " << fac->description() << std::endl;
318  fac->CallBuild(coarseLevel);
319  // Coordinates transfer is marginally different from all other operations
320  // because it is *optional*, and not required. For instance, we may need
321  // coordinates only on level 4 if we start repartitioning from that level,
322  // but we don't need them on level 1,2,3. As our current Hierarchy setup
323  // assumes propagation of dependencies only through three levels, this
324  // means that we need to rely on other methods to propagate optional data.
325  //
326  // The method currently used is through RAP transfer factories, which are
327  // simply factories which are called at the end of RAP with a single goal:
328  // transfer some fine data to coarser level. Because these factories are
329  // kind of outside of the mainline factories, they behave different. In
330  // particular, we call their Build method explicitly, rather than through
331  // Get calls. This difference is significant, as the Get call is smart
332  // enough to know when to release all factory dependencies, and Build is
333  // dumb. This led to the following CoordinatesTransferFactory sequence:
334  // 1. Request level 0
335  // 2. Request level 1
336  // 3. Request level 0
337  // 4. Release level 0
338  // 5. Release level 1
339  //
340  // The problem is missing "6. Release level 0". Because it was missing,
341  // we had outstanding request on "Coordinates", "Aggregates" and
342  // "CoarseMap" on level 0.
343  //
344  // This was fixed by explicitly calling Release on transfer factories in
345  // RAPFactory. I am still unsure how exactly it works, but now we have
346  // clear data requests for all levels.
347  coarseLevel.Release(*fac);
348  }
349  }
350 }
351 
352 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354  // check if it's a TwoLevelFactoryBase based transfer factory
355  TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::rcp_dynamic_cast<const TwoLevelFactoryBase>(factory) == Teuchos::null, Exceptions::BadCast,
356  "MueLu::RAPFactory::AddTransferFactory: Transfer factory is not derived from TwoLevelFactoryBase. "
357  "This is very strange. (Note: you can remove this exception if there's a good reason for)");
358  TEUCHOS_TEST_FOR_EXCEPTION(hasDeclaredInput_, Exceptions::RuntimeError, "MueLu::RAPFactory::AddTransferFactory: Factory is being added after we have already declared input");
359  transferFacts_.push_back(factory);
360 }
361 
362 } // namespace MueLu
363 
364 #define MUELU_RAPFACTORY_SHORT
365 #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)
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