46 #ifndef MUELU_RAPFACTORY_DEF_HPP
47 #define MUELU_RAPFACTORY_DEF_HPP
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>
64 #include "MueLu_PerfUtils.hpp"
70 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
72 : hasDeclaredInput_(false) { }
74 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
76 RCP<ParameterList> validParamList = rcp(
new ParameterList());
78 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
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");
88 validParamList->set<
bool > (
"CheckMainDiagonal",
false,
"Check main diagonal for zeros");
89 validParamList->set<
bool > (
"RepairMainDiagonal",
false,
"Repair zeros on main diagonal");
92 ParameterList norecurse;
93 norecurse.disableRecursiveValidation();
94 validParamList->set<ParameterList> (
"matrixmatrix: kernel params", norecurse,
"MatrixMatrix kernel parameters");
96 return validParamList;
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");
105 Input(fineLevel,
"A");
106 Input(coarseLevel,
"P");
109 for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it)
110 (*it)->CallDeclareInput(coarseLevel);
112 hasDeclaredInput_ =
true;
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;
122 std::ostringstream levelstr;
127 "MueLu::RAPFactory::Build(): CallDeclareInput has not been called before Build!");
129 const Teuchos::ParameterList& pL = GetParameterList();
130 RCP<Matrix> A = Get< RCP<Matrix> >(fineLevel,
"A");
131 RCP<Matrix> P = Get< RCP<Matrix> >(coarseLevel,
"P"), AP, Ac;
133 bool isEpetra = A->getRowMap()->lib() == Xpetra::UseEpetra;
134 #ifdef KOKKOS_ENABLE_CUDA
135 bool isCuda =
typeid(
Node).name() ==
typeid(Kokkos::Compat::KokkosCudaWrapperNode).name();
140 if (pL.get<
bool>(
"rap: triple product") ==
false || isEpetra || isCuda) {
141 if (pL.get<
bool>(
"rap: triple product") && isEpetra)
142 GetOStream(
Warnings1) <<
"Switching from triple product to R x (A x P) since triple product has not been implemented for Epetra.\n";
143 #ifdef KOKKOS_ENABLE_CUDA
144 if (pL.get<
bool>(
"rap: triple product") && isCuda)
145 GetOStream(
Warnings1) <<
"Switching from triple product to R x (A x P) since triple product has not been implemented for Cuda.\n";
149 RCP<ParameterList> APparams = rcp(
new ParameterList);
150 if(pL.isSublist(
"matrixmatrix: kernel params"))
151 APparams->sublist(
"matrixmatrix: kernel params") = pL.sublist(
"matrixmatrix: kernel params");
154 APparams->set(
"compute global constants: temporaries",APparams->get(
"compute global constants: temporaries",
false));
155 APparams->set(
"compute global constants",APparams->get(
"compute global constants",
false));
157 if (coarseLevel.IsAvailable(
"AP reuse data",
this)) {
158 GetOStream(
static_cast<MsgType>(
Runtime0 |
Test)) <<
"Reusing previous AP data" << std::endl;
160 APparams = coarseLevel.Get< RCP<ParameterList> >(
"AP reuse data",
this);
162 if (APparams->isParameter(
"graph"))
163 AP = APparams->get< RCP<Matrix> >(
"graph");
169 AP = MatrixMatrix::Multiply(*A, !doTranspose, *P, !doTranspose, AP, GetOStream(
Statistics2),
170 doFillComplete, doOptimizeStorage, labelstr+std::string(
"MueLu::A*P-")+levelstr.str(), APparams);
174 RCP<ParameterList> RAPparams = rcp(
new ParameterList);
175 if(pL.isSublist(
"matrixmatrix: kernel params"))
176 RAPparams->sublist(
"matrixmatrix: kernel params") = pL.sublist(
"matrixmatrix: kernel params");
178 if (coarseLevel.IsAvailable(
"RAP reuse data",
this)) {
179 GetOStream(
static_cast<MsgType>(
Runtime0 |
Test)) <<
"Reusing previous RAP data" << std::endl;
181 RAPparams = coarseLevel.Get< RCP<ParameterList> >(
"RAP reuse data",
this);
183 if (RAPparams->isParameter(
"graph"))
184 Ac = RAPparams->get< RCP<Matrix> >(
"graph");
188 Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
192 RAPparams->set(
"compute global constants: temporaries",RAPparams->get(
"compute global constants: temporaries",
false));
193 RAPparams->set(
"compute global constants",
true);
199 if (pL.get<
bool>(
"transpose: use implicit") ==
true) {
202 Ac = MatrixMatrix::Multiply(*P, doTranspose, *AP, !doTranspose, Ac, GetOStream(
Statistics2),
203 doFillComplete, doOptimizeStorage, labelstr+std::string(
"MueLu::R*(AP)-implicit-")+levelstr.str(), RAPparams);
206 RCP<Matrix> R = Get< RCP<Matrix> >(coarseLevel,
"R");
210 Ac = MatrixMatrix::Multiply(*R, !doTranspose, *AP, !doTranspose, Ac, GetOStream(
Statistics2),
211 doFillComplete, doOptimizeStorage, labelstr+std::string(
"MueLu::R*(AP)-explicit-")+levelstr.str(), RAPparams);
214 Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >(
"rap: relative diagonal floor")();
215 if(relativeFloor.size() > 0) {
216 Xpetra::MatrixUtils<SC,LO,GO,NO>::RelativeDiagonalBoost(Ac, relativeFloor,GetOStream(
Statistics2));
219 bool repairZeroDiagonals = pL.get<
bool>(
"RepairMainDiagonal") || pL.get<
bool>(
"rap: fix zero diagonals");
220 bool checkAc = pL.get<
bool>(
"CheckMainDiagonal")|| pL.get<
bool>(
"rap: fix zero diagonals"); ;
221 if (checkAc || repairZeroDiagonals)
222 Xpetra::MatrixUtils<SC,LO,GO,NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(
Warnings1));
225 RCP<ParameterList> params = rcp(
new ParameterList());;
226 params->set(
"printLoadBalancingInfo",
true);
227 params->set(
"printCommInfo",
true);
231 if(!Ac.is_null()) {std::ostringstream oss; oss <<
"A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
232 Set(coarseLevel,
"A", Ac);
234 APparams->set(
"graph", AP);
235 Set(coarseLevel,
"AP reuse data", APparams);
236 RAPparams->set(
"graph", Ac);
237 Set(coarseLevel,
"RAP reuse data", RAPparams);
239 RCP<ParameterList> RAPparams = rcp(
new ParameterList);
240 if(pL.isSublist(
"matrixmatrix: kernel params"))
241 RAPparams->sublist(
"matrixmatrix: kernel params") = pL.sublist(
"matrixmatrix: kernel params");
243 if (coarseLevel.IsAvailable(
"RAP reuse data",
this)) {
244 GetOStream(
static_cast<MsgType>(
Runtime0 |
Test)) <<
"Reusing previous RAP data" << std::endl;
246 RAPparams = coarseLevel.Get< RCP<ParameterList> >(
"RAP reuse data",
this);
248 if (RAPparams->isParameter(
"graph"))
249 Ac = RAPparams->get< RCP<Matrix> >(
"graph");
253 Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
257 RAPparams->set(
"compute global constants: temporaries",RAPparams->get(
"compute global constants: temporaries",
false));
258 RAPparams->set(
"compute global constants",
true);
260 if (pL.get<
bool>(
"transpose: use implicit") ==
true) {
262 Ac = MatrixFactory::Build(P->getDomainMap(), Teuchos::as<LO>(0));
266 Xpetra::TripleMatrixMultiply<SC,LO,GO,NO>::
267 MultiplyRAP(*P, doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
268 doOptimizeStorage, labelstr+std::string(
"MueLu::R*A*P-implicit-")+levelstr.str(),
272 RCP<Matrix> R = Get< RCP<Matrix> >(coarseLevel,
"R");
273 Ac = MatrixFactory::Build(R->getRowMap(), Teuchos::as<LO>(0));
277 Xpetra::TripleMatrixMultiply<SC,LO,GO,NO>::
278 MultiplyRAP(*R, !doTranspose, *A, !doTranspose, *P, !doTranspose, *Ac, doFillComplete,
279 doOptimizeStorage, labelstr+std::string(
"MueLu::R*A*P-explicit-")+levelstr.str(),
283 Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >(
"rap: relative diagonal floor")();
284 if(relativeFloor.size() > 0) {
285 Xpetra::MatrixUtils<SC,LO,GO,NO>::RelativeDiagonalBoost(Ac, relativeFloor,GetOStream(
Statistics2));
288 bool repairZeroDiagonals = pL.get<
bool>(
"RepairMainDiagonal") || pL.get<
bool>(
"rap: fix zero diagonals");
289 bool checkAc = pL.get<
bool>(
"CheckMainDiagonal")|| pL.get<
bool>(
"rap: fix zero diagonals"); ;
290 if (checkAc || repairZeroDiagonals)
291 Xpetra::MatrixUtils<SC,LO,GO,NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(
Warnings1));
296 RCP<ParameterList> params = rcp(
new ParameterList());;
297 params->set(
"printLoadBalancingInfo",
true);
298 params->set(
"printCommInfo",
true);
302 if(!Ac.is_null()) {std::ostringstream oss; oss <<
"A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
303 Set(coarseLevel,
"A", Ac);
305 RAPparams->set(
"graph", Ac);
306 Set(coarseLevel,
"RAP reuse data", RAPparams);
312 if (transferFacts_.begin() != transferFacts_.end()) {
316 for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it) {
317 RCP<const FactoryBase> fac = *it;
318 GetOStream(
Runtime0) <<
"RAPFactory: call transfer factory: " << fac->description() << std::endl;
319 fac->CallBuild(coarseLevel);
354 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
357 TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::rcp_dynamic_cast<const TwoLevelFactoryBase>(factory) == Teuchos::null,
Exceptions::BadCast,
358 "MueLu::RAPFactory::AddTransferFactory: Transfer factory is not derived from TwoLevelFactoryBase. "
359 "This is very strange. (Note: you can remove this exception if there's a good reason for)");
360 TEUCHOS_TEST_FOR_EXCEPTION(hasDeclaredInput_,
Exceptions::RuntimeError,
"MueLu::RAPFactory::AddTransferFactory: Factory is being added after we have already declared input");
361 transferFacts_.push_back(factory);
366 #define MUELU_RAPFACTORY_SHORT
#define SET_VALID_ENTRY(name)
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.
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.