MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_Hierarchy_decl.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_HIERARCHY_DECL_HPP
11 #define MUELU_HIERARCHY_DECL_HPP
12 
14 #include <Teuchos_Ptr.hpp>
15 
16 #include <Xpetra_ConfigDefs.hpp> // global_size_t
17 #include <Xpetra_Matrix_fwd.hpp>
20 #include <Xpetra_Operator_fwd.hpp>
21 
22 #include <MueLu_SmootherCloner.hpp>
23 #include "MueLu_ConfigDefs.hpp"
24 #include "MueLu_BaseClass.hpp"
25 #include "MueLu_Hierarchy_fwd.hpp"
26 
27 #include "MueLu_Types.hpp"
28 
30 #include "MueLu_FactoryManager.hpp" // no fwd declaration because constructor of FactoryManager is used as a default parameter of Setup()
31 #include "MueLu_KeepType.hpp"
32 #include "MueLu_Level_fwd.hpp"
33 #include "MueLu_MasterList.hpp"
34 #include "MueLu_NoFactory.hpp"
35 #include "MueLu_PerfUtils_fwd.hpp"
36 #include "MueLu_PFactory_fwd.hpp"
39 #include "MueLu_Utilities_fwd.hpp"
40 
41 namespace MueLu {
42 
43 enum class ConvergenceStatus {
44  Converged,
46  Undefined
47 };
48 
59 template <class Scalar = DefaultScalar,
62  class Node = DefaultNode>
63 class Hierarchy : public BaseClass {
64 #undef MUELU_HIERARCHY_SHORT
65 #include "MueLu_UseShortNames.hpp"
66 
69 
71  struct ConvData {
73  : maxIts_(1)
74  , tol_(-STS::magnitude(STS::one())) {}
75  ConvData(LO maxIts)
76  : maxIts_(maxIts)
77  , tol_(-STS::magnitude(STS::one())) {}
79  : maxIts_(10000)
80  , tol_(tol) {}
81  ConvData(std::pair<LO, MagnitudeType> p)
82  : maxIts_(p.first)
83  , tol_(p.second) {}
84 
87  };
88 
89  public:
91 
92 
94  Hierarchy();
96  Hierarchy(const std::string& label);
97 
99  Hierarchy(const RCP<Matrix>& A);
100 
102  Hierarchy(const RCP<Matrix>& A, const std::string& label);
103 
105  virtual ~Hierarchy();
106 
108 
110 
111 
113  static CycleType GetDefaultCycle() { return MasterList::getDefault<std::string>("cycle type") == "V" ? VCYCLE : WCYCLE; }
114  static int GetDefaultCycleStartLevel() { return MasterList::getDefault<int>("W cycle start level"); }
115  static bool GetDefaultImplicitTranspose() { return MasterList::getDefault<bool>("transpose: use implicit"); }
116  static bool GetDefaultFuseProlongationAndUpdate() { return MasterList::getDefault<bool>("fuse prolongation and update"); }
117  static Xpetra::global_size_t GetDefaultMaxCoarseSize() { return MasterList::getDefault<int>("coarse: max size"); }
118  static int GetDefaultMaxLevels() { return MasterList::getDefault<int>("max levels"); }
119  static bool GetDefaultPRrebalance() { return MasterList::getDefault<bool>("repartition: rebalance P and R"); }
120 
122  bool GetImplicitTranspose() const { return implicitTranspose_; }
124 
125  void SetMaxCoarseSize(Xpetra::global_size_t maxCoarseSize) { maxCoarseSize_ = maxCoarseSize; }
126  void SetPRrebalance(bool doPRrebalance) { doPRrebalance_ = doPRrebalance; }
127  void SetPRViaCopyrebalance(bool doPRViaCopyrebalance) { doPRViaCopyrebalance_ = doPRViaCopyrebalance; }
128  void SetImplicitTranspose(const bool& implicit) { implicitTranspose_ = implicit; }
129  void SetFuseProlongationAndUpdate(const bool& fuse) { fuseProlongationAndUpdate_ = fuse; }
130 
132 
134 
135  template <class S2, class LO2, class GO2, class N2>
136  friend class Hierarchy;
137 
138  int LastLevelID() const { return Levels_.size() - 1; }
139 
140  private:
141  void DumpCurrentGraph(int level) const;
142 
143  public:
145  void AddLevel(const RCP<Level>& level);
146 
148  void AddNewLevel();
149 
151  RCP<Level>& GetLevel(const int levelID = 0);
152 
153  int GetNumLevels() const;
154  int GetGlobalNumLevels() const;
155 
156  MagnitudeType GetRate() const { return rate_; }
157 
158  // This function is global
159  double GetOperatorComplexity() const;
160 
161  // This function is global
162  double GetSmootherComplexity() const;
163 
165  void CheckLevel(Level& level, int levelID);
166 
167  void SetMatvecParams(RCP<ParameterList> matvecParams);
168 
170 
207  bool Setup(int coarseLevelID, const RCP<const FactoryManagerBase> fineLevelManager /* = Teuchos::null */, const RCP<const FactoryManagerBase> coarseLevelManager,
208  const RCP<const FactoryManagerBase> nextLevelManager = Teuchos::null);
209 
211  void Setup(const FactoryManagerBase& manager = FactoryManager(), int startLevel = 0, int numDesiredLevels = GetDefaultMaxLevels());
212 
213  void SetupRe();
214 
216 
218  void Clear(int startLevel = 0);
219  void ExpertClear();
220 
222  CycleType GetCycle() const { return Cycle_; }
223 
225  void SetCycle(CycleType Cycle) { Cycle_ = Cycle; }
226 
227  void SetCycleStartLevel(int cycleStart) { WCycleStartLevel_ = cycleStart; }
228 
230  void SetProlongatorScalingFactor(double scalingFactor) { scalingFactor_ = scalingFactor; }
231 
244  ConvergenceStatus Iterate(const MultiVector& B, MultiVector& X, ConvData conv = ConvData(),
245  bool InitialGuessIsZero = false, LO startLevel = 0);
246 
256  void Write(const LO& start = -1, const LO& end = -1, const std::string& suffix = "");
257 
259 
261 
262 
264  void Keep(const std::string& ename, const FactoryBase* factory = NoFactory::get());
265 
267  void Delete(const std::string& ename, const FactoryBase* factory = NoFactory::get());
268 
270  void AddKeepFlag(const std::string& ename, const FactoryBase* factory = NoFactory::get(), KeepType keep = MueLu::Keep);
271 
273  void RemoveKeepFlag(const std::string& ename, const FactoryBase* factory, KeepType keep = MueLu::All);
274 
276 
278 
279 
281  std::string description() const;
282 
288  void describe(Teuchos::FancyOStream& out, const VerbLevel verbLevel = Default) const;
290 
292  void print(std::ostream& out = std::cout, const VerbLevel verbLevel = (MueLu::Parameters | MueLu::Statistics0)) const;
293 
298  void IsPreconditioner(const bool flag);
299 
301 
302  void EnableGraphDumping(const std::string& filename, int levelID = 1) {
303  isDumpingEnabled_ = true;
304  dumpLevel_ = levelID;
305  dumpFile_ = filename;
306  }
307 
308  void setlib(Xpetra::UnderlyingLib inlib) { lib_ = inlib; }
310 
313  description_ = "";
314  }
315 
316  void AllocateLevelMultiVectors(int numvecs, bool forceMapCheck = false);
318 
319  protected:
320  const RCP<const FactoryManagerBase>& GetLevelManager(const int levelID) const {
321  return levelManagers_[levelID];
322  }
323 
324  private:
326  Hierarchy(const Hierarchy& h);
327 
329  bool IsCalculationOfResidualRequired(const LO startLevel, const ConvData& conv) const;
330 
339  const MagnitudeType convergenceTolerance) const;
340 
342  void PrintResidualHistory(const LO iteration,
343  const Teuchos::Array<MagnitudeType>& residualNorm) const;
344 
346  ConvergenceStatus ComputeResidualAndPrintHistory(const Operator& A, const MultiVector& X,
347  const MultiVector& B, const LO iteration,
348  const LO startLevel, const ConvData& conv, MagnitudeType& previousResidualNorm);
349 
352 
357  void ReplaceCoordinateMap(Level& level);
358 
362 
366 
370 
374  bool doPRViaCopyrebalance_; // fully explicit, needed for CombinePFactory
375 
378 
381 
384 
387 
390 
392  mutable std::string description_ = ""; // mutable so that we can lazily initialize in description(), which is declared const
393 
400  // -1 = dump all levels, -2 = dump nothing
402  std::string dumpFile_;
403 
406 
409 
413 
414 }; // class Hierarchy
415 
416 } // namespace MueLu
417 
418 #define MUELU_HIERARCHY_SHORT
419 #endif // MUELU_HIERARCHY_DECL_HPP
void IsPreconditioner(const bool flag)
Array< RCP< MultiVector > > coarseExport_
This class specifies the default factory that should generate some data on a Level if the data does n...
MueLu::DefaultLocalOrdinal LocalOrdinal
RCP< Level > & GetLevel(const int levelID=0)
Retrieve a certain level from hierarchy.
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
Hierarchy()
Default constructor.
Array< RCP< MultiVector > > coarseImport_
static bool GetDefaultFuseProlongationAndUpdate()
Print parameters.
void AddLevel(const RCP< Level > &level)
Add a level at the end of the hierarchy.
void CheckForEmptySmoothersAndCoarseSolve()
void PrintResidualHistory(const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const
Print residualNorm for this iteration to the screen.
short KeepType
int sizeOfAllocatedLevelMultiVectors_
Caching (Multi)Vectors used in Hierarchy::Iterate()
int WCycleStartLevel_
Level at which to start W-cycle.
void SetMaxCoarseSize(Xpetra::global_size_t maxCoarseSize)
void DumpCurrentGraph(int level) const
Xpetra::global_size_t GetMaxCoarseSize() const
void AddNewLevel()
Add a new level at the end of the hierarchy.
LocalOrdinal LO
static Xpetra::global_size_t GetDefaultMaxCoarseSize()
void SetCycle(CycleType Cycle)
Supports VCYCLE and WCYCLE types.
void Clear(int startLevel=0)
Clear impermanent data from previous setup.
void SetFuseProlongationAndUpdate(const bool &fuse)
bool GetFuseProlongationAndUpdate() const
bool GetImplicitTranspose() const
const RCP< const FactoryManagerBase > & GetLevelManager(const int levelID) const
Teuchos::ScalarTraits< SC > STS
MueLu::DefaultNode Node
std::string description_
cache description to avoid recreating in each call to description() - use ResetDescription() to force...
MagnitudeType rate_
Convergece rate.
static const NoFactory * get()
ConvergenceStatus Iterate(const MultiVector &B, MultiVector &X, ConvData conv=ConvData(), bool InitialGuessIsZero=false, LO startLevel=0)
Apply the multigrid preconditioner.
void ReplaceCoordinateMap(Level &level)
void AddKeepFlag(const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep)
Call Level::AddKeepFlag for each level of the Hierarchy.
void print(std::ostream &out=std::cout, const VerbLevel verbLevel=(MueLu::Parameters|MueLu::Statistics0)) const
Hierarchy::print is local hierarchy function, thus the statistics can be different from global ones...
Xpetra::UnderlyingLib lib()
double GetSmootherComplexity() const
Base class for factories (e.g., R, P, and A_coarse).
Print statistics that do not involve significant additional computation.
ConvData(std::pair< LO, MagnitudeType > p)
Xpetra::UnderlyingLib lib_
Epetra/Tpetra mode.
static CycleType GetDefaultCycle()
MueLu::DefaultScalar Scalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
MueLu::DefaultGlobalOrdinal GlobalOrdinal
void Write(const LO &start=-1, const LO &end=-1, const std::string &suffix="")
Print matrices in the multigrid hierarchy to file.
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
std::string description() const
Return a simple one-line description of this object.
void SetCycleStartLevel(int cycleStart)
Class that provides default factories within Needs class.
void ResetDescription()
force recreation of cached description_ next time description() is called:
void RemoveKeepFlag(const std::string &ename, const FactoryBase *factory, KeepType keep=MueLu::All)
Call Level::RemoveKeepFlag for each level of the Hierarchy.
void CheckLevel(Level &level, int levelID)
Helper function.
Always keep data, even accross run. This flag is set by Level::Keep(). This flag is propagated to coa...
static bool GetDefaultImplicitTranspose()
void AllocateLevelMultiVectors(int numvecs, bool forceMapCheck=false)
void SetPRrebalance(bool doPRrebalance)
size_t global_size_t
void SetMatvecParams(RCP< ParameterList > matvecParams)
int GetGlobalNumLevels() const
void SetPRViaCopyrebalance(bool doPRViaCopyrebalance)
static int GetDefaultMaxLevels()
CycleType GetCycle() const
Returns multigrid cycle type (supports VCYCLE and WCYCLE)
Base class for MueLu classes.
void SetImplicitTranspose(const bool &implicit)
CycleType Cycle_
V- or W-cycle.
Data struct for defining stopping criteria of multigrid iteration.
Array< RCP< MultiVector > > coarseRhs_
STS::magnitudeType MagnitudeType
void Delete(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Delete(ename, factory) for each level of the Hierarchy.
Array< RCP< Level > > Levels_
Container for Level objects.
bool isPreconditioner_
Hierarchy may be used in a standalone mode, or as a preconditioner.
static bool GetDefaultPRrebalance()
bool IsCalculationOfResidualRequired(const LO startLevel, const ConvData &conv) const
Decide if the residual needs to be computed.
MagnitudeType GetRate() const
static int GetDefaultCycleStartLevel()
void SetProlongatorScalingFactor(double scalingFactor)
Specify damping factor alpha such that x = x + alpha*P*c, where c is the coarse grid correction...
virtual ~Hierarchy()
Destructor.
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
double scalingFactor_
Scaling factor to be applied to coarse grid correction.
Array< RCP< MultiVector > > residual_
void Keep(const std::string &ename, const FactoryBase *factory=NoFactory::get())
Call Level::Keep(ename, factory) for each level of the Hierarchy.
ConvergenceStatus IsConverged(const Teuchos::Array< MagnitudeType > &residualNorm, const MagnitudeType convergenceTolerance) const
Decide if the multigrid iteration is converged.
bool isDumpingEnabled_
Graph dumping.
void setlib(Xpetra::UnderlyingLib inlib)
double GetOperatorComplexity() const
bool Setup(int coarseLevelID, const RCP< const FactoryManagerBase > fineLevelManager, const RCP< const FactoryManagerBase > coarseLevelManager, const RCP< const FactoryManagerBase > nextLevelManager=Teuchos::null)
Multi-level setup phase: build a new level of the hierarchy.
Array< RCP< MultiVector > > coarseX_
Array< RCP< MultiVector > > correction_
void EnableGraphDumping(const std::string &filename, int levelID=1)
ConvergenceStatus ComputeResidualAndPrintHistory(const Operator &A, const MultiVector &X, const MultiVector &B, const LO iteration, const LO startLevel, const ConvData &conv, MagnitudeType &previousResidualNorm)
Compute the residual norm and print it depending on the verbosity level.
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Array< RCP< const FactoryManagerBase > > levelManagers_
Level managers used during the Setup.
Xpetra::global_size_t maxCoarseSize_