10#ifndef MUELU_EMINPFACTORY_DEF_HPP
11#define MUELU_EMINPFACTORY_DEF_HPP
13#include <Xpetra_Matrix.hpp>
14#include <Xpetra_StridedMapFactory.hpp>
18#include "MueLu_CGSolver.hpp"
19#include "MueLu_Constraint.hpp"
20#include "MueLu_GMRESSolver.hpp"
23#include "MueLu_PerfUtils.hpp"
25#include "MueLu_SteepestDescentSolver.hpp"
29template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
31 RCP<ParameterList> validParamList = rcp(
new ParameterList());
33#define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
38 validParamList->getEntry(
"emin: iterative method").setValidator(rcp(
new Teuchos::StringValidator(Teuchos::tuple<std::string>(
"cg",
"sd",
"gmres"))));
42 validParamList->set<RCP<const FactoryBase> >(
"A", Teuchos::null,
"Generating factory for the matrix A used during internal iterations");
43 validParamList->set<RCP<const FactoryBase> >(
"P", Teuchos::null,
"Generating factory for the initial guess");
44 validParamList->set<RCP<const FactoryBase> >(
"Constraint", Teuchos::null,
"Generating factory for constraints");
46 validParamList->set<RCP<Matrix> >(
"P0", Teuchos::null,
"Initial guess at P");
47 validParamList->set<
bool>(
"Keep P0",
false,
"Keep an initial P0 (for reuse)");
49 validParamList->set<RCP<Constraint> >(
"Constraint0", Teuchos::null,
"Initial Constraint");
50 validParamList->set<
bool>(
"Keep Constraint0",
false,
"Keep an initial Constraint (for reuse)");
52 return validParamList;
55template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
57 Input(fineLevel,
"A");
59 static bool isAvailableP0 =
false;
60 static bool isAvailableConstraint0 =
false;
74 isAvailableP0 = coarseLevel.
IsAvailable(
"P0",
this);
75 isAvailableConstraint0 = coarseLevel.
IsAvailable(
"Constraint0",
this);
78 if (isAvailableP0 ==
false)
79 Input(coarseLevel,
"P");
81 if (isAvailableConstraint0 ==
false)
82 Input(coarseLevel,
"Constraint");
85template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
87 BuildP(fineLevel, coarseLevel);
90template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
110 P0 = coarseLevel.
Get<RCP<Matrix> >(
"P0",
this);
111 numIts = pL.get<
int>(
"emin: num reuse iterations");
117 numIts = pL.get<
int>(
"emin: num iterations");
125 if (coarseLevel.
IsAvailable(
"Constraint0",
this)) {
127 X = coarseLevel.
Get<RCP<Constraint> >(
"Constraint0",
this);
136 std::string solverType = pL.get<std::string>(
"emin: iterative method");
137 RCP<SolverBase> solver;
138 if (solverType ==
"cg")
140 else if (solverType ==
"sd")
142 else if (solverType ==
"gmres")
146 solver->Iterate(*A, *X, *P0, P);
149 if (!P->IsView(
"stridedMaps")) {
150 if (A->IsView(
"stridedMaps") ==
true) {
157 std::vector<size_t> stridingInfo(1, 1);
158 RCP<const StridedMap> dMap = StridedMapFactory::Build(X->GetPattern()->getDomainMap(), stridingInfo);
160 P->CreateView(
"stridedMaps", A->getRowMap(
"stridedMaps"), dMap);
163 P->CreateView(
"stridedMaps", P->getRangeMap(), P->getDomainMap());
170 Set(coarseLevel,
"P", P);
172 if (pL.get<
bool>(
"Keep P0")) {
176 coarseLevel.
Keep(
"P0",
this);
177 Set(coarseLevel,
"P0", P);
179 if (pL.get<
bool>(
"Keep Constraint0")) {
183 coarseLevel.
Keep(
"Constraint0",
this);
184 Set(coarseLevel,
"Constraint0", X);
188 RCP<ParameterList> params = rcp(
new ParameterList());
189 params->set(
"printLoadBalancingInfo",
true);
190 params->set(
"printCommInfo",
true);
203 Set(coarseLevel,
"R", R);
206 RCP<ParameterList> params = rcp(
new ParameterList());
207 params->set(
"printLoadBalancingInfo",
true);
208 params->set(
"printCommInfo",
true);
#define SET_VALID_ENTRY(name)
Implements conjugate gradient algorithm for energy-minimization.
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
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.
Timer to be used in factories. Similar to Monitor but with additional timers.
void Input(Level &level, const std::string &varName) const
T Get(Level &level, const std::string &varName) const
void Set(Level &level, const std::string &varName, const T &data) const
Implements conjugate gradient algorithm for energy-minimization.
Class that holds all level-specific information.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need's value has been saved.
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access)....
RequestMode GetRequestMode() const
void Keep(const std::string &ename, const FactoryBase *factory)
Request to keep variable 'ename' generated by 'factory' after the setup phase.
virtual const Teuchos::ParameterList & GetParameterList() const
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
Implements steepest descent algorithm for energy-minimization.
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
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)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
bool IsPrint(MsgType type, int thisProcRankOnly=-1) const
Find out whether we need to print out information for a specific message type.
Namespace for MueLu classes and methods.
@ Statistics2
Print even more statistics.
@ Runtime0
One-liner description of what is happening.
@ Runtime1
Description of what is happening (more verbose).