MueLu  Version of the Day
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  return validParamList;
80  }
81 
82  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
84  Input(fineLevel, "A");
85 
86  // Get default tentative prolongator factory
87  // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
88  RCP<const FactoryBase> initialPFact = GetFactory("P");
89  if (initialPFact == Teuchos::null) { initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent"); }
90  coarseLevel.DeclareInput("P", initialPFact.get(), this); // --
91  }
92 
93  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
95  return BuildP(fineLevel, coarseLevel);
96  }
97 
98  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
100  FactoryMonitor m(*this, "Prolongator smoothing", coarseLevel);
101  std::ostringstream levelstr;
102  levelstr << coarseLevel.GetLevelID();
103 
104  typedef typename Teuchos::ScalarTraits<SC>::magnitudeType Magnitude;
105 
106  // Get default tentative prolongator factory
107  // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
108  // -- Warning: Do not use directly initialPFact_. Use initialPFact instead everywhere!
109  RCP<const FactoryBase> initialPFact = GetFactory("P");
110  if (initialPFact == Teuchos::null) { initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent"); }
111 
112  // Level Get
113  RCP<Matrix> A = Get< RCP<Matrix> >(fineLevel, "A");
114  RCP<Matrix> Ptent = coarseLevel.Get< RCP<Matrix> >("P", initialPFact.get());
115 
116  if(restrictionMode_) {
117  SubFactoryMonitor m2(*this, "Transpose A", coarseLevel);
118  A = Utilities::Transpose(*A, true); // build transpose of A explicitely
119  }
120 
121  //Build final prolongator
122  RCP<Matrix> finalP; // output
123 
124  const ParameterList & pL = GetParameterList();
125  Scalar dampingFactor = as<Scalar>(pL.get<double>("sa: damping factor"));
126  LO maxEigenIterations = as<LO>(pL.get<int>("sa: eigenvalue estimate num iterations"));
127  bool estimateMaxEigen = pL.get<bool>("sa: calculate eigenvalue estimate");
128  if (dampingFactor != Teuchos::ScalarTraits<Scalar>::zero()) {
129 
130  Scalar lambdaMax;
131  {
132  SubFactoryMonitor m2(*this, "Eigenvalue estimate", coarseLevel);
133  lambdaMax = A->GetMaxEigenvalueEstimate();
134  if (lambdaMax == -Teuchos::ScalarTraits<SC>::one() || estimateMaxEigen) {
135  GetOStream(Statistics1) << "Calculating max eigenvalue estimate now (max iters = "<< maxEigenIterations << ")" << std::endl;
136  Magnitude stopTol = 1e-4;
137  lambdaMax = Utilities::PowerMethod(*A, true, maxEigenIterations, stopTol);
138  A->SetMaxEigenvalueEstimate(lambdaMax);
139  } else {
140  GetOStream(Statistics1) << "Using cached max eigenvalue estimate" << std::endl;
141  }
142  GetOStream(Statistics0) << "Prolongator damping factor = " << dampingFactor/lambdaMax << " (" << dampingFactor << " / " << lambdaMax << ")" << std::endl;
143  }
144 
145  {
146  SubFactoryMonitor m2(*this, "Fused (I-omega*D^{-1} A)*Ptent", coarseLevel);
147  Teuchos::RCP<Vector> invDiag = Utilities::GetMatrixDiagonalInverse(*A);
148 
149  SC omega = dampingFactor / lambdaMax;
150 
151  // finalP = Ptent + (I - \omega D^{-1}A) Ptent
152  finalP = Xpetra::IteratorOps<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Jacobi(omega, *invDiag, *A, *Ptent, finalP, GetOStream(Statistics2),std::string("MueLu::SaP-")+levelstr.str());
153  }
154 
155  } else {
156  finalP = Ptent;
157  }
158 
159  // Level Set
160  if (!restrictionMode_) {
161  // prolongation factory is in prolongation mode
162  Set(coarseLevel, "P", finalP);
163 
164  // NOTE: EXPERIMENTAL
165  if (Ptent->IsView("stridedMaps"))
166  finalP->CreateView("stridedMaps", Ptent);
167 
168  } else {
169  // prolongation factory is in restriction mode
170  RCP<Matrix> R = Utilities::Transpose(*finalP, true);
171  Set(coarseLevel, "R", R);
172 
173  // NOTE: EXPERIMENTAL
174  if (Ptent->IsView("stridedMaps"))
175  R->CreateView("stridedMaps", Ptent, true);
176  }
177 
178  if (IsPrint(Statistics1)) {
179  RCP<ParameterList> params = rcp(new ParameterList());
180  params->set("printLoadBalancingInfo", true);
181  params->set("printCommInfo", true);
182  GetOStream(Statistics1) << PerfUtils::PrintMatrixInfo(*finalP, (!restrictionMode_ ? "P" : "R"), params);
183  }
184 
185  } //Build()
186 
187 } //namespace MueLu
188 
189 #endif // MUELU_SAPFACTORY_DEF_HPP
190 
191 //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->Get< RCP<Matrix> >("A", factory) if factory == NULL => use default 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())
Timer to be used in factories. Similar to Monitor but with additional timers.
Print more statistics.
Namespace for MueLu classes and methods.
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.
Print statistics that do not involve significant additional computation.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:99
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 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.
const RCP< const FactoryManagerBase > GetFactoryManager()
returns the current factory manager
Definition: MueLu_Level.cpp:96
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
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)