10 #ifndef MUELU_HIERARCHY_DECL_HPP
11 #define MUELU_HIERARCHY_DECL_HPP
14 #include <Teuchos_Ptr.hpp>
30 #include "MueLu_FactoryManager.hpp"
64 #undef MUELU_HIERARCHY_SHORT
135 template <
class S2,
class LO2,
class GO2,
class N2>
218 void Clear(
int startLevel = 0);
245 bool InitialGuessIsZero =
false,
LO startLevel = 0);
256 void Write(
const LO& start = -1,
const LO& end = -1,
const std::string& suffix =
"");
347 const MultiVector& B,
const LO iteration,
348 const LO startLevel,
const ConvData& conv,
MagnitudeType& previousResidualNorm);
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.
void DeleteLevelMultiVectors()
Array< RCP< MultiVector > > coarseImport_
static bool GetDefaultFuseProlongationAndUpdate()
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.
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.
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
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.
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...
ConvData(MagnitudeType tol)
static bool GetDefaultImplicitTranspose()
void AllocateLevelMultiVectors(int numvecs, bool forceMapCheck=false)
void SetPRrebalance(bool doPRrebalance)
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.
bool fuseProlongationAndUpdate_
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)
bool doPRViaCopyrebalance_
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_