MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_SaPFactory_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // MueLu: A package for multigrid based preconditioning
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 #ifndef MUELU_SAPFACTORY_DEF_HPP
47 #define MUELU_SAPFACTORY_DEF_HPP
48 
49 #include <Xpetra_Matrix.hpp>
50 #include <Xpetra_IteratorOps.hpp> // containing routines to generate Jacobi iterator
51 #include <sstream>
52 
54 
56 #include "MueLu_Level.hpp"
57 #include "MueLu_MasterList.hpp"
58 #include "MueLu_Monitor.hpp"
59 #include "MueLu_PerfUtils.hpp"
61 #include "MueLu_TentativePFactory.hpp"
62 #include "MueLu_Utilities.hpp"
63 
64 namespace MueLu {
65 
66  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68  RCP<ParameterList> validParamList = rcp(new ParameterList());
69 
70 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
71  SET_VALID_ENTRY("sa: damping factor");
72  SET_VALID_ENTRY("sa: calculate eigenvalue estimate");
73  SET_VALID_ENTRY("sa: eigenvalue estimate num iterations");
74 #undef SET_VALID_ENTRY
75 
76  validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A used during the prolongator smoothing process");
77  validParamList->set< RCP<const FactoryBase> >("P", Teuchos::null, "Tentative prolongator factory");
78 
79  // Make sure we don't recursively validate options for the matrixmatrix kernels
80  ParameterList norecurse;
81  norecurse.disableRecursiveValidation();
82  validParamList->set<ParameterList> ("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
83 
84 
85  return validParamList;
86  }
87 
88  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
90  Input(fineLevel, "A");
91 
92  // Get default tentative prolongator factory
93  // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
94  RCP<const FactoryBase> initialPFact = GetFactory("P");
95  if (initialPFact == Teuchos::null) { initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent"); }
96  coarseLevel.DeclareInput("P", initialPFact.get(), this); // --
97  }
98 
99  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
101  return BuildP(fineLevel, coarseLevel);
102  }
103 
104  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
106  FactoryMonitor m(*this, "Prolongator smoothing", coarseLevel);
107 
108  std::string levelIDs = toString(coarseLevel.GetLevelID());
109 
110  const std::string prefix = "MueLu::SaPFactory(" + levelIDs + "): ";
111 
112  typedef typename Teuchos::ScalarTraits<SC>::coordinateType Coordinate;
113 
114  // Get default tentative prolongator factory
115  // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
116  // -- Warning: Do not use directly initialPFact_. Use initialPFact instead everywhere!
117  RCP<const FactoryBase> initialPFact = GetFactory("P");
118  if (initialPFact == Teuchos::null) { initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent"); }
119  const ParameterList& pL = GetParameterList();
120 
121  // Level Get
122  RCP<Matrix> A = Get< RCP<Matrix> >(fineLevel, "A");
123  RCP<Matrix> Ptent = coarseLevel.Get< RCP<Matrix> >("P", initialPFact.get());
124 
125  if (restrictionMode_) {
126  SubFactoryMonitor m2(*this, "Transpose A", coarseLevel);
127  A = Utilities::Transpose(*A, true); // build transpose of A explicitely
128  }
129 
130  // Build final prolongator
131  RCP<Matrix> finalP;
132 
133  // Reuse pattern if available
134  RCP<ParameterList> APparams = rcp(new ParameterList);;
135  if(pL.isSublist("matrixmatrix: kernel params"))
136  APparams->sublist("matrixmatrix: kernel params") = pL.sublist("matrixmatrix: kernel params");
137 
138  if (coarseLevel.IsAvailable("AP reuse data", this)) {
139  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous AP data" << std::endl;
140 
141  APparams = coarseLevel.Get< RCP<ParameterList> >("AP reuse data", this);
142 
143  if (APparams->isParameter("graph"))
144  finalP = APparams->get< RCP<Matrix> >("graph");
145  }
146  // By default, we don't need global constants for SaP
147  APparams->set("compute global constants: temporaries",APparams->get("compute global constants: temporaries",false));
148  APparams->set("compute global constants",APparams->get("compute global constants",false));
149 
150  SC dampingFactor = as<SC>(pL.get<double>("sa: damping factor"));
151  LO maxEigenIterations = as<LO>(pL.get<int> ("sa: eigenvalue estimate num iterations"));
152  bool estimateMaxEigen = pL.get<bool> ("sa: calculate eigenvalue estimate");
153  if (dampingFactor != Teuchos::ScalarTraits<SC>::zero()) {
154 
155  Scalar lambdaMax;
156  {
157  SubFactoryMonitor m2(*this, "Eigenvalue estimate", coarseLevel);
158  lambdaMax = A->GetMaxEigenvalueEstimate();
159  if (lambdaMax == -Teuchos::ScalarTraits<SC>::one() || estimateMaxEigen) {
160  GetOStream(Statistics1) << "Calculating max eigenvalue estimate now (max iters = "<< maxEigenIterations << ")" << std::endl;
161  Coordinate stopTol = 1e-4;
162  lambdaMax = Utilities::PowerMethod(*A, true, maxEigenIterations, stopTol);
163  A->SetMaxEigenvalueEstimate(lambdaMax);
164  } else {
165  GetOStream(Statistics1) << "Using cached max eigenvalue estimate" << std::endl;
166  }
167  GetOStream(Statistics1) << "Prolongator damping factor = " << dampingFactor/lambdaMax << " (" << dampingFactor << " / " << lambdaMax << ")" << std::endl;
168  }
169 
170  {
171  SubFactoryMonitor m2(*this, "Fused (I-omega*D^{-1} A)*Ptent", coarseLevel);
173 
174  SC omega = dampingFactor / lambdaMax;
175  TEUCHOS_TEST_FOR_EXCEPTION(!std::isfinite(Teuchos::ScalarTraits<SC>::magnitude(omega)), Exceptions::RuntimeError, "Prolongator damping factor needs to be finite.");
176 
177  // finalP = Ptent + (I - \omega D^{-1}A) Ptent
178  finalP = Xpetra::IteratorOps<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Jacobi(omega, *invDiag, *A, *Ptent, finalP,
179  GetOStream(Statistics2), std::string("MueLu::SaP-")+levelIDs, APparams);
180  }
181 
182  } else {
183  finalP = Ptent;
184  }
185 
186  // Level Set
187  if (!restrictionMode_) {
188  // The factory is in prolongation mode
189  if(!finalP.is_null()) {std::ostringstream oss; oss << "P_" << coarseLevel.GetLevelID(); finalP->setObjectLabel(oss.str());}
190  Set(coarseLevel, "P", finalP);
191 
192  APparams->set("graph", finalP);
193  Set(coarseLevel, "AP reuse data", APparams);
194 
195  // NOTE: EXPERIMENTAL
196  if (Ptent->IsView("stridedMaps"))
197  finalP->CreateView("stridedMaps", Ptent);
198 
199  if (IsPrint(Statistics2)) {
200  RCP<ParameterList> params = rcp(new ParameterList());
201  params->set("printLoadBalancingInfo", true);
202  params->set("printCommInfo", true);
203  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*finalP, "P", params);
204  }
205 
206  } else {
207  // The factory is in restriction mode
208  RCP<Matrix> R;
209  {
210  SubFactoryMonitor m2(*this, "Transpose P", coarseLevel);
211  R = Utilities::Transpose(*finalP, true);
212  if(!R.is_null()) {std::ostringstream oss; oss << "R_" << coarseLevel.GetLevelID(); R->setObjectLabel(oss.str());}
213  }
214 
215  Set(coarseLevel, "R", R);
216 
217  // NOTE: EXPERIMENTAL
218  if (Ptent->IsView("stridedMaps"))
219  R->CreateView("stridedMaps", Ptent, true/*transposeA*/);
220 
221  if (IsPrint(Statistics2)) {
222  RCP<ParameterList> params = rcp(new ParameterList());
223  params->set("printLoadBalancingInfo", true);
224  params->set("printCommInfo", true);
225  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*R, "R", params);
226  }
227  }
228 
229  } //Build()
230 
231 } //namespace MueLu
232 
233 #endif // MUELU_SAPFACTORY_DEF_HPP
234 
235 //TODO: restrictionMode_ should use the parameter list.
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.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Print more statistics.
LocalOrdinal LO
One-liner description of what is happening.
T * get() const
static RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetMatrixDiagonalInverse(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps()*100)
Print even more statistics.
bool isParameter(const std::string &name) const
static RCP< Matrix > Jacobi(SC omega, const Vector &Dinv, const Matrix &A, const Matrix &B, RCP< Matrix > C_in, Teuchos::FancyOStream &fos, const std::string &label, RCP< ParameterList > &params)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
bool isSublist(const std::string &name) const
Timer to be used in factories. Similar to SubMonitor but adds a timer level by level.
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
#define SET_VALID_ENTRY(name)
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 BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Scalar SC
const RCP< const FactoryManagerBase > GetFactoryManager()
returns the current factory manager
Definition: MueLu_Level.cpp:96
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
Exception throws to report errors in the internal logical of the program.
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
static Scalar PowerMethod(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, bool scaleByDiag=true, LocalOrdinal niters=10, Magnitude tolerance=1e-2, bool verbose=false, unsigned int seed=123)
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.
bool is_null() const