MueLu
Version of the Day
|
Provides methods to build a multigrid hierarchy and apply multigrid cycles. More...
#include <MueLu_Hierarchy_decl.hpp>
Classes | |
struct | ConvData |
Data struct for defining stopping criteria of multigrid iteration. More... | |
Public Member Functions | |
int | LastLevelID () const |
void | AddLevel (const RCP< Level > &level) |
Add a level at the end of the hierarchy. More... | |
void | AddNewLevel () |
Add a new level at the end of the hierarchy. More... | |
RCP< Level > & | GetLevel (const int levelID=0) |
Retrieve a certain level from hierarchy. More... | |
int | GetNumLevels () const |
int | GetGlobalNumLevels () const |
MagnitudeType | GetRate () const |
double | GetOperatorComplexity () const |
double | GetSmootherComplexity () const |
void | CheckLevel (Level &level, int levelID) |
Helper function. More... | |
void | SetMatvecParams (RCP< ParameterList > matvecParams) |
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. More... | |
void | Setup (const FactoryManagerBase &manager=FactoryManager(), int startLevel=0, int numDesiredLevels=GetDefaultMaxLevels()) |
void | SetupRe () |
void | CheckForEmptySmoothersAndCoarseSolve () |
void | Clear (int startLevel=0) |
Clear impermanent data from previous setup. More... | |
void | ExpertClear () |
CycleType | GetCycle () const |
Returns multigrid cycle type (supports VCYCLE and WCYCLE) More... | |
void | SetCycle (CycleType Cycle) |
Supports VCYCLE and WCYCLE types. More... | |
void | SetCycleStartLevel (int cycleStart) |
void | SetProlongatorScalingFactor (double scalingFactor) |
Specify damping factor alpha such that x = x + alpha*P*c, where c is the coarse grid correction. More... | |
ConvergenceStatus | Iterate (const MultiVector &B, MultiVector &X, ConvData conv=ConvData(), bool InitialGuessIsZero=false, LO startLevel=0) |
Apply the multigrid preconditioner. More... | |
void | Write (const LO &start=-1, const LO &end=-1, const std::string &suffix="") |
Print matrices in the multigrid hierarchy to file. More... | |
void | EnableGraphDumping (const std::string &filename, int levelID=1) |
void | setlib (Xpetra::UnderlyingLib inlib) |
Xpetra::UnderlyingLib | lib () |
void | ResetDescription () |
force recreation of cached description_ next time description() is called: More... | |
void | AllocateLevelMultiVectors (int numvecs, bool forceMapCheck=false) |
void | DeleteLevelMultiVectors () |
Public Member Functions inherited from MueLu::BaseClass | |
virtual | ~BaseClass () |
Destructor. More... | |
Public Member Functions inherited from MueLu::VerboseObject | |
VerbLevel | GetVerbLevel () const |
Get the verbosity level. More... | |
void | SetVerbLevel (const VerbLevel verbLevel) |
Set the verbosity level of this object. More... | |
int | GetProcRankVerbose () const |
Get proc rank used for printing. Do not use this information for any other purpose. More... | |
int | SetProcRankVerbose (int procRank) const |
Set proc rank used for printing. More... | |
bool | IsPrint (MsgType type, int thisProcRankOnly=-1) const |
Find out whether we need to print out information for a specific message type. More... | |
Teuchos::FancyOStream & | GetOStream (MsgType type, int thisProcRankOnly=0) const |
Get an output stream for outputting the input message type. More... | |
Teuchos::FancyOStream & | GetBlackHole () const |
VerboseObject () | |
virtual | ~VerboseObject () |
Destructor. More... | |
Public Member Functions inherited from Teuchos::VerboseObject< VerboseObject > | |
VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObject & | setVerbLevel (const EVerbosityLevel verbLevel) const |
virtual const VerboseObject & | setOverridingVerbLevel (const EVerbosityLevel verbLevel) const |
virtual EVerbosityLevel | getVerbLevel () const |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT RCP< const ParameterList > | getValidVerboseObjectSublist () |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | setupVerboseObjectSublist (ParameterList *paramList) |
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT void | readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel) |
void | readVerboseObjectSublist (ParameterList *paramList, VerboseObject< ObjectType > *verboseObject) |
Public Member Functions inherited from Teuchos::VerboseObjectBase | |
virtual | ~VerboseObjectBase () |
VerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) | |
virtual const VerboseObjectBase & | setOStream (const RCP< FancyOStream > &oStream) const |
virtual const VerboseObjectBase & | setOverridingOStream (const RCP< FancyOStream > &oStream) const |
virtual VerboseObjectBase & | setLinePrefix (const std::string &linePrefix) |
virtual RCP< FancyOStream > | getOStream () const |
virtual RCP< FancyOStream > | getOverridingOStream () const |
virtual std::string | getLinePrefix () const |
virtual OSTab | getOSTab (const int tabs=1, const std::string &linePrefix="") const |
Public Member Functions inherited from MueLu::Describable | |
virtual | ~Describable () |
Destructor. More... | |
virtual std::string | ShortClassName () const |
Return the class name of the object, without template parameters and without namespace. More... | |
Public Member Functions inherited from Teuchos::Describable | |
void | describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const |
LabeledObject () | |
virtual | ~LabeledObject () |
virtual void | setObjectLabel (const std::string &objectLabel) |
virtual std::string | getObjectLabel () const |
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
Protected Member Functions | |
const RCP< const FactoryManagerBase > & | GetLevelManager (const int levelID) const |
Protected Member Functions inherited from Teuchos::VerboseObject< VerboseObject > | |
void | initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) |
Protected Member Functions inherited from Teuchos::VerboseObjectBase | |
void | initializeVerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) |
virtual void | informUpdatedVerbosityState () const |
Private Types | |
typedef Teuchos::ScalarTraits< SC > | STS |
typedef STS::magnitudeType | MagnitudeType |
Private Member Functions | |
void | DumpCurrentGraph (int level) const |
Hierarchy (const Hierarchy &h) | |
Copy constructor is not implemented. More... | |
bool | IsCalculationOfResidualRequired (const LO startLevel, const ConvData &conv) const |
Decide if the residual needs to be computed. More... | |
ConvergenceStatus | IsConverged (const Teuchos::Array< MagnitudeType > &residualNorm, const MagnitudeType convergenceTolerance) const |
Decide if the multigrid iteration is converged. More... | |
void | PrintResidualHistory (const LO iteration, const Teuchos::Array< MagnitudeType > &residualNorm) const |
Print residualNorm for this iteration to the screen. More... | |
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. More... | |
void | ReplaceCoordinateMap (Level &level) |
Private Attributes | |
Array< RCP< Level > > | Levels_ |
Container for Level objects. More... | |
Xpetra::global_size_t | maxCoarseSize_ |
bool | implicitTranspose_ |
bool | fuseProlongationAndUpdate_ |
bool | doPRrebalance_ |
bool | doPRViaCopyrebalance_ |
bool | isPreconditioner_ |
Hierarchy may be used in a standalone mode, or as a preconditioner. More... | |
CycleType | Cycle_ |
V- or W-cycle. More... | |
int | WCycleStartLevel_ |
Level at which to start W-cycle. More... | |
double | scalingFactor_ |
Scaling factor to be applied to coarse grid correction. More... | |
Xpetra::UnderlyingLib | lib_ |
Epetra/Tpetra mode. More... | |
std::string | description_ = "" |
cache description to avoid recreating in each call to description() - use ResetDescription() to force recreation in Setup, SetupRe, etc. More... | |
bool | isDumpingEnabled_ |
Graph dumping. More... | |
int | dumpLevel_ |
std::string | dumpFile_ |
MagnitudeType | rate_ |
Convergece rate. More... | |
Array< RCP< const FactoryManagerBase > > | levelManagers_ |
Level managers used during the Setup. More... | |
int | sizeOfAllocatedLevelMultiVectors_ |
Caching (Multi)Vectors used in Hierarchy::Iterate() More... | |
Array< RCP< MultiVector > > | residual_ |
Array< RCP< MultiVector > > | coarseRhs_ |
Array< RCP< MultiVector > > | coarseX_ |
Array< RCP< MultiVector > > | coarseImport_ |
Array< RCP< MultiVector > > | coarseExport_ |
Array< RCP< MultiVector > > | correction_ |
Friends | |
template<class S2 , class LO2 , class GO2 , class N2 > | |
class | Hierarchy |
Constructors/Destructors | |
Hierarchy () | |
Default constructor. More... | |
Hierarchy (const std::string &label) | |
Constructor that labels the hierarchy. More... | |
Hierarchy (const RCP< Matrix > &A) | |
Constructor. More... | |
Hierarchy (const RCP< Matrix > &A, const std::string &label) | |
Constructor. More... | |
virtual | ~Hierarchy () |
Destructor. More... | |
Set/Get Methods. | |
Xpetra::global_size_t | GetMaxCoarseSize () const |
bool | GetImplicitTranspose () const |
bool | GetFuseProlongationAndUpdate () const |
void | SetMaxCoarseSize (Xpetra::global_size_t maxCoarseSize) |
void | SetPRrebalance (bool doPRrebalance) |
void | SetPRViaCopyrebalance (bool doPRViaCopyrebalance) |
void | SetImplicitTranspose (const bool &implicit) |
void | SetFuseProlongationAndUpdate (const bool &fuse) |
static CycleType | GetDefaultCycle () |
static int | GetDefaultCycleStartLevel () |
static bool | GetDefaultImplicitTranspose () |
static bool | GetDefaultFuseProlongationAndUpdate () |
static Xpetra::global_size_t | GetDefaultMaxCoarseSize () |
static int | GetDefaultMaxLevels () |
static bool | GetDefaultPRrebalance () |
Permanent storage | |
void | Keep (const std::string &ename, const FactoryBase *factory=NoFactory::get()) |
Call Level::Keep(ename, factory) for each level of the Hierarchy. More... | |
void | Delete (const std::string &ename, const FactoryBase *factory=NoFactory::get()) |
Call Level::Delete(ename, factory) for each level of the Hierarchy. More... | |
void | AddKeepFlag (const std::string &ename, const FactoryBase *factory=NoFactory::get(), KeepType keep=MueLu::Keep) |
Call Level::AddKeepFlag for each level of the Hierarchy. More... | |
void | RemoveKeepFlag (const std::string &ename, const FactoryBase *factory, KeepType keep=MueLu::All) |
Call Level::RemoveKeepFlag for each level of the Hierarchy. More... | |
Overridden from Teuchos::Describable | |
std::string | description () const |
Return a simple one-line description of this object. More... | |
void | describe (Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const |
Print the Hierarchy with some verbosity level to a FancyOStream object. More... | |
void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::VERB_HIGH) const |
Print the object with some verbosity level to an FancyOStream object. More... | |
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. More... | |
void | IsPreconditioner (const bool flag) |
Additional Inherited Members | |
Static Public Member Functions inherited from MueLu::VerboseObject | |
static void | SetMueLuOStream (const Teuchos::RCP< Teuchos::FancyOStream > &mueluOStream) |
static void | SetMueLuOFileStream (const std::string &filename) |
static Teuchos::RCP < Teuchos::FancyOStream > | GetMueLuOStream () |
static void | SetDefaultVerbLevel (const VerbLevel defaultVerbLevel) |
Set the default (global) verbosity level. More... | |
static VerbLevel | GetDefaultVerbLevel () |
Get the default (global) verbosity level. More... | |
Static Public Member Functions inherited from Teuchos::VerboseObject< VerboseObject > | |
static void | setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel) |
static EVerbosityLevel | getDefaultVerbLevel () |
Static Public Member Functions inherited from Teuchos::VerboseObjectBase | |
static void | setDefaultOStream (const RCP< FancyOStream > &defaultOStream) |
static RCP< FancyOStream > | getDefaultOStream () |
Static Public Attributes inherited from Teuchos::Describable | |
static const EVerbosityLevel | verbLevel_default |
Provides methods to build a multigrid hierarchy and apply multigrid cycles.
Allows users to manually populate operators at different levels within a multigrid method and push them into the hierarchy via SetLevel() and/or to supply factories for automatically generating prolongators, restrictors, and coarse level discretizations. Additionally, this class contains an apply method that supports V and W cycles.
Definition at line 63 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 67 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 68 of file MueLu_Hierarchy_decl.hpp.
MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Hierarchy | ( | ) |
Default constructor.
Definition at line 41 of file MueLu_Hierarchy_def.hpp.
MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Hierarchy | ( | const std::string & | label | ) |
Constructor that labels the hierarchy.
Definition at line 60 of file MueLu_Hierarchy_def.hpp.
MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Hierarchy | ( | const RCP< Matrix > & | A | ) |
Constructor.
Definition at line 67 of file MueLu_Hierarchy_def.hpp.
MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Hierarchy | ( | const RCP< Matrix > & | A, |
const std::string & | label | ||
) |
Constructor.
Definition at line 90 of file MueLu_Hierarchy_def.hpp.
|
virtualdefault |
Destructor.
|
private |
Copy constructor is not implemented.
|
inlinestatic |
Definition at line 113 of file MueLu_Hierarchy_decl.hpp.
|
inlinestatic |
Definition at line 114 of file MueLu_Hierarchy_decl.hpp.
|
inlinestatic |
Definition at line 115 of file MueLu_Hierarchy_decl.hpp.
|
inlinestatic |
Definition at line 116 of file MueLu_Hierarchy_decl.hpp.
|
inlinestatic |
Definition at line 117 of file MueLu_Hierarchy_decl.hpp.
|
inlinestatic |
Definition at line 118 of file MueLu_Hierarchy_decl.hpp.
|
inlinestatic |
Definition at line 119 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 121 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 122 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 123 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 125 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 126 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 127 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 128 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 129 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 138 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 1398 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::AddLevel | ( | const RCP< Level > & | level | ) |
Add a level at the end of the hierarchy.
Definition at line 100 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::AddNewLevel | ( | ) |
Add a new level at the end of the hierarchy.
Definition at line 116 of file MueLu_Hierarchy_def.hpp.
RCP< Level > & MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetLevel | ( | const int | levelID = 0 | ) |
Retrieve a certain level from hierarchy.
Definition at line 123 of file MueLu_Hierarchy_def.hpp.
int MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetNumLevels | ( | ) | const |
Definition at line 130 of file MueLu_Hierarchy_def.hpp.
int MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetGlobalNumLevels | ( | ) | const |
Definition at line 135 of file MueLu_Hierarchy_def.hpp.
|
inline |
Definition at line 156 of file MueLu_Hierarchy_decl.hpp.
double MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetOperatorComplexity | ( | ) | const |
Definition at line 147 of file MueLu_Hierarchy_def.hpp.
double MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::GetSmootherComplexity | ( | ) | const |
Definition at line 170 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CheckLevel | ( | Level & | level, |
int | levelID | ||
) |
Helper function.
Definition at line 212 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetMatvecParams | ( | RCP< ParameterList > | matvecParams | ) |
Definition at line 222 of file MueLu_Hierarchy_def.hpp.
bool MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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.
This method is aimed to be used in a loop building the hierarchy level by level. See Hierarchy::Setup(manager, startLevel, numDesiredLevels) for an example of usage. @param coarseLevelID ID of the level to be built. @param fineLevelManager defines how to build missing data of the fineLevel (example: aggregates) @param coarseLevelManager defines how to build the level @param nextLevelManager defines how the next coarse level will be built. This is used to post corresponding request before building the coarse level to keep useful data.
CoarseLevel is considered to be the last level if:
Pre-condition: FineLevel:
Post-condition: FineLevel:
Definition at line 266 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Setup | ( | const FactoryManagerBase & | manager = FactoryManager() , |
int | startLevel = 0 , |
||
int | numDesiredLevels = GetDefaultMaxLevels() |
||
) |
Definition at line 575 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetupRe | ( | ) |
Definition at line 533 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CheckForEmptySmoothersAndCoarseSolve | ( | ) |
Definition at line 649 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Clear | ( | int | startLevel = 0 | ) |
Clear impermanent data from previous setup.
Definition at line 658 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ExpertClear | ( | ) |
Definition at line 667 of file MueLu_Hierarchy_def.hpp.
|
inline |
Returns multigrid cycle type (supports VCYCLE and WCYCLE)
Definition at line 222 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Supports VCYCLE and WCYCLE types.
Definition at line 225 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 227 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Specify damping factor alpha such that x = x + alpha*P*c, where c is the coarse grid correction.
Definition at line 230 of file MueLu_Hierarchy_decl.hpp.
ConvergenceStatus MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Iterate | ( | const MultiVector & | B, |
MultiVector & | X, | ||
ConvData | conv = ConvData() , |
||
bool | InitialGuessIsZero = false , |
||
LO | startLevel = 0 |
||
) |
Apply the multigrid preconditioner.
In theory, more general cycle types than just V- and W-cycles are possible. However, the enumerated type CycleType would have to be extended.
B | right-hand side of linear problem |
X | initial and final (approximate) solution of linear problem |
ConvData | struct which stores convergence criteria (maximum number of multigrid iterations or stopping tolerance) |
InitialGuessIsZero | Indicates whether the initial guess is zero |
startLevel | index of starting level to build multigrid hierarchy (default = 0) |
Definition at line 871 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Write | ( | const LO & | start = -1 , |
const LO & | end = -1 , |
||
const std::string & | suffix = "" |
||
) |
Print matrices in the multigrid hierarchy to file.
[in] | start | start level |
[in] | end | end level |
Default behavior is to print system and transfer matrices from the entire hierarchy. Files are named "A_0.m", "P_1.m", "R_1.m", etc, and are in matrix market coordinate format.
Definition at line 1141 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Keep | ( | const std::string & | ename, |
const FactoryBase * | factory = NoFactory::get() |
||
) |
Call Level::Keep(ename, factory) for each level of the Hierarchy.
Definition at line 1170 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Delete | ( | const std::string & | ename, |
const FactoryBase * | factory = NoFactory::get() |
||
) |
Call Level::Delete(ename, factory) for each level of the Hierarchy.
Definition at line 1176 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::AddKeepFlag | ( | const std::string & | ename, |
const FactoryBase * | factory = NoFactory::get() , |
||
KeepType | keep = MueLu::Keep |
||
) |
Call Level::AddKeepFlag for each level of the Hierarchy.
Definition at line 1182 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::RemoveKeepFlag | ( | const std::string & | ename, |
const FactoryBase * | factory, | ||
KeepType | keep = MueLu::All |
||
) |
Call Level::RemoveKeepFlag for each level of the Hierarchy.
Definition at line 1188 of file MueLu_Hierarchy_def.hpp.
|
virtual |
Return a simple one-line description of this object.
Reimplemented from MueLu::Describable.
Definition at line 1194 of file MueLu_Hierarchy_def.hpp.
|
virtual |
Print the Hierarchy with some verbosity level to a FancyOStream object.
[in] | out | The Teuchos::FancyOstream. |
[in] | verbLevel | Controls amount of output. |
Reimplemented from MueLu::Describable.
Definition at line 1210 of file MueLu_Hierarchy_def.hpp.
|
virtual |
Print the object with some verbosity level to an FancyOStream object.
Reimplemented from MueLu::Describable.
Definition at line 1205 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::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.
Definition at line 1386 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::IsPreconditioner | ( | const bool | flag | ) |
Indicate whether the multigrid method is a preconditioner or a solver.
This is used in conjunction with the verbosity level to determine whether the residuals can be printed.
Definition at line 1393 of file MueLu_Hierarchy_def.hpp.
|
inline |
Definition at line 302 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 308 of file MueLu_Hierarchy_decl.hpp.
|
inline |
Definition at line 309 of file MueLu_Hierarchy_decl.hpp.
|
inline |
force recreation of cached description_ next time description() is called:
Definition at line 312 of file MueLu_Hierarchy_decl.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::AllocateLevelMultiVectors | ( | int | numvecs, |
bool | forceMapCheck = false |
||
) |
Definition at line 1525 of file MueLu_Hierarchy_def.hpp.
void MueLu::Hierarchy< Scalar, LocalOrdinal, GlobalOrdinal, Node >::DeleteLevelMultiVectors | ( | ) |
Definition at line 1601 of file MueLu_Hierarchy_def.hpp.
|
inlineprotected |
Definition at line 320 of file MueLu_Hierarchy_decl.hpp.
|
private |
Decide if the residual needs to be computed.
Definition at line 1613 of file MueLu_Hierarchy_def.hpp.
|
private |
Decide if the multigrid iteration is converged.
We judge convergence by comparing the current residualNorm
to the user given convergenceTolerance
and then return the appropriate ConvergenceStatus
Definition at line 1619 of file MueLu_Hierarchy_def.hpp.
|
private |
Print residualNorm
for this iteration
to the screen.
Definition at line 1639 of file MueLu_Hierarchy_def.hpp.
|
private |
Compute the residual norm and print it depending on the verbosity level.
Definition at line 1650 of file MueLu_Hierarchy_def.hpp.
|
private |
We replace coordinates GIDs to make them consistent with matrix GIDs, even if user does not do that. Ideally, though, we should completely remove any notion of coordinate GIDs, and deal only with LIDs, assuming that they are consistent with matrix block IDs
Definition at line 1452 of file MueLu_Hierarchy_def.hpp.
|
friend |
Definition at line 136 of file MueLu_Hierarchy_decl.hpp.
|
private |
Container for Level objects.
Definition at line 351 of file MueLu_Hierarchy_decl.hpp.
|
private |
Minimum size of a matrix on any level. If we fall below that, we stop the coarsening
Definition at line 361 of file MueLu_Hierarchy_decl.hpp.
|
private |
Potential speed up of the setup by skipping R construction, and using transpose matrix-matrix product for RAP
Definition at line 365 of file MueLu_Hierarchy_decl.hpp.
|
private |
Potential speed up of the solve by fusing prolongation and update steps. This can lead to more iterations to round-off error accumulation.
Definition at line 369 of file MueLu_Hierarchy_decl.hpp.
|
private |
Potential speed up of the setup by skipping rebalancing of P and R, and doing extra import during solve
Definition at line 373 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 374 of file MueLu_Hierarchy_decl.hpp.
|
private |
Hierarchy may be used in a standalone mode, or as a preconditioner.
Definition at line 377 of file MueLu_Hierarchy_decl.hpp.
|
private |
V- or W-cycle.
Definition at line 380 of file MueLu_Hierarchy_decl.hpp.
|
private |
Level at which to start W-cycle.
Definition at line 383 of file MueLu_Hierarchy_decl.hpp.
|
private |
Scaling factor to be applied to coarse grid correction.
Definition at line 386 of file MueLu_Hierarchy_decl.hpp.
|
private |
Epetra/Tpetra mode.
Definition at line 389 of file MueLu_Hierarchy_decl.hpp.
|
mutableprivate |
cache description to avoid recreating in each call to description() - use ResetDescription() to force recreation in Setup, SetupRe, etc.
Definition at line 392 of file MueLu_Hierarchy_decl.hpp.
|
private |
Graph dumping.
If enabled, we dump the graph on a specified level into a specified file
Definition at line 399 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 401 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 402 of file MueLu_Hierarchy_decl.hpp.
|
private |
Convergece rate.
Definition at line 405 of file MueLu_Hierarchy_decl.hpp.
|
private |
Level managers used during the Setup.
Definition at line 408 of file MueLu_Hierarchy_decl.hpp.
|
private |
Caching (Multi)Vectors used in Hierarchy::Iterate()
Definition at line 411 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 412 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 412 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 412 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 412 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 412 of file MueLu_Hierarchy_decl.hpp.
|
private |
Definition at line 412 of file MueLu_Hierarchy_decl.hpp.