46 #ifndef MUELU_TRANSPFACTORY_DEF_HPP
47 #define MUELU_TRANSPFACTORY_DEF_HPP
49 #include <Teuchos_ParameterList.hpp>
50 #include <Teuchos_Time.hpp>
52 #include <Xpetra_Matrix.hpp>
60 #include "MueLu_PerfUtils.hpp"
61 #include "MueLu_Utilities.hpp"
65 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
67 RCP<ParameterList> validParamList = rcp(
new ParameterList());
68 validParamList->set< RCP<const FactoryBase> >(
"P", Teuchos::null,
"Generating factory of the matrix P");
71 ParameterList norecurse;
72 norecurse.disableRecursiveValidation();
73 validParamList->set<ParameterList> (
"matrixmatrix: kernel params", norecurse,
"MatrixMatrix kernel parameters");
75 return validParamList;
78 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
80 Input(coarseLevel,
"P");
83 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
88 RCP<Matrix> P = Get< RCP<Matrix> >(coarseLevel,
"P");
89 const Teuchos::ParameterList& pL = GetParameterList();
92 RCP<ParameterList> Tparams;
93 if(pL.isSublist(
"matrixmatrix: kernel params"))
94 Tparams=rcp(
new ParameterList(pL.sublist(
"matrixmatrix: kernel params")));
96 Tparams= rcp(
new ParameterList);
99 Tparams->set(
"compute global constants: temporaries",Tparams->get(
"compute global constants: temporaries",
false));
100 Tparams->set(
"compute global constants",Tparams->get(
"compute global constants",
false));
105 RCP<ParameterList> params = rcp(
new ParameterList());
106 params->set(
"printLoadBalancingInfo",
true);
107 params->set(
"printCommInfo",
true);
111 Set(coarseLevel,
"R", R);
114 if (P->IsView(
"stridedMaps"))
115 R->CreateView(
"stridedMaps", P,
true);
Timer to be used in factories. Similar to Monitor but with additional timers.
Class that holds all level-specific information.
int GetLevelID() const
Return level number.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
static RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Transpose(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op, bool optimizeTranspose=false, const std::string &label=std::string(), const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.