MueLu  Version of the Day
MueLu_RAPFactory_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_RAPFACTORY_DEF_HPP
47 #define MUELU_RAPFACTORY_DEF_HPP
48 
49 
50 #include <sstream>
51 
52 #include <Xpetra_Matrix.hpp>
53 #include <Xpetra_MatrixFactory.hpp>
54 #include <Xpetra_MatrixMatrix.hpp>
55 #include <Xpetra_MatrixUtils.hpp>
56 #include <Xpetra_TripleMatrixMultiply.hpp>
57 #include <Xpetra_Vector.hpp>
58 #include <Xpetra_VectorFactory.hpp>
59 
61 
62 #include "MueLu_MasterList.hpp"
63 #include "MueLu_Monitor.hpp"
64 #include "MueLu_PerfUtils.hpp"
66 //#include "MueLu_Utilities.hpp"
67 
68 namespace MueLu {
69 
70  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
72  : hasDeclaredInput_(false) { }
73 
74  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
76  RCP<ParameterList> validParamList = rcp(new ParameterList());
77 
78 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
79  SET_VALID_ENTRY("transpose: use implicit");
80  SET_VALID_ENTRY("rap: triple product");
81  SET_VALID_ENTRY("rap: fix zero diagonals");
82  SET_VALID_ENTRY("rap: relative diagonal floor");
83 #undef SET_VALID_ENTRY
84  validParamList->set< RCP<const FactoryBase> >("A", null, "Generating factory of the matrix A used during the prolongator smoothing process");
85  validParamList->set< RCP<const FactoryBase> >("P", null, "Prolongator factory");
86  validParamList->set< RCP<const FactoryBase> >("R", null, "Restrictor factory");
87 
88  validParamList->set< bool > ("CheckMainDiagonal", false, "Check main diagonal for zeros");
89  validParamList->set< bool > ("RepairMainDiagonal", false, "Repair zeros on main diagonal");
90 
91  // Make sure we don't recursively validate options for the matrixmatrix kernels
92  ParameterList norecurse;
93  norecurse.disableRecursiveValidation();
94  validParamList->set<ParameterList> ("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
95 
96  return validParamList;
97  }
98 
99  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
101  const Teuchos::ParameterList& pL = GetParameterList();
102  if (pL.get<bool>("transpose: use implicit") == false)
103  Input(coarseLevel, "R");
104 
105  Input(fineLevel, "A");
106  Input(coarseLevel, "P");
107 
108  // call DeclareInput of all user-given transfer factories
109  for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it)
110  (*it)->CallDeclareInput(coarseLevel);
111 
112  hasDeclaredInput_ = true;
113  }
114 
115  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
117  const bool doTranspose = true;
118  const bool doFillComplete = true;
119  const bool doOptimizeStorage = true;
120  RCP<Matrix> Ac;
121  {
122  FactoryMonitor m(*this, "Computing Ac", coarseLevel);
123  std::ostringstream levelstr;
124  levelstr << coarseLevel.GetLevelID();
125  std::string labelstr = FormattingHelper::getColonLabel(coarseLevel.getObjectLabel());
126 
127  TEUCHOS_TEST_FOR_EXCEPTION(hasDeclaredInput_ == false, Exceptions::RuntimeError,
128  "MueLu::RAPFactory::Build(): CallDeclareInput has not been called before Build!");
129 
130  const Teuchos::ParameterList& pL = GetParameterList();
131  RCP<Matrix> A = Get< RCP<Matrix> >(fineLevel, "A");
132  RCP<Matrix> P = Get< RCP<Matrix> >(coarseLevel, "P"), AP;
133 
134  bool isEpetra = A->getRowMap()->lib() == Xpetra::UseEpetra;
135 #ifdef KOKKOS_ENABLE_CUDA
136  bool isCuda = typeid(Node).name() == typeid(Kokkos::Compat::KokkosCudaWrapperNode).name();
137 #else
138  bool isCuda = false;
139 #endif
140 
141  if (pL.get<bool>("rap: triple product") == false || isEpetra || isCuda) {
142  if (pL.get<bool>("rap: triple product") && isEpetra)
143  GetOStream(Warnings1) << "Switching from triple product to R x (A x P) since triple product has not been implemented for Epetra.\n";
144 #ifdef KOKKOS_ENABLE_CUDA
145  if (pL.get<bool>("rap: triple product") && isCuda)
146  GetOStream(Warnings1) << "Switching from triple product to R x (A x P) since triple product has not been implemented for Cuda.\n";
147 #endif
148 
149  // Reuse pattern if available (multiple solve)
150  RCP<ParameterList> APparams = rcp(new ParameterList);
151  if(pL.isSublist("matrixmatrix: kernel params"))
152  APparams->sublist("matrixmatrix: kernel params") = pL.sublist("matrixmatrix: kernel params");
153 
154  // By default, we don't need global constants for A*P
155  APparams->set("compute global constants: temporaries",APparams->get("compute global constants: temporaries",false));
156  APparams->set("compute global constants",APparams->get("compute global constants",false));
157 
158  if (coarseLevel.IsAvailable("AP reuse data", this)) {
159  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous AP data" << std::endl;
160 
161  APparams = coarseLevel.Get< RCP<ParameterList> >("AP reuse data", this);
162 
163  if (APparams->isParameter("graph"))
164  AP = APparams->get< RCP<Matrix> >("graph");
165  }
166 
167  {
168  SubFactoryMonitor subM(*this, "MxM: A x P", coarseLevel);
169 
170  AP = MatrixMatrix::Multiply(*A, !doTranspose, *P, !doTranspose, AP, GetOStream(Statistics2),
171  doFillComplete, doOptimizeStorage, labelstr+std::string("MueLu::A*P-")+levelstr.str(), APparams);
172  }
173 
174  // Reuse coarse matrix memory if available (multiple solve)
175  RCP<ParameterList> RAPparams = rcp(new ParameterList);
176  if(pL.isSublist("matrixmatrix: kernel params"))
177  RAPparams->sublist("matrixmatrix: kernel params") = pL.sublist("matrixmatrix: kernel params");
178 
179  if (coarseLevel.IsAvailable("RAP reuse data", this)) {
180  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous RAP data" << std::endl;
181 
182  RAPparams = coarseLevel.Get< RCP<ParameterList> >("RAP reuse data", this);
183 
184  if (RAPparams->isParameter("graph"))
185  Ac = RAPparams->get< RCP<Matrix> >("graph");
186 
187  // Some eigenvalue may have been cached with the matrix in the previous run.
188  // As the matrix values will be updated, we need to reset the eigenvalue.
189  Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
190  }
191 
192  // We *always* need global constants for the RAP, but not for the temps
193  RAPparams->set("compute global constants: temporaries",RAPparams->get("compute global constants: temporaries",false));
194  RAPparams->set("compute global constants",true);
195 
196  // Allow optimization of storage.
197  // This is necessary for new faster Epetra MM kernels.
198  // Seems to work with matrix modifications to repair diagonal entries.
199 
200  if (pL.get<bool>("transpose: use implicit") == true) {
201  SubFactoryMonitor m2(*this, "MxM: P' x (AP) (implicit)", coarseLevel);
202 
203  Ac = MatrixMatrix::Multiply(*P, doTranspose, *AP, !doTranspose, Ac, GetOStream(Statistics2),
204  doFillComplete, doOptimizeStorage, labelstr+std::string("MueLu::R*(AP)-implicit-")+levelstr.str(), RAPparams);
205 
206  } else {
207  RCP<Matrix> R = Get< RCP<Matrix> >(coarseLevel, "R");
208 
209  SubFactoryMonitor m2(*this, "MxM: R x (AP) (explicit)", coarseLevel);
210 
211  Ac = MatrixMatrix::Multiply(*R, !doTranspose, *AP, !doTranspose, Ac, GetOStream(Statistics2),
212  doFillComplete, doOptimizeStorage, labelstr+std::string("MueLu::R*(AP)-explicit-")+levelstr.str(), RAPparams);
213  }
214 
215  Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >("rap: relative diagonal floor")();
216  if(relativeFloor.size() > 0) {
217  Xpetra::MatrixUtils<SC,LO,GO,NO>::RelativeDiagonalBoost(Ac, relativeFloor,GetOStream(Statistics2));
218  }
219 
220  bool repairZeroDiagonals = pL.get<bool>("RepairMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
221  bool checkAc = pL.get<bool>("CheckMainDiagonal")|| pL.get<bool>("rap: fix zero diagonals"); ;
222  if (checkAc || repairZeroDiagonals)
223  Xpetra::MatrixUtils<SC,LO,GO,NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(Warnings1));
224 
225  if (IsPrint(Statistics2)) {
226  RCP<ParameterList> params = rcp(new ParameterList());;
227  params->set("printLoadBalancingInfo", true);
228  params->set("printCommInfo", true);
229  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
230  }
231 
232  if(!Ac.is_null()) {std::ostringstream oss; oss << "A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
233  Set(coarseLevel, "A", Ac);
234 
235  APparams->set("graph", AP);
236  Set(coarseLevel, "AP reuse data", APparams);
237  RAPparams->set("graph", Ac);
238  Set(coarseLevel, "RAP reuse data", RAPparams);
239  } else {
240  RCP<ParameterList> RAPparams = rcp(new ParameterList);
241  if(pL.isSublist("matrixmatrix: kernel params"))
242  RAPparams->sublist("matrixmatrix: kernel params") = pL.sublist("matrixmatrix: kernel params");
243 
244  if (coarseLevel.IsAvailable("RAP reuse data", this)) {
245  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous RAP data" << std::endl;
246 
247  RAPparams = coarseLevel.Get< RCP<ParameterList> >("RAP reuse data", this);
248 
249  if (RAPparams->isParameter("graph"))
250  Ac = RAPparams->get< RCP<Matrix> >("graph");
251 
252  // Some eigenvalue may have been cached with the matrix in the previous run.
253  // As the matrix values will be updated, we need to reset the eigenvalue.
254  Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
255  }
256 
257  // We *always* need global constants for the RAP, but not for the temps
258  RAPparams->set("compute global constants: temporaries",RAPparams->get("compute global constants: temporaries",false));
259  RAPparams->set("compute global constants",true);
260 
261  if (pL.get<bool>("transpose: use implicit") == true) {
262 
263  Ac = MatrixFactory::Build(P->getDomainMap(), Teuchos::as<LO>(0));
264 
265  SubFactoryMonitor m2(*this, "MxMxM: R x A x P (implicit)", coarseLevel);
266 
267  Xpetra::TripleMatrixMultiply<SC,LO,GO,NO>::
268  MultiplyRAP(*P, doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
269  doOptimizeStorage, labelstr+std::string("MueLu::R*A*P-implicit-")+levelstr.str(),
270  RAPparams);
271 
272  } else {
273  RCP<Matrix> R = Get< RCP<Matrix> >(coarseLevel, "R");
274  Ac = MatrixFactory::Build(R->getRowMap(), Teuchos::as<LO>(0));
275 
276  SubFactoryMonitor m2(*this, "MxMxM: R x A x P (explicit)", coarseLevel);
277 
278  Xpetra::TripleMatrixMultiply<SC,LO,GO,NO>::
279  MultiplyRAP(*R, !doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
280  doOptimizeStorage, labelstr+std::string("MueLu::R*A*P-explicit-")+levelstr.str(),
281  RAPparams);
282  }
283 
284  Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >("rap: relative diagonal floor")();
285  if(relativeFloor.size() > 0) {
286  Xpetra::MatrixUtils<SC,LO,GO,NO>::RelativeDiagonalBoost(Ac, relativeFloor,GetOStream(Statistics2));
287  }
288 
289  bool repairZeroDiagonals = pL.get<bool>("RepairMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
290  bool checkAc = pL.get<bool>("CheckMainDiagonal")|| pL.get<bool>("rap: fix zero diagonals"); ;
291  if (checkAc || repairZeroDiagonals)
292  Xpetra::MatrixUtils<SC,LO,GO,NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(Warnings1));
293 
294 
295 
296  if (IsPrint(Statistics2)) {
297  RCP<ParameterList> params = rcp(new ParameterList());;
298  params->set("printLoadBalancingInfo", true);
299  params->set("printCommInfo", true);
300  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
301  }
302 
303  if(!Ac.is_null()) {std::ostringstream oss; oss << "A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
304  Set(coarseLevel, "A", Ac);
305 
306  RAPparams->set("graph", Ac);
307  Set(coarseLevel, "RAP reuse data", RAPparams);
308  }
309 
310 
311  }
312 
313 #ifdef HAVE_MUELU_DEBUG
314  MatrixUtils::checkLocalRowMapMatchesColMap(*Ac);
315 #endif // HAVE_MUELU_DEBUG
316 
317  if (transferFacts_.begin() != transferFacts_.end()) {
318  SubFactoryMonitor m(*this, "Projections", coarseLevel);
319 
320  // call Build of all user-given transfer factories
321  for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it) {
322  RCP<const FactoryBase> fac = *it;
323  GetOStream(Runtime0) << "RAPFactory: call transfer factory: " << fac->description() << std::endl;
324  fac->CallBuild(coarseLevel);
325  // Coordinates transfer is marginally different from all other operations
326  // because it is *optional*, and not required. For instance, we may need
327  // coordinates only on level 4 if we start repartitioning from that level,
328  // but we don't need them on level 1,2,3. As our current Hierarchy setup
329  // assumes propagation of dependencies only through three levels, this
330  // means that we need to rely on other methods to propagate optional data.
331  //
332  // The method currently used is through RAP transfer factories, which are
333  // simply factories which are called at the end of RAP with a single goal:
334  // transfer some fine data to coarser level. Because these factories are
335  // kind of outside of the mainline factories, they behave different. In
336  // particular, we call their Build method explicitly, rather than through
337  // Get calls. This difference is significant, as the Get call is smart
338  // enough to know when to release all factory dependencies, and Build is
339  // dumb. This led to the following CoordinatesTransferFactory sequence:
340  // 1. Request level 0
341  // 2. Request level 1
342  // 3. Request level 0
343  // 4. Release level 0
344  // 5. Release level 1
345  //
346  // The problem is missing "6. Release level 0". Because it was missing,
347  // we had outstanding request on "Coordinates", "Aggregates" and
348  // "CoarseMap" on level 0.
349  //
350  // This was fixed by explicitly calling Release on transfer factories in
351  // RAPFactory. I am still unsure how exactly it works, but now we have
352  // clear data requests for all levels.
353  coarseLevel.Release(*fac);
354  }
355  }
356 
357  }
358 
359  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
361  // check if it's a TwoLevelFactoryBase based transfer factory
362  TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::rcp_dynamic_cast<const TwoLevelFactoryBase>(factory) == Teuchos::null, Exceptions::BadCast,
363  "MueLu::RAPFactory::AddTransferFactory: Transfer factory is not derived from TwoLevelFactoryBase. "
364  "This is very strange. (Note: you can remove this exception if there's a good reason for)");
365  TEUCHOS_TEST_FOR_EXCEPTION(hasDeclaredInput_, Exceptions::RuntimeError, "MueLu::RAPFactory::AddTransferFactory: Factory is being added after we have already declared input");
366  transferFacts_.push_back(factory);
367  }
368 
369 } //namespace MueLu
370 
371 #define MUELU_RAPFACTORY_SHORT
372 #endif // MUELU_RAPFACTORY_DEF_HPP
#define SET_VALID_ENTRY(name)
MueLu::DefaultNode Node
Exception indicating invalid cast attempted.
Exception throws to report errors in the internal logical of the program.
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
void Release(const FactoryBase &factory)
Decrement the storage counter for all the inputs of a factory.
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Runtime0
One-liner description of what is happening.
@ Warnings1
Additional warnings.
static std::string getColonLabel(const std::string &label)
Helper function for object label.