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;
123 std::ostringstream levelstr;
128 "MueLu::RAPFactory::Build(): CallDeclareInput has not been called before Build!");
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;
134 bool isEpetra = A->getRowMap()->lib() == Xpetra::UseEpetra;
135 #ifdef KOKKOS_ENABLE_CUDA
136 bool isCuda =
typeid(
Node).name() ==
typeid(Kokkos::Compat::KokkosCudaWrapperNode).name();
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";
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");
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));
158 if (coarseLevel.IsAvailable(
"AP reuse data",
this)) {
159 GetOStream(
static_cast<MsgType>(
Runtime0 |
Test)) <<
"Reusing previous AP data" << std::endl;
161 APparams = coarseLevel.Get< RCP<ParameterList> >(
"AP reuse data",
this);
163 if (APparams->isParameter(
"graph"))
164 AP = APparams->get< RCP<Matrix> >(
"graph");
170 AP = MatrixMatrix::Multiply(*A, !doTranspose, *P, !doTranspose, AP, GetOStream(
Statistics2),
171 doFillComplete, doOptimizeStorage, labelstr+std::string(
"MueLu::A*P-")+levelstr.str(), APparams);
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");
179 if (coarseLevel.IsAvailable(
"RAP reuse data",
this)) {
180 GetOStream(
static_cast<MsgType>(
Runtime0 |
Test)) <<
"Reusing previous RAP data" << std::endl;
182 RAPparams = coarseLevel.Get< RCP<ParameterList> >(
"RAP reuse data",
this);
184 if (RAPparams->isParameter(
"graph"))
185 Ac = RAPparams->get< RCP<Matrix> >(
"graph");
189 Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
193 RAPparams->set(
"compute global constants: temporaries",RAPparams->get(
"compute global constants: temporaries",
false));
194 RAPparams->set(
"compute global constants",
true);
200 if (pL.get<
bool>(
"transpose: use implicit") ==
true) {
203 Ac = MatrixMatrix::Multiply(*P, doTranspose, *AP, !doTranspose, Ac, GetOStream(
Statistics2),
204 doFillComplete, doOptimizeStorage, labelstr+std::string(
"MueLu::R*(AP)-implicit-")+levelstr.str(), RAPparams);
207 RCP<Matrix> R = Get< RCP<Matrix> >(coarseLevel,
"R");
211 Ac = MatrixMatrix::Multiply(*R, !doTranspose, *AP, !doTranspose, Ac, GetOStream(
Statistics2),
212 doFillComplete, doOptimizeStorage, labelstr+std::string(
"MueLu::R*(AP)-explicit-")+levelstr.str(), RAPparams);
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));
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));
226 RCP<ParameterList> params = rcp(
new ParameterList());;
227 params->set(
"printLoadBalancingInfo",
true);
228 params->set(
"printCommInfo",
true);
232 if(!Ac.is_null()) {std::ostringstream oss; oss <<
"A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
233 Set(coarseLevel,
"A", Ac);
235 APparams->set(
"graph", AP);
236 Set(coarseLevel,
"AP reuse data", APparams);
237 RAPparams->set(
"graph", Ac);
238 Set(coarseLevel,
"RAP reuse data", RAPparams);
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");
244 if (coarseLevel.IsAvailable(
"RAP reuse data",
this)) {
245 GetOStream(
static_cast<MsgType>(
Runtime0 |
Test)) <<
"Reusing previous RAP data" << std::endl;
247 RAPparams = coarseLevel.Get< RCP<ParameterList> >(
"RAP reuse data",
this);
249 if (RAPparams->isParameter(
"graph"))
250 Ac = RAPparams->get< RCP<Matrix> >(
"graph");
254 Ac->SetMaxEigenvalueEstimate(-Teuchos::ScalarTraits<SC>::one());
258 RAPparams->set(
"compute global constants: temporaries",RAPparams->get(
"compute global constants: temporaries",
false));
259 RAPparams->set(
"compute global constants",
true);
261 if (pL.get<
bool>(
"transpose: use implicit") ==
true) {
263 Ac = MatrixFactory::Build(P->getDomainMap(), Teuchos::as<LO>(0));
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(),
273 RCP<Matrix> R = Get< RCP<Matrix> >(coarseLevel,
"R");
274 Ac = MatrixFactory::Build(R->getRowMap(), Teuchos::as<LO>(0));
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(),
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));
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));
297 RCP<ParameterList> params = rcp(
new ParameterList());;
298 params->set(
"printLoadBalancingInfo",
true);
299 params->set(
"printCommInfo",
true);
303 if(!Ac.is_null()) {std::ostringstream oss; oss <<
"A_" << coarseLevel.GetLevelID(); Ac->setObjectLabel(oss.str());}
304 Set(coarseLevel,
"A", Ac);
306 RAPparams->set(
"graph", Ac);
307 Set(coarseLevel,
"RAP reuse data", RAPparams);
313 #ifdef HAVE_MUELU_DEBUG
314 MatrixUtils::checkLocalRowMapMatchesColMap(*Ac);
317 if (transferFacts_.begin() != transferFacts_.end()) {
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);
359 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
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);
371 #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.