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 <Xpetra_IO.hpp>
52 #include <sstream>
53 
55 
57 #include "MueLu_Level.hpp"
58 #include "MueLu_MasterList.hpp"
59 #include "MueLu_Monitor.hpp"
60 #include "MueLu_PerfUtils.hpp"
62 #include "MueLu_TentativePFactory.hpp"
63 #include "MueLu_Utilities.hpp"
64 
65 namespace MueLu {
66 
67  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
69  RCP<ParameterList> validParamList = rcp(new ParameterList());
70 
71 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
72  SET_VALID_ENTRY("sa: damping factor");
73  SET_VALID_ENTRY("sa: calculate eigenvalue estimate");
74  SET_VALID_ENTRY("sa: eigenvalue estimate num iterations");
75  SET_VALID_ENTRY("sa: use rowsumabs diagonal scaling");
76  SET_VALID_ENTRY("sa: enforce constraints");
77  SET_VALID_ENTRY("tentative: calculate qr");
78  SET_VALID_ENTRY("sa: max eigenvalue");
79  SET_VALID_ENTRY("sa: rowsumabs diagonal replacement tolerance");
80  SET_VALID_ENTRY("sa: rowsumabs diagonal replacement value");
81 #undef SET_VALID_ENTRY
82 
83  validParamList->set< RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A used during the prolongator smoothing process");
84  validParamList->set< RCP<const FactoryBase> >("P", Teuchos::null, "Tentative prolongator factory");
85 
86  // Make sure we don't recursively validate options for the matrixmatrix kernels
87  ParameterList norecurse;
88  norecurse.disableRecursiveValidation();
89  validParamList->set<ParameterList> ("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
90 
91 
92  return validParamList;
93  }
94 
95  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
97  Input(fineLevel, "A");
98 
99  // Get default tentative prolongator factory
100  // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
101  RCP<const FactoryBase> initialPFact = GetFactory("P");
102  if (initialPFact == Teuchos::null) { initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent"); }
103  coarseLevel.DeclareInput("P", initialPFact.get(), this); // --
104  }
105 
106  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
108  return BuildP(fineLevel, coarseLevel);
109  }
110 
111  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
113  FactoryMonitor m(*this, "Prolongator smoothing", coarseLevel);
114 
115  std::string levelIDs = toString(coarseLevel.GetLevelID());
116 
117  const std::string prefix = "MueLu::SaPFactory(" + levelIDs + "): ";
118 
119  typedef typename Teuchos::ScalarTraits<SC>::coordinateType Coordinate;
120 
121  // Get default tentative prolongator factory
122  // Getting it that way ensure that the same factory instance will be used for both SaPFactory and NullspaceFactory.
123  // -- Warning: Do not use directly initialPFact_. Use initialPFact instead everywhere!
124  RCP<const FactoryBase> initialPFact = GetFactory("P");
125  if (initialPFact == Teuchos::null) { initialPFact = coarseLevel.GetFactoryManager()->GetFactory("Ptent"); }
126  const ParameterList& pL = GetParameterList();
127 
128  // Level Get
129  RCP<Matrix> A = Get< RCP<Matrix> >(fineLevel, "A");
130  RCP<Matrix> Ptent = coarseLevel.Get< RCP<Matrix> >("P", initialPFact.get());
131 
132  if (restrictionMode_) {
133  SubFactoryMonitor m2(*this, "Transpose A", coarseLevel);
134  A = Utilities::Transpose(*A, true); // build transpose of A explicitely
135  }
136 
137  // Build final prolongator
138  RCP<Matrix> finalP;
139 
140  // Reuse pattern if available
141  RCP<ParameterList> APparams = rcp(new ParameterList);;
142  if(pL.isSublist("matrixmatrix: kernel params"))
143  APparams->sublist("matrixmatrix: kernel params") = pL.sublist("matrixmatrix: kernel params");
144 
145  if (coarseLevel.IsAvailable("AP reuse data", this)) {
146  GetOStream(static_cast<MsgType>(Runtime0 | Test)) << "Reusing previous AP data" << std::endl;
147 
148  APparams = coarseLevel.Get< RCP<ParameterList> >("AP reuse data", this);
149 
150  if (APparams->isParameter("graph"))
151  finalP = APparams->get< RCP<Matrix> >("graph");
152  }
153  // By default, we don't need global constants for SaP
154  APparams->set("compute global constants: temporaries",APparams->get("compute global constants: temporaries",false));
155  APparams->set("compute global constants",APparams->get("compute global constants",false));
156 
157  const SC dampingFactor = as<SC>(pL.get<double>("sa: damping factor"));
158  const LO maxEigenIterations = as<LO>(pL.get<int> ("sa: eigenvalue estimate num iterations"));
159  const bool estimateMaxEigen = pL.get<bool> ("sa: calculate eigenvalue estimate");
160  const bool useAbsValueRowSum= pL.get<bool> ("sa: use rowsumabs diagonal scaling");
161  const bool doQRStep = pL.get<bool>("tentative: calculate qr");
162  const bool enforceConstraints= pL.get<bool>("sa: enforce constraints");
163  const SC userDefinedMaxEigen = as<SC>(pL.get<double>("sa: max eigenvalue"));
164  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType Magnitude;
165  double dTol = pL.get<double>("sa: rowsumabs diagonal replacement tolerance");
166  const Magnitude diagonalReplacementTolerance = (dTol == as<double>(-1) ? Teuchos::ScalarTraits<Scalar>::eps()*100 : as<Magnitude>(pL.get<double>("sa: rowsumabs diagonal replacement tolerance")));
167  const SC diagonalReplacementValue = as<SC>(pL.get<double>("sa: rowsumabs diagonal replacement value"));
168 
169  // Sanity checking
170  TEUCHOS_TEST_FOR_EXCEPTION(doQRStep && enforceConstraints,Exceptions::RuntimeError,
171  "MueLu::TentativePFactory::MakeTentative: cannot use 'enforce constraints' and 'calculate qr' at the same time");
172 
173  if (dampingFactor != Teuchos::ScalarTraits<SC>::zero()) {
174 
175  Scalar lambdaMax;
176  Teuchos::RCP<Vector> invDiag;
177  if (userDefinedMaxEigen == -1.)
178  {
179  SubFactoryMonitor m2(*this, "Eigenvalue estimate", coarseLevel);
180  lambdaMax = A->GetMaxEigenvalueEstimate();
181  if (lambdaMax == -Teuchos::ScalarTraits<SC>::one() || estimateMaxEigen) {
182  GetOStream(Statistics1) << "Calculating max eigenvalue estimate now (max iters = "<< maxEigenIterations <<
183  ( (useAbsValueRowSum) ? ", use rowSumAbs diagonal)" : ", use point diagonal)") << std::endl;
184  Coordinate stopTol = 1e-4;
185  if (useAbsValueRowSum) {
186  const bool returnReciprocal=true;
187  const bool replaceSingleEntryRowWithZero=true;
188  invDiag = Utilities::GetLumpedMatrixDiagonal(*A,returnReciprocal,
189  diagonalReplacementTolerance,
190  diagonalReplacementValue,
191  replaceSingleEntryRowWithZero);
192  TEUCHOS_TEST_FOR_EXCEPTION(invDiag.is_null(), Exceptions::RuntimeError,
193  "SaPFactory: eigenvalue estimate: diagonal reciprocal is null.");
194  lambdaMax = Utilities::PowerMethod(*A, invDiag, maxEigenIterations, stopTol);
195  } else
196  lambdaMax = Utilities::PowerMethod(*A, true, maxEigenIterations, stopTol);
197  A->SetMaxEigenvalueEstimate(lambdaMax);
198  } else {
199  GetOStream(Statistics1) << "Using cached max eigenvalue estimate" << std::endl;
200  }
201  }
202  else {
203  lambdaMax = userDefinedMaxEigen;
204  A->SetMaxEigenvalueEstimate(lambdaMax);
205  }
206  GetOStream(Statistics1) << "Prolongator damping factor = " << dampingFactor/lambdaMax << " (" << dampingFactor << " / " << lambdaMax << ")" << std::endl;
207 
208  {
209  SubFactoryMonitor m2(*this, "Fused (I-omega*D^{-1} A)*Ptent", coarseLevel);
210  if (!useAbsValueRowSum)
211  invDiag = Utilities::GetMatrixDiagonalInverse(*A); //default
212  else if (invDiag == Teuchos::null) {
213  GetOStream(Runtime0) << "Using rowsumabs diagonal" << std::endl;
214  const bool returnReciprocal=true;
215  const bool replaceSingleEntryRowWithZero=true;
216  invDiag = Utilities::GetLumpedMatrixDiagonal(*A,returnReciprocal,
217  diagonalReplacementTolerance,
218  diagonalReplacementValue,
219  replaceSingleEntryRowWithZero);
220  TEUCHOS_TEST_FOR_EXCEPTION(invDiag.is_null(), Exceptions::RuntimeError, "SaPFactory: diagonal reciprocal is null.");
221  }
222 
223  SC omega = dampingFactor / lambdaMax;
224  TEUCHOS_TEST_FOR_EXCEPTION(!std::isfinite(Teuchos::ScalarTraits<SC>::magnitude(omega)), Exceptions::RuntimeError, "Prolongator damping factor needs to be finite.");
225 
226  // finalP = Ptent + (I - \omega D^{-1}A) Ptent
227  finalP = Xpetra::IteratorOps<Scalar,LocalOrdinal,GlobalOrdinal,Node>::Jacobi(omega, *invDiag, *A, *Ptent, finalP,
228  GetOStream(Statistics2), std::string("MueLu::SaP-")+levelIDs, APparams);
229  if (enforceConstraints) SatisfyPConstraints( A, finalP);
230  }
231 
232  } else {
233  finalP = Ptent;
234  }
235 
236  // Level Set
237  if (!restrictionMode_) {
238  // The factory is in prolongation mode
239  if(!finalP.is_null()) {std::ostringstream oss; oss << "P_" << coarseLevel.GetLevelID(); finalP->setObjectLabel(oss.str());}
240  Set(coarseLevel, "P", finalP);
241 
242  APparams->set("graph", finalP);
243  Set(coarseLevel, "AP reuse data", APparams);
244 
245  // NOTE: EXPERIMENTAL
246  if (Ptent->IsView("stridedMaps"))
247  finalP->CreateView("stridedMaps", Ptent);
248 
249  if (IsPrint(Statistics2)) {
250  RCP<ParameterList> params = rcp(new ParameterList());
251  params->set("printLoadBalancingInfo", true);
252  params->set("printCommInfo", true);
253  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*finalP, "P", params);
254  }
255 
256  } else {
257  // The factory is in restriction mode
258  RCP<Matrix> R;
259  {
260  SubFactoryMonitor m2(*this, "Transpose P", coarseLevel);
261  R = Utilities::Transpose(*finalP, true);
262  if(!R.is_null()) {std::ostringstream oss; oss << "R_" << coarseLevel.GetLevelID(); R->setObjectLabel(oss.str());}
263  }
264 
265  Set(coarseLevel, "R", R);
266 
267  // NOTE: EXPERIMENTAL
268  if (Ptent->IsView("stridedMaps"))
269  R->CreateView("stridedMaps", Ptent, true/*transposeA*/);
270 
271  if (IsPrint(Statistics2)) {
272  RCP<ParameterList> params = rcp(new ParameterList());
273  params->set("printLoadBalancingInfo", true);
274  params->set("printCommInfo", true);
275  GetOStream(Statistics2) << PerfUtils::PrintMatrixInfo(*R, "R", params);
276  }
277  }
278 
279  } //Build()
280 
281  // Analyze the grid transfer produced by smoothed aggregation and make
282  // modifications if it does not look right. In particular, if there are
283  // negative entries or entries larger than 1, modify P's rows.
284  //
285  // Note: this kind of evaluation probably only makes sense if not doing QR
286  // when constructing tentative P.
287  //
288  // For entries that do not satisfy the two constraints (>= 0 or <=1) we set
289  // these entries to the constraint value and modify the rest of the row
290  // so that the row sum remains the same as before by adding an equal
291  // amount to each remaining entry. However, if the original row sum value
292  // violates the constraints, we set the row sum back to 1 (the row sum of
293  // tentative P). After doing the modification to a row, we need to check
294  // again the entire row to make sure that the modified row does not violate
295  // the constraints.
296 
297  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
298  void SaPFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::SatisfyPConstraints(const RCP<Matrix> A, RCP<Matrix>& P) const {
299 
300  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
301  const Scalar one = Teuchos::ScalarTraits<Scalar>::one();
302  LO nPDEs = A->GetFixedBlockSize();
303  Teuchos::ArrayRCP<Scalar> ConstraintViolationSum(nPDEs);
304  Teuchos::ArrayRCP<Scalar> Rsum(nPDEs);
305  Teuchos::ArrayRCP<size_t> nPositive(nPDEs);
306  for (size_t k=0; k < (size_t) nPDEs; k++) ConstraintViolationSum[k] = zero;
307  for (size_t k=0; k < (size_t) nPDEs; k++) Rsum[k] = zero;
308  for (size_t k=0; k < (size_t) nPDEs; k++) nPositive[k] = 0;
309 
310 
311  for (size_t i = 0; i < as<size_t>(P->getRowMap()->getNodeNumElements()); i++) {
312 
313  Teuchos::ArrayView<const LocalOrdinal> indices;
314  Teuchos::ArrayView<const Scalar> vals1;
315  Teuchos::ArrayView< Scalar> vals;
316  P->getLocalRowView((LO) i, indices, vals1);
317  size_t nnz = indices.size();
318  if (nnz == 0) continue;
319 
320  vals = ArrayView<Scalar>(const_cast<SC*>(vals1.getRawPtr()), nnz);
321 
322 
323  bool checkRow = true;
324 
325  while (checkRow) {
326 
327  // check constraints and compute the row sum
328 
329  for (LO j = 0; j < indices.size(); j++) {
330  Rsum[ j%nPDEs ] += vals[j];
331  if (Teuchos::ScalarTraits<SC>::real(vals[j]) < Teuchos::ScalarTraits<SC>::real(zero)) {
332  ConstraintViolationSum[ j%nPDEs ] += vals[j];
333  vals[j] = zero;
334  }
335  else {
336  if (Teuchos::ScalarTraits<SC>::real(vals[j]) != Teuchos::ScalarTraits<SC>::real(zero))
337  (nPositive[ j%nPDEs])++;
338 
339  if (Teuchos::ScalarTraits<SC>::real(vals[j]) > Teuchos::ScalarTraits<SC>::real(1.00001 )) {
340  ConstraintViolationSum[ j%nPDEs ] += (vals[j] - one);
341  vals[j] = one;
342  }
343  }
344  }
345 
346  checkRow = false;
347 
348  // take into account any row sum that violates the contraints
349 
350  for (size_t k=0; k < (size_t) nPDEs; k++) {
351 
352  if (Teuchos::ScalarTraits<SC>::real(Rsum[ k ]) < Teuchos::ScalarTraits<SC>::magnitude(zero)) {
353  ConstraintViolationSum[k] += (-Rsum[k]); // rstumin
354  }
355  else if (Teuchos::ScalarTraits<SC>::real(Rsum[ k ]) > Teuchos::ScalarTraits<SC>::magnitude(1.00001)) {
356  ConstraintViolationSum[k] += (one - Rsum[k]); // rstumin
357  }
358  }
359 
360  // check if row need modification
361  for (size_t k=0; k < (size_t) nPDEs; k++) {
362  if (Teuchos::ScalarTraits<SC>::magnitude(ConstraintViolationSum[ k ]) != Teuchos::ScalarTraits<SC>::magnitude(zero))
363  checkRow = true;
364  }
365  // modify row
366  if (checkRow) {
367  for (LO j = 0; j < indices.size(); j++) {
368  if (Teuchos::ScalarTraits<SC>::real(vals[j]) > Teuchos::ScalarTraits<SC>::real(zero)) {
369  vals[j] += (ConstraintViolationSum[j%nPDEs]/ (as<Scalar>(nPositive[j%nPDEs])));
370  }
371  }
372  for (size_t k=0; k < (size_t) nPDEs; k++) ConstraintViolationSum[k] = zero;
373  }
374  for (size_t k=0; k < (size_t) nPDEs; k++) Rsum[k] = zero;
375  for (size_t k=0; k < (size_t) nPDEs; k++) nPositive[k] = 0;
376  } // while (checkRow) ...
377  } // for (size_t i = 0; i < as<size_t>(P->getRowMap()->getNumNodeElements()); i++) ...
378  } //SatsifyPConstraints()
379 
380 } //namespace MueLu
381 
382 #endif // MUELU_SAPFACTORY_DEF_HPP
383 
384 //TODO: restrictionMode_ should use the parameter list.
void SatisfyPConstraints(RCP< Matrix > A, RCP< Matrix > &P) const
Enforce constraints on prolongator.
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.
std::string toString(const T &what)
Little helper function to convert non-string types to strings.
scalar_type dampingFactor
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
Timer to be used in factories. Similar to Monitor but with additional timers.
Print more statistics.
One-liner description of what is happening.
Namespace for MueLu classes and methods.
static Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetLumpedMatrixDiagonal(Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > const &A, const bool doReciprocal=false, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar tolReplacement=Teuchos::ScalarTraits< Scalar >::zero(), const bool replaceSingleEntryRowWithZero=false)
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
Print even more statistics.
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:76
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
MueLu::DefaultScalar Scalar
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.
#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)
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
const RCP< const FactoryManagerBase > GetFactoryManager()
returns the current factory manager
Definition: MueLu_Level.cpp:96
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.
Exception throws to report errors in the internal logical of the program.
static RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > GetMatrixDiagonalInverse(const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > &A, Magnitude tol=Teuchos::ScalarTraits< Scalar >::eps() *100, Scalar valReplacement=Teuchos::ScalarTraits< Scalar >::zero())
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
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)