MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_RAPShiftFactory_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_RAPSHIFTFACTORY_DEF_HPP
11 #define MUELU_RAPSHIFTFACTORY_DEF_HPP
12 
13 #include <sstream>
14 
15 #include <Xpetra_Matrix.hpp>
16 #include <Xpetra_MatrixMatrix.hpp>
17 #include <Xpetra_MatrixUtils.hpp>
18 #include <Xpetra_Vector.hpp>
19 #include <Xpetra_VectorFactory.hpp>
20 #include <Xpetra_MatrixFactory2.hpp>
21 
23 #include "MueLu_MasterList.hpp"
24 #include "MueLu_Monitor.hpp"
25 #include "MueLu_PerfUtils.hpp"
26 
27 namespace MueLu {
28 
29 /*********************************************************************************************************/
30 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
32  : implicitTranspose_(false) {}
33 
34 /*********************************************************************************************************/
35 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
37  RCP<ParameterList> validParamList = rcp(new ParameterList());
38 
39 #define SET_VALID_ENTRY(name) validParamList->setEntry(name, MasterList::getEntry(name))
40  SET_VALID_ENTRY("transpose: use implicit");
41  SET_VALID_ENTRY("rap: fix zero diagonals");
42  SET_VALID_ENTRY("rap: shift");
43  SET_VALID_ENTRY("rap: shift array");
44  SET_VALID_ENTRY("rap: cfl array");
45  SET_VALID_ENTRY("rap: shift diagonal M");
46  SET_VALID_ENTRY("rap: shift low storage");
47  SET_VALID_ENTRY("rap: relative diagonal floor");
48 #undef SET_VALID_ENTRY
49 
50  validParamList->set<RCP<const FactoryBase> >("A", Teuchos::null, "Generating factory of the matrix A used during the prolongator smoothing process");
51  validParamList->set<RCP<const FactoryBase> >("M", Teuchos::null, "Generating factory of the matrix M used during the non-Galerkin RAP");
52  validParamList->set<RCP<const FactoryBase> >("Mdiag", Teuchos::null, "Generating factory of the matrix Mdiag used during the non-Galerkin RAP");
53  validParamList->set<RCP<const FactoryBase> >("K", Teuchos::null, "Generating factory of the matrix K used during the non-Galerkin RAP");
54  validParamList->set<RCP<const FactoryBase> >("P", Teuchos::null, "Prolongator factory");
55  validParamList->set<RCP<const FactoryBase> >("R", Teuchos::null, "Restrictor factory");
56 
57  validParamList->set<bool>("CheckMainDiagonal", false, "Check main diagonal for zeros");
58  validParamList->set<bool>("RepairMainDiagonal", false, "Repair zeros on main diagonal");
59 
60  validParamList->set<RCP<const FactoryBase> >("deltaT", Teuchos::null, "user deltaT");
61  validParamList->set<RCP<const FactoryBase> >("cfl", Teuchos::null, "user cfl");
62  validParamList->set<RCP<const FactoryBase> >("cfl-based shift array", Teuchos::null, "MueLu-generated shift array for CFL-based shifting");
63 
64  // Make sure we don't recursively validate options for the matrixmatrix kernels
65  ParameterList norecurse;
66  norecurse.disableRecursiveValidation();
67  validParamList->set<ParameterList>("matrixmatrix: kernel params", norecurse, "MatrixMatrix kernel parameters");
68 
69  return validParamList;
70 }
71 
72 /*********************************************************************************************************/
73 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
75  const Teuchos::ParameterList &pL = GetParameterList();
76 
77  bool use_mdiag = false;
78  if (pL.isParameter("rap: shift diagonal M"))
79  use_mdiag = pL.get<bool>("rap: shift diagonal M");
80 
81  // The low storage version requires mdiag
82  bool use_low_storage = false;
83  if (pL.isParameter("rap: shift low storage")) {
84  use_low_storage = pL.get<bool>("rap: shift low storage");
85  use_mdiag = use_low_storage ? true : use_mdiag;
86  }
87 
88  if (implicitTranspose_ == false) {
89  Input(coarseLevel, "R");
90  }
91 
92  if (!use_low_storage)
93  Input(fineLevel, "K");
94  else
95  Input(fineLevel, "A");
96  Input(coarseLevel, "P");
97 
98  if (!use_mdiag)
99  Input(fineLevel, "M");
100  else
101  Input(fineLevel, "Mdiag");
102 
103  // CFL array stuff
104  if (pL.isParameter("rap: cfl array") && pL.get<Teuchos::Array<double> >("rap: cfl array").size() > 0) {
105  if (fineLevel.GetLevelID() == 0) {
106  if (fineLevel.IsAvailable("deltaT", NoFactory::get())) {
107  fineLevel.DeclareInput("deltaT", NoFactory::get(), this);
108  } else {
109  TEUCHOS_TEST_FOR_EXCEPTION(fineLevel.IsAvailable("fine deltaT", NoFactory::get()),
111  "deltaT was not provided by the user on level0!");
112  }
113 
114  if (fineLevel.IsAvailable("cfl", NoFactory::get())) {
115  fineLevel.DeclareInput("cfl", NoFactory::get(), this);
116  } else {
117  TEUCHOS_TEST_FOR_EXCEPTION(fineLevel.IsAvailable("fine cfl", NoFactory::get()),
119  "cfl was not provided by the user on level0!");
120  }
121  } else {
122  Input(fineLevel, "cfl-based shift array");
123  }
124  }
125 
126  // call DeclareInput of all user-given transfer factories
127  for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it) {
128  (*it)->CallDeclareInput(coarseLevel);
129  }
130 }
131 
132 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
133 void RAPShiftFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Build(Level &fineLevel, Level &coarseLevel) const { // FIXME make fineLevel const
134  {
135  FactoryMonitor m(*this, "Computing Ac", coarseLevel);
136  const Teuchos::ParameterList &pL = GetParameterList();
137 
138  bool M_is_diagonal = false;
139  if (pL.isParameter("rap: shift diagonal M"))
140  M_is_diagonal = pL.get<bool>("rap: shift diagonal M");
141 
142  // The low storage version requires mdiag
143  bool use_low_storage = false;
144  if (pL.isParameter("rap: shift low storage")) {
145  use_low_storage = pL.get<bool>("rap: shift low storage");
146  M_is_diagonal = use_low_storage ? true : M_is_diagonal;
147  }
148 
150  Teuchos::ArrayRCP<double> myshifts;
151  if (pL.isParameter("rap: shift array") && pL.get<Teuchos::Array<double> >("rap: shift array").size() > 0) {
152  // Do we have an array of shifts? If so, we set doubleShifts_
153  doubleShifts = pL.get<Teuchos::Array<double> >("rap: shift array")();
154  }
155  if (pL.isParameter("rap: cfl array") && pL.get<Teuchos::Array<double> >("rap: cfl array").size() > 0) {
156  // Do we have an array of CFLs? If so, we calculated the shifts from them.
157  Teuchos::ArrayView<const double> CFLs = pL.get<Teuchos::Array<double> >("rap: cfl array")();
158  if (fineLevel.GetLevelID() == 0) {
159  double dt = Get<double>(fineLevel, "deltaT");
160  double cfl = Get<double>(fineLevel, "cfl");
161  double ts_at_cfl1 = dt / cfl;
162  myshifts.resize(CFLs.size());
163  Teuchos::Array<double> myCFLs(CFLs.size());
164  myCFLs[0] = cfl;
165 
166  // Never make the CFL bigger
167  for (int i = 1; i < (int)CFLs.size(); i++)
168  myCFLs[i] = (CFLs[i] > cfl) ? cfl : CFLs[i];
169 
170  {
171  std::ostringstream ofs;
172  ofs << "RAPShiftFactory: CFL schedule = ";
173  for (int i = 0; i < (int)CFLs.size(); i++)
174  ofs << " " << myCFLs[i];
175  GetOStream(Statistics0) << ofs.str() << std::endl;
176  }
177  GetOStream(Statistics0) << "RAPShiftFactory: Timestep at CFL=1 is " << ts_at_cfl1 << " " << std::endl;
178 
179  // The shift array needs to be 1/dt
180  for (int i = 0; i < (int)myshifts.size(); i++)
181  myshifts[i] = 1.0 / (ts_at_cfl1 * myCFLs[i]);
182  doubleShifts = myshifts();
183 
184  {
185  std::ostringstream ofs;
186  ofs << "RAPShiftFactory: shift schedule = ";
187  for (int i = 0; i < (int)doubleShifts.size(); i++)
188  ofs << " " << doubleShifts[i];
189  GetOStream(Statistics0) << ofs.str() << std::endl;
190  }
191  Set(coarseLevel, "cfl-based shift array", myshifts);
192  } else {
193  myshifts = Get<Teuchos::ArrayRCP<double> >(fineLevel, "cfl-based shift array");
194  doubleShifts = myshifts();
195  Set(coarseLevel, "cfl-based shift array", myshifts);
196  // NOTE: If we're not on level zero, then we should have a shift array
197  }
198  }
199 
200  // Inputs: K, M, P
201  // Note: In the low-storage case we do not keep a separate "K", we just use A
202  RCP<Matrix> K;
203  RCP<Matrix> M;
204  RCP<Vector> Mdiag;
205 
206  if (use_low_storage)
207  K = Get<RCP<Matrix> >(fineLevel, "A");
208  else
209  K = Get<RCP<Matrix> >(fineLevel, "K");
210  if (!M_is_diagonal)
211  M = Get<RCP<Matrix> >(fineLevel, "M");
212  else
213  Mdiag = Get<RCP<Vector> >(fineLevel, "Mdiag");
214 
215  RCP<Matrix> P = Get<RCP<Matrix> >(coarseLevel, "P");
216 
217  // Build Kc = RKP, Mc = RMP
218  RCP<Matrix> KP, MP;
219 
220  // Reuse pattern if available (multiple solve)
221  // FIXME: Old style reuse doesn't work any more
222  // if (IsAvailable(coarseLevel, "AP Pattern")) {
223  // KP = Get< RCP<Matrix> >(coarseLevel, "AP Pattern");
224  // MP = Get< RCP<Matrix> >(coarseLevel, "AP Pattern");
225  // }
226 
227  {
228  SubFactoryMonitor subM(*this, "MxM: K x P", coarseLevel);
230  if (!M_is_diagonal) {
232  } else {
234  MP->leftScale(*Mdiag);
235  }
236 
237  Set(coarseLevel, "AP Pattern", KP);
238  }
239 
240  bool doOptimizedStorage = true;
241 
242  RCP<Matrix> Ac, Kc, Mc;
243 
244  // Reuse pattern if available (multiple solve)
245  // if (IsAvailable(coarseLevel, "RAP Pattern"))
246  // Ac = Get< RCP<Matrix> >(coarseLevel, "RAP Pattern");
247 
248  bool doFillComplete = true;
249  if (implicitTranspose_) {
250  SubFactoryMonitor m2(*this, "MxM: P' x (KP) (implicit)", coarseLevel);
251  Kc = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*P, true, *KP, false, Kc, GetOStream(Statistics2), doFillComplete, doOptimizedStorage);
252  Mc = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*P, true, *MP, false, Mc, GetOStream(Statistics2), doFillComplete, doOptimizedStorage);
253  } else {
254  RCP<Matrix> R = Get<RCP<Matrix> >(coarseLevel, "R");
255  SubFactoryMonitor m2(*this, "MxM: R x (KP) (explicit)", coarseLevel);
256  Kc = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*R, false, *KP, false, Kc, GetOStream(Statistics2), doFillComplete, doOptimizedStorage);
257  Mc = Xpetra::MatrixMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::Multiply(*R, false, *MP, false, Mc, GetOStream(Statistics2), doFillComplete, doOptimizedStorage);
258  }
259 
260  // Get the shift
261  // FIXME - We should really get rid of the shifts array and drive this the same way everything else works
262  // If we're using the recursive "low storage" version, we need to shift by ( \prod_{i=1}^k shift[i] - \prod_{i=1}^{k-1} shift[i]) to
263  // get the recursive relationships correct
264  int level = coarseLevel.GetLevelID();
266  if (!use_low_storage) {
267  // High Storage version
268  if (level < (int)shifts_.size())
269  shift = shifts_[level];
270  else
271  shift = Teuchos::as<Scalar>(pL.get<double>("rap: shift"));
272  } else {
273  // Low Storage Version
274  if (level < (int)shifts_.size()) {
275  if (level == 1)
276  shift = shifts_[level];
277  else {
279  for (int i = 1; i < level - 1; i++) {
280  prod1 *= shifts_[i];
281  }
282  shift = (prod1 * shifts_[level] - prod1);
283  }
284  } else if (doubleShifts.size() != 0) {
285  double d_shift = 0.0;
286  if (level < doubleShifts.size())
287  d_shift = doubleShifts[level] - doubleShifts[level - 1];
288 
289  if (d_shift < 0.0)
290  GetOStream(Warnings1) << "WARNING: RAPShiftFactory has detected a negative shift... This implies a less stable coarse grid." << std::endl;
291  shift = Teuchos::as<Scalar>(d_shift);
292  } else {
293  double base_shift = pL.get<double>("rap: shift");
294  if (level == 1)
295  shift = Teuchos::as<Scalar>(base_shift);
296  else
297  shift = Teuchos::as<Scalar>(pow(base_shift, level) - pow(base_shift, level - 1));
298  }
299  }
300  GetOStream(Runtime0) << "RAPShiftFactory: Using shift " << shift << std::endl;
301 
302  // recombine to get K+shift*M
303  {
304  SubFactoryMonitor m2(*this, "Add: RKP + s*RMP", coarseLevel);
306  Ac->fillComplete();
307  }
308 
309  Teuchos::ArrayView<const double> relativeFloor = pL.get<Teuchos::Array<double> >("rap: relative diagonal floor")();
310  if (relativeFloor.size() > 0)
312 
313  bool repairZeroDiagonals = pL.get<bool>("RepairMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
314  bool checkAc = pL.get<bool>("CheckMainDiagonal") || pL.get<bool>("rap: fix zero diagonals");
315  ;
316  if (checkAc || repairZeroDiagonals)
317  Xpetra::MatrixUtils<SC, LO, GO, NO>::CheckRepairMainDiagonal(Ac, repairZeroDiagonals, GetOStream(Warnings1));
318 
319  RCP<ParameterList> params = rcp(new ParameterList());
320  ;
321  params->set("printLoadBalancingInfo", true);
322  GetOStream(Statistics0) << PerfUtils::PrintMatrixInfo(*Ac, "Ac", params);
323 
324  Set(coarseLevel, "A", Ac);
325  // We only need K in the 'high storage' mode
326  if (!use_low_storage)
327  Set(coarseLevel, "K", Kc);
328 
329  if (!M_is_diagonal) {
330  Set(coarseLevel, "M", Mc);
331  } else {
332  // If M is diagonal, then we only pass that part down the hierarchy
333  // NOTE: Should we be doing some kind of rowsum instead?
335  Mc->getLocalDiagCopy(*Mcv);
336  Set(coarseLevel, "Mdiag", Mcv);
337  }
338 
339  // Set(coarseLevel, "RAP Pattern", Ac);
340  }
341 
342  if (transferFacts_.begin() != transferFacts_.end()) {
343  SubFactoryMonitor m(*this, "Projections", coarseLevel);
344 
345  // call Build of all user-given transfer factories
346  for (std::vector<RCP<const FactoryBase> >::const_iterator it = transferFacts_.begin(); it != transferFacts_.end(); ++it) {
347  RCP<const FactoryBase> fac = *it;
348  GetOStream(Runtime0) << "RAPShiftFactory: call transfer factory: " << fac->description() << std::endl;
349  fac->CallBuild(coarseLevel);
350  // AP (11/11/13): I am not sure exactly why we need to call Release, but we do need it to get rid
351  // of dangling data for CoordinatesTransferFactory
352  coarseLevel.Release(*fac);
353  }
354  }
355 }
356 
357 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
359  // check if it's a TwoLevelFactoryBase based transfer factory
360  TEUCHOS_TEST_FOR_EXCEPTION(Teuchos::rcp_dynamic_cast<const TwoLevelFactoryBase>(factory) == Teuchos::null, Exceptions::BadCast, "MueLu::RAPShiftFactory::AddTransferFactory: Transfer factory is not derived from TwoLevelFactoryBase. This is very strange. (Note: you can remove this exception if there's a good reason for)");
361  transferFacts_.push_back(factory);
362 }
363 
364 } // namespace MueLu
365 
366 #define MUELU_RAPSHIFTFACTORY_SHORT
367 #endif // MUELU_RAPSHIFTFACTORY_DEF_HPP
Exception indicating invalid cast attempted.
void Build(Level &fineLevel, Level &coarseLevel) const
Build an object with this factory.
virtual void CallBuild(Level &requestedLevel) const =0
void Release(const FactoryBase &factory)
Decrement the storage counter for all the inputs of a factory.
T & get(const std::string &name, T def_value)
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
size_type size() const
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 ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
size_type size() const
static void CheckRepairMainDiagonal(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> &Ac, bool const &repairZeroDiagonals, Teuchos::FancyOStream &fos, const typename Teuchos::ScalarTraits< Scalar >::magnitudeType threshold=Teuchos::ScalarTraits< typename Teuchos::ScalarTraits< Scalar >::magnitudeType >::zero())
static const NoFactory * get()
Print even more statistics.
Additional warnings.
void resize(const size_type n, const T &val=T())
bool isParameter(const std::string &name) const
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
Print statistics that do not involve significant additional computation.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
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 > > BuildCopy(const RCP< const Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > A)
static std::string PrintMatrixInfo(const Matrix &A, const std::string &msgTag, RCP< const Teuchos::ParameterList > params=Teuchos::null)
static RCP< Vector > Build(const Teuchos::RCP< const Map > &map, bool zeroOut=true)
static void TwoMatrixAdd(const Matrix &A, bool transposeA, SC alpha, Matrix &B, SC beta)
void AddTransferFactory(const RCP< const FactoryBase > &factory)
Add transfer factory in the end of list of transfer factories in RepartitionAcFactory.
static void Multiply(const Matrix &A, bool transposeA, const Matrix &B, bool transposeB, Matrix &C, bool call_FillComplete_on_result=true, bool doOptimizeStorage=true, const std::string &label=std::string(), const RCP< ParameterList > &params=null)
int GetLevelID() const
Return level number.
Definition: MueLu_Level.cpp:51
Exception throws to report errors in the internal logical of the program.
void DeclareInput(const std::string &ename, const FactoryBase *factory, const FactoryBase *requestedBy=NoFactory::get())
Callback from FactoryBase::CallDeclareInput() and FactoryBase::DeclareInput()
virtual std::string description() const
Return a simple one-line description of this object.
bool IsAvailable(const std::string &ename, const FactoryBase *factory=NoFactory::get()) const
Test whether a need&#39;s value has been saved.
#define SET_VALID_ENTRY(name)