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