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 // ***********************************************************************
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_HIERARCHY_DECL_HPP
47 #define MUELU_HIERARCHY_DECL_HPP
48 
50 #include <Teuchos_Ptr.hpp>
51 
52 #include <Xpetra_ConfigDefs.hpp> // global_size_t
53 #include <Xpetra_Matrix_fwd.hpp>
54 #include <Xpetra_MultiVector_fwd.hpp>
55 #include <Xpetra_MultiVectorFactory_fwd.hpp>
56 #include <Xpetra_Operator_fwd.hpp>
57 
58 #include <MueLu_SmootherCloner.hpp>
59 #include "MueLu_ConfigDefs.hpp"
60 #include "MueLu_BaseClass.hpp"
61 #include "MueLu_Hierarchy_fwd.hpp"
62 
63 #include "MueLu_Types.hpp"
64 
66 #include "MueLu_FactoryManager.hpp" // no fwd declaration because constructor of FactoryManager is used as a default parameter of Setup()
68 #include "MueLu_KeepType.hpp"
69 #include "MueLu_Level_fwd.hpp"
70 #include "MueLu_MasterList.hpp"
71 #include "MueLu_NoFactory.hpp"
72 #include "MueLu_PerfUtils_fwd.hpp"
73 #include "MueLu_PFactory_fwd.hpp"
74 #include "MueLu_RFactory_fwd.hpp"
79 #include "MueLu_Utilities_fwd.hpp"
80 
81 namespace MueLu {
82 
83  enum ReturnType {
87  };
88 
99  template <class Scalar = DefaultScalar,
102  class Node = DefaultNode>
103  class Hierarchy : public BaseClass {
104 #undef MUELU_HIERARCHY_SHORT
105 #include "MueLu_UseShortNames.hpp"
106 
109 
111  struct ConvData {
112  ConvData() : maxIts_(1), tol_(-STS::magnitude(STS::one())) { }
113  ConvData(LO maxIts) : maxIts_(maxIts), tol_(-STS::magnitude(STS::one())) { }
114  ConvData(MagnitudeType tol) : maxIts_(10000), tol_(tol) { }
115  ConvData(std::pair<LO,MagnitudeType> p) : maxIts_(p.first), tol_(p.second) { }
116 
119  };
120 
121  public:
122 
124 
125 
127  Hierarchy();
129  Hierarchy(const std::string& label);
130 
132  Hierarchy(const RCP<Matrix> & A);
133 
135  Hierarchy(const RCP<Matrix> & A, const std::string& label);
136 
138  virtual ~Hierarchy() { }
139 
141 
143 
144 
146  static CycleType GetDefaultCycle() { return MasterList::getDefault<std::string>("cycle type") == "V" ? VCYCLE : WCYCLE; }
147  static int GetDefaultCycleStartLevel() { return MasterList::getDefault<int>("W cycle start level"); }
148  static bool GetDefaultImplicitTranspose() { return MasterList::getDefault<bool>("transpose: use implicit"); }
149  static bool GetDefaultFuseProlongationAndUpdate() { return MasterList::getDefault<bool>("fuse prolongation and update"); }
150  static Xpetra::global_size_t GetDefaultMaxCoarseSize() { return MasterList::getDefault<int>("coarse: max size"); }
151  static int GetDefaultMaxLevels() { return MasterList::getDefault<int>("max levels"); }
152  static bool GetDefaultPRrebalance() { return MasterList::getDefault<bool>("repartition: rebalance P and R"); }
153 
154  Xpetra::global_size_t GetMaxCoarseSize() const { return maxCoarseSize_; }
155  bool GetImplicitTranspose() const { return implicitTranspose_; }
157 
158  void SetMaxCoarseSize(Xpetra::global_size_t maxCoarseSize) { maxCoarseSize_ = maxCoarseSize; }
159  void SetPRrebalance(bool doPRrebalance) { doPRrebalance_ = doPRrebalance; }
160  void SetImplicitTranspose(const bool& implicit) { implicitTranspose_ = implicit; }
161  void SetFuseProlongationAndUpdate(const bool& fuse) { fuseProlongationAndUpdate_ = fuse; }
162 
164 
166 
167  template<class S2, class LO2, class GO2, class N2>
168  friend class Hierarchy;
169 
170  private:
171  int LastLevelID() const { return Levels_.size() - 1; }
172  void DumpCurrentGraph() const;
173 
174  public:
175 
177  void AddLevel(const RCP<Level> & level);
178 
180  void AddNewLevel();
181 
183  RCP<Level> & GetLevel(const int levelID = 0);
184 
185  int GetNumLevels() const;
186  int GetGlobalNumLevels() const;
187 
188  MagnitudeType GetRate() const { return rate_; }
189 
190  // This function is global
191  double GetOperatorComplexity() const;
192 
193  // This function is global
194  double GetSmootherComplexity() const;
195 
197  void CheckLevel(Level& level, int levelID);
198 
199  void SetMatvecParams(RCP<ParameterList> matvecParams);
200 
202 
239  bool Setup(int coarseLevelID, const RCP<const FactoryManagerBase> fineLevelManager /* = Teuchos::null */, const RCP<const FactoryManagerBase> coarseLevelManager,
240  const RCP<const FactoryManagerBase> nextLevelManager = Teuchos::null);
241 
243  void Setup(const FactoryManagerBase& manager = FactoryManager(), int startLevel = 0, int numDesiredLevels = GetDefaultMaxLevels());
244 
245  void SetupRe();
246 
248  void Clear(int startLevel = 0);
249  void ExpertClear();
250 
252  CycleType GetCycle() const { return Cycle_; }
253 
255  void SetCycle(CycleType Cycle) { Cycle_ = Cycle; }
256 
257  void SetCycleStartLevel(int cycleStart) { WCycleStartLevel_ = cycleStart; }
258 
260  void SetProlongatorScalingFactor(double scalingFactor) { scalingFactor_ = scalingFactor; }
261 
274  ReturnType Iterate(const MultiVector& B, MultiVector& X, ConvData conv = ConvData(),
275  bool InitialGuessIsZero = false, LO startLevel = 0);
276 
286  void Write(const LO &start=-1, const LO &end=-1, const std::string &suffix="");
287 
289 
291 
292 
294  void Keep(const std::string & ename, const FactoryBase* factory = NoFactory::get());
295 
297  void Delete(const std::string& ename, const FactoryBase* factory = NoFactory::get());
298 
300  void AddKeepFlag(const std::string & ename, const FactoryBase* factory = NoFactory::get(), KeepType keep = MueLu::Keep);
301 
303  void RemoveKeepFlag(const std::string & ename, const FactoryBase* factory, KeepType keep = MueLu::All);
304 
306 
308 
309 
311  std::string description() const;
312 
318  void describe(Teuchos::FancyOStream& out, const VerbLevel verbLevel = Default) const;
320 
322  void print(std::ostream& out = std::cout, const VerbLevel verbLevel = (MueLu::Parameters | MueLu::Statistics0)) const;
323 
328  void IsPreconditioner(const bool flag);
329 
331 
332  void EnableGraphDumping(const std::string& filename, int levelID = 1) {
333  isDumpingEnabled_ = true;
334  dumpLevel_ = levelID;
335  dumpFile_ = filename;
336  }
337 
338  void setlib(Xpetra::UnderlyingLib inlib) { lib_ = inlib; }
339  Xpetra::UnderlyingLib lib() { return lib_; }
340 
343  description_ = "";
344  }
345 
346  void AllocateLevelMultiVectors(int numvecs);
348 
349  protected:
350  const RCP<const FactoryManagerBase>& GetLevelManager(const int levelID) const {
351  return levelManagers_[levelID];
352  }
353 
354  private:
356  Hierarchy(const Hierarchy &h);
357 
360 
365  void ReplaceCoordinateMap(Level& level);
366 
369  Xpetra::global_size_t maxCoarseSize_;
370 
374 
378 
382 
385 
388 
391 
394 
396  Xpetra::UnderlyingLib lib_;
397 
399  mutable std::string description_ = ""; // mutable so that we can lazily initialize in description(), which is declared const
400 
408  std::string dumpFile_;
409 
412 
415 
419 
420 
421  }; //class Hierarchy
422 
423 } //namespace MueLu
424 
425 #define MUELU_HIERARCHY_SHORT
426 #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
KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
RCP< Level > & GetLevel(const int levelID=0)
Retrieve a certain level from hierarchy.
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.
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)
Xpetra::global_size_t GetMaxCoarseSize() const
void AddNewLevel()
Add a new level at the end of the hierarchy.
static Xpetra::global_size_t GetDefaultMaxCoarseSize()
void DumpCurrentGraph() const
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()
void ReplaceCoordinateMap(Level &level)
void AllocateLevelMultiVectors(int numvecs)
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
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:99
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 SetPRrebalance(bool doPRrebalance)
void SetMatvecParams(RCP< ParameterList > matvecParams)
int GetGlobalNumLevels() const
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.
virtual ~Hierarchy()
Destructor.
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()
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...
void describe(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the Hierarchy with some verbosity level to a FancyOStream object.
ReturnType Iterate(const MultiVector &B, MultiVector &X, ConvData conv=ConvData(), bool InitialGuessIsZero=false, LO startLevel=0)
Apply the multigrid preconditioner.
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.
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)
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_