10#ifndef MUELU_TRANSPFACTORY_DEF_HPP
11#define MUELU_TRANSPFACTORY_DEF_HPP
13#include <Teuchos_ParameterList.hpp>
14#include <Teuchos_Time.hpp>
16#include <Xpetra_Matrix.hpp>
21#include "MueLu_PerfUtils.hpp"
22#include "MueLu_Utilities.hpp"
26template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
28 RCP<ParameterList> validParamList = rcp(
new ParameterList());
29 validParamList->set<RCP<const FactoryBase> >(
"P", Teuchos::null,
"Generating factory of the matrix P");
32 ParameterList norecurse;
33 norecurse.disableRecursiveValidation();
34 validParamList->set<ParameterList>(
"matrixmatrix: kernel params", norecurse,
"MatrixMatrix kernel parameters");
36 return validParamList;
39template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
41 Input(coarseLevel,
"P");
44template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
47 std::string label =
"MueLu::TransP-" + Teuchos::toString(coarseLevel.
GetLevelID());
49 RCP<Matrix> P = Get<RCP<Matrix> >(coarseLevel,
"P");
52 if (P == Teuchos::null)
return;
54 const Teuchos::ParameterList& pL = GetParameterList();
57 RCP<ParameterList> Tparams;
58 if (pL.isSublist(
"matrixmatrix: kernel params"))
59 Tparams = rcp(
new ParameterList(pL.sublist(
"matrixmatrix: kernel params")));
61 Tparams = rcp(
new ParameterList);
64 Tparams->set(
"compute global constants: temporaries", Tparams->get(
"compute global constants: temporaries",
false));
65 Tparams->set(
"compute global constants", Tparams->get(
"compute global constants",
false));
70 RCP<ParameterList> params = rcp(
new ParameterList());
71 params->set(
"printLoadBalancingInfo",
true);
72 params->set(
"printCommInfo",
true);
76 Set(coarseLevel,
"R", R);
79 if (P->IsView(
"stridedMaps"))
80 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.