MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_EminPFactory_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef MUELU_EMINPFACTORY_DEF_HPP
11 #define MUELU_EMINPFACTORY_DEF_HPP
12 
13 #include <Xpetra_Matrix.hpp>
14 #include <Xpetra_StridedMapFactory.hpp>
15 
17 
18 #include "MueLu_CGSolver.hpp"
19 #include "MueLu_Constraint.hpp"
20 #include "MueLu_GMRESSolver.hpp"
21 #include "MueLu_MasterList.hpp"
22 #include "MueLu_Monitor.hpp"
23 #include "MueLu_PerfUtils.hpp"
24 #include "MueLu_SolverBase.hpp"
25 #include "MueLu_SteepestDescentSolver.hpp"
26 
27 namespace MueLu {
28 
29 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31  RCP<ParameterList> validParamList = rcp(new ParameterList());
32 
33 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
34  SET_VALID_ENTRY("emin: num iterations");
35  SET_VALID_ENTRY("emin: num reuse iterations");
36  SET_VALID_ENTRY("emin: iterative method");
37  {
38  validParamList->getEntry("emin: iterative method").setValidator(rcp(new Teuchos::StringValidator(Teuchos::tuple<std::string>("cg", "sd", "gmres"))));
39  }
40 #undef SET_VALID_ENTRY
41 
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");
45 
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)");
48 
49  validParamList->set<RCP<Constraint> >("Constraint0", Teuchos::null, "Initial Constraint");
50  validParamList->set<bool>("Keep Constraint0", false, "Keep an initial Constraint (for reuse)");
51 
52  return validParamList;
53 }
54 
55 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
57  Input(fineLevel, "A");
58 
59  static bool isAvailableP0 = false;
60  static bool isAvailableConstraint0 = false;
61 
62  // Here is a tricky little piece of code
63  // We don't want to request (aka call Input) when we reuse and P0 is available
64  // However, we cannot run something simple like this:
65  // if (!coarseLevel.IsAvailable("P0", this))
66  // Input(coarseLevel, "P");
67  // The reason is that it works fine during the request stage, but fails
68  // in the release stage as we _construct_ P0 during Build process. Therefore,
69  // we need to understand whether we are in Request or Release mode
70  // NOTE: This is a very unique situation, please try not to propagate the
71  // mode check any further
72 
73  if (coarseLevel.GetRequestMode() == Level::REQUEST) {
74  isAvailableP0 = coarseLevel.IsAvailable("P0", this);
75  isAvailableConstraint0 = coarseLevel.IsAvailable("Constraint0", this);
76  }
77 
78  if (isAvailableP0 == false)
79  Input(coarseLevel, "P");
80 
81  if (isAvailableConstraint0 == false)
82  Input(coarseLevel, "Constraint");
83 }
84 
85 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
87  BuildP(fineLevel, coarseLevel);
88 }
89 
90 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
92  FactoryMonitor m(*this, "Prolongator minimization", coarseLevel);
93 
94  const ParameterList& pL = GetParameterList();
95 
96  // Get the matrix
97  RCP<Matrix> A = Get<RCP<Matrix> >(fineLevel, "A");
98 
99  if (restrictionMode_) {
100  SubFactoryMonitor m2(*this, "Transpose A", coarseLevel);
101 
102  A = Utilities::Transpose(*A, true);
103  }
104 
105  // Get/make initial guess
106  RCP<Matrix> P0;
107  int numIts;
108  if (coarseLevel.IsAvailable("P0", this)) {
109  // Reuse data
110  P0 = coarseLevel.Get<RCP<Matrix> >("P0", this);
111  numIts = pL.get<int>("emin: num reuse iterations");
112  GetOStream(Runtime0) << "Reusing P0" << std::endl;
113 
114  } else {
115  // Construct data
116  P0 = Get<RCP<Matrix> >(coarseLevel, "P");
117  numIts = pL.get<int>("emin: num iterations");
118  }
119  // NOTE: the main assumption here that P0 satisfies both constraints:
120  // - nonzero pattern
121  // - nullspace preservation
122 
123  // Get/make constraint operator
124  RCP<Constraint> X;
125  if (coarseLevel.IsAvailable("Constraint0", this)) {
126  // Reuse data
127  X = coarseLevel.Get<RCP<Constraint> >("Constraint0", this);
128  GetOStream(Runtime0) << "Reusing Constraint0" << std::endl;
129 
130  } else {
131  // Construct data
132  X = Get<RCP<Constraint> >(coarseLevel, "Constraint");
133  }
134  GetOStream(Runtime0) << "Number of emin iterations = " << numIts << std::endl;
135 
136  std::string solverType = pL.get<std::string>("emin: iterative method");
137  RCP<SolverBase> solver;
138  if (solverType == "cg")
139  solver = rcp(new CGSolver(numIts));
140  else if (solverType == "sd")
141  solver = rcp(new SteepestDescentSolver(numIts));
142  else if (solverType == "gmres")
143  solver = rcp(new GMRESSolver(numIts));
144 
145  RCP<Matrix> P;
146  solver->Iterate(*A, *X, *P0, P);
147 
148  // NOTE: EXPERIMENTAL and FRAGILE
149  if (!P->IsView("stridedMaps")) {
150  if (A->IsView("stridedMaps") == true) {
151  GetOStream(Runtime1) << "Using A to fillComplete P" << std::endl;
152 
153  // FIXME: X->GetPattern() actually returns a CrsGraph.
154  // CrsGraph has no knowledge of Xpetra's sup/Matrix views. As such,
155  // it has no idea about strided maps. We create one, which is
156  // most likely incorrect for many use cases.
157  std::vector<size_t> stridingInfo(1, 1);
158  RCP<const StridedMap> dMap = StridedMapFactory::Build(X->GetPattern()->getDomainMap(), stridingInfo);
159 
160  P->CreateView("stridedMaps", A->getRowMap("stridedMaps"), dMap);
161 
162  } else {
163  P->CreateView("stridedMaps", P->getRangeMap(), P->getDomainMap());
164  }
165  }
166 
167  // Level Set
168  if (!restrictionMode_) {
169  // The factory is in prolongation mode
170  Set(coarseLevel, "P", P);
171 
172  if (pL.get<bool>("Keep P0")) {
173  // NOTE: we must do Keep _before_ set as the Needs class only sets if
174  // a) data has been requested (which is not the case here), or
175  // b) data has some keep flag
176  coarseLevel.Keep("P0", this);
177  Set(coarseLevel, "P0", P);
178  }
179  if (pL.get<bool>("Keep Constraint0")) {
180  // NOTE: we must do Keep _before_ set as the Needs class only sets if
181  // a) data has been requested (which is not the case here), or
182  // b) data has some keep flag
183  coarseLevel.Keep("Constraint0", this);
184  Set(coarseLevel, "Constraint0", X);
185  }
186 
187  if (IsPrint(Statistics2)) {
188  RCP<ParameterList> params = rcp(new ParameterList());
189  params->set("printLoadBalancingInfo", true);
190  params->set("printCommInfo", true);
191  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*P, "P", params);
192  }
193 
194  } else {
195  // The factory is in restriction mode
196  RCP<Matrix> R;
197  {
198  SubFactoryMonitor m2(*this, "Transpose P", coarseLevel);
199 
200  R = Utilities::Transpose(*P, true);
201  }
202 
203  Set(coarseLevel, "R", R);
204 
205  if (IsPrint(Statistics2)) {
206  RCP<ParameterList> params = rcp(new ParameterList());
207  params->set("printLoadBalancingInfo", true);
208  params->set("printCommInfo", true);
209  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*R, "R", params);
210  }
211  }
212 }
213 
214 } // namespace MueLu
215 
216 #endif // MUELU_EMINPFACTORY_DEF_HPP
void Keep(const std::string &ename, const FactoryBase *factory)
Request to keep variable &#39;ename&#39; generated by &#39;factory&#39; after the setup phase.
#define SET_VALID_ENTRY(name)
T & Get(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Get data without decrementing associated storage counter (i.e., read-only access). Usage: Level-&gt;Get&lt; RCP&lt;Matrix&gt; &gt;(&quot;A&quot;, factory) if factory == NULL =&gt; use default factory.
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
void setValidator(RCP< const ParameterEntryValidator > const &validator)
T & get(const std::string &name, T def_value)
Implements conjugate gradient algorithm for energy-minimization.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Timer to be used in factories. Similar to Monitor but with additional timers.
One-liner description of what is happening.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
RCP< const CrsGraph > GetPattern() const
Print even more statistics.
RequestMode GetRequestMode() const
Implements conjugate gradient algorithm for energy-minimization.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
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 > &params=Teuchos::null)
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.
Implements steepest descent algorithm for energy-minimization.
Description of what is happening (more verbose)
ParameterEntry & getEntry(const std::string &name)
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.