MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node > Class Template Reference

#include <MueLu_ParameterListInterpreter_decl.hpp>

Inheritance diagram for MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >:
MueLu::HierarchyManager< Scalar, LocalOrdinal, GlobalOrdinal, Node > MueLu::HierarchyFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > MueLu::BaseClass MueLu::VerboseObject MueLu::Describable Teuchos::VerboseObject< VerboseObject > Teuchos::Describable Teuchos::VerboseObjectBase Teuchos::LabeledObject

Public Member Functions

void SetParameterList (const Teuchos::ParameterList &paramList)
 Set parameter list for Parameter list interpreter. More...
 
void SetupHierarchy (Hierarchy &H) const
 Call the SetupHierarchy routine from the HiearchyManager object. More...
 
- Public Member Functions inherited from MueLu::HierarchyManager< Scalar, LocalOrdinal, GlobalOrdinal, Node >
 HierarchyManager (int numDesiredLevel=MasterList::getDefault< int >("max levels"))
 
virtual ~HierarchyManager ()
 
void AddFactoryManager (int startLevel, int numDesiredLevel, RCP< FactoryManagerBase > manager)
 
RCP< FactoryManagerBaseGetFactoryManager (int levelID) const
 
size_t getNumFactoryManagers () const
 returns number of factory managers stored in levelManagers_ vector. More...
 
void CheckConfig ()
 
virtual RCP< HierarchyCreateHierarchy () const
 Create an empty Hierarchy object. More...
 
virtual RCP< HierarchyCreateHierarchy (const std::string &label) const
 Create a labeled empty Hierarchy object. More...
 
- Public Member Functions inherited from MueLu::HierarchyFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node >
virtual ~HierarchyFactory ()
 Destructor. More...
 
- 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::FancyOStreamGetOStream (MsgType type, int thisProcRankOnly=0) const
 Get an output stream for outputting the input message type. More...
 
Teuchos::FancyOStreamGetBlackHole () 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 VerboseObjectsetVerbLevel (const EVerbosityLevel verbLevel) const
 
virtual const VerboseObjectsetOverridingVerbLevel (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 VerboseObjectBasesetOStream (const RCP< FancyOStream > &oStream) const
 
virtual const VerboseObjectBasesetOverridingOStream (const RCP< FancyOStream > &oStream) const
 
virtual VerboseObjectBasesetLinePrefix (const std::string &linePrefix)
 
virtual RCP< FancyOStreamgetOStream () const
 
virtual RCP< FancyOStreamgetOverridingOStream () 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...
 
virtual void describe (Teuchos::FancyOStream &out_arg, const VerbLevel verbLevel=Default) const
 
virtual std::string description () const
 Return a simple one-line description of this object. More...
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
 Print the object with some verbosity level to an FancyOStream object. 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)
 

Private Types

typedef std::pair< std::string,
const FactoryBase * > 
keep_pair
 

Private Member Functions

virtual void SetupOperator (Operator &A) const
 Setup Operator object. More...
 

Private Attributes

int blockSize_
 block size of matrix (fixed block size) More...
 
CycleType Cycle_
 multigrid cycle type (V-cycle or W-cycle) More...
 
int WCycleStartLevel_
 in case of W-cycle, level on which cycle should start More...
 
double scalingFactor_
 prolongator scaling factor More...
 
GlobalOrdinal dofOffset_
 global offset variable describing offset of DOFs in operator More...
 

Constructors/Destructors

 ParameterListInterpreter ()
 Empty constructor. More...
 
 ParameterListInterpreter (Teuchos::ParameterList &paramList, Teuchos::RCP< const Teuchos::Comm< int > > comm=Teuchos::null, Teuchos::RCP< FactoryFactory > factFact=Teuchos::null, Teuchos::RCP< FacadeClassFactory > facadeFact=Teuchos::null)
 Constructor that accepts a user-provided ParameterList. More...
 
 ParameterListInterpreter (const std::string &xmlFileName, const Teuchos::Comm< int > &comm, Teuchos::RCP< FactoryFactory > factFact=Teuchos::null, Teuchos::RCP< FacadeClassFactory > facadeFact=Teuchos::null)
 Constructor that reads parameters from an XML file. More...
 
virtual ~ParameterListInterpreter ()
 Destructor. More...
 
bool changedPRrebalance_
 Easy interpreter stuff. More...
 
bool changedImplicitTranspose_
 
bool useCoordinates_
 
bool useKokkos_
 
void SetEasyParameterList (const Teuchos::ParameterList &paramList)
 
void Validate (const Teuchos::ParameterList &paramList) const
 
void UpdateFactoryManager (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_Smoothers (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_CoarseSolvers (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_Aggregation_TentativeP (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_Restriction (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_RAP (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_Coordinates (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_Repartition (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
 
void UpdateFactoryManager_Nullspace (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps, RCP< Factory > &nullSpaceFactory) const
 
void UpdateFactoryManager_SemiCoarsen (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_PCoarsen (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_SA (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_Emin (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_PG (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
void UpdateFactoryManager_Matlab (Teuchos::ParameterList &paramList, const Teuchos::ParameterList &defaultList, FactoryManager &manager, int levelID, std::vector< keep_pair > &keeps) const
 
typedef std::map< std::string,
RCP< const FactoryBase > > 
FactoryMap
 
typedef std::map< std::string,
RCP< FactoryManagerBase > > 
FactoryManagerMap
 
Teuchos::RCP< FactoryFactoryfactFact_
 Internal factory for factories. More...
 
Teuchos::RCP
< MueLu::FacadeClassFactory
< Scalar, LocalOrdinal,
GlobalOrdinal, Node > > 
facadeFact_
 FacadeClass factory. More...
 
void SetFactoryParameterList (const Teuchos::ParameterList &paramList)
 Factory interpreter stuff. More...
 
void BuildFactoryMap (const Teuchos::ParameterList &paramList, const FactoryMap &factoryMapIn, FactoryMap &factoryMapOut, FactoryManagerMap &factoryManagers) const
 Interpret "Factories" sublist. More...
 

Additional Inherited Members

- Public Types inherited from MueLu::HierarchyManager< Scalar, LocalOrdinal, GlobalOrdinal, Node >
typedef std::map< std::string,
RCP< const FactoryBase > > 
FactoryMap
 
- 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< FancyOStreamgetDefaultOStream ()
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 
- Protected Member Functions inherited from MueLu::HierarchyManager< Scalar, LocalOrdinal, GlobalOrdinal, Node >
virtual void SetupExtra (Hierarchy &) const
 Setup extra data. More...
 
Teuchos::RCP< FactoryManagerBaseLvlMngr (int levelID, int lastLevelID) 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
 
- Protected Attributes inherited from MueLu::HierarchyManager< Scalar, LocalOrdinal, GlobalOrdinal, Node >
int numDesiredLevel_
 
Xpetra::global_size_t maxCoarseSize_
 
MsgType verbosity_
 
bool doPRrebalance_
 
bool implicitTranspose_
 
bool fuseProlongationAndUpdate_
 
int sizeOfMultiVectors_
 
int graphOutputLevel_
 
Teuchos::Array< int > matricesToPrint_
 
Teuchos::Array< int > prolongatorsToPrint_
 
Teuchos::Array< int > restrictorsToPrint_
 
Teuchos::Array< int > nullspaceToPrint_
 
Teuchos::Array< int > coordinatesToPrint_
 
Teuchos::Array< int > elementToNodeMapsToPrint_
 
Teuchos::RCP
< Teuchos::ParameterList
matvecParams_
 
std::map< int, std::vector
< keep_pair > > 
keep_
 

Detailed Description

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
class MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >

Definition at line 118 of file MueLu_ParameterListInterpreter_decl.hpp.

Member Typedef Documentation

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
typedef std::pair<std::string, const FactoryBase*> MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::keep_pair
private

Definition at line 122 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
typedef std::map<std::string, RCP<const FactoryBase> > MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::FactoryMap
private

Definition at line 260 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
typedef std::map<std::string, RCP<FactoryManagerBase> > MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::FactoryManagerMap
private

Definition at line 261 of file MueLu_ParameterListInterpreter_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ParameterListInterpreter ( )
inlineprotected

Empty constructor.

Constructor for derived classes

Definition at line 133 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ParameterListInterpreter ( Teuchos::ParameterList paramList,
Teuchos::RCP< const Teuchos::Comm< int > >  comm = Teuchos::null,
Teuchos::RCP< FactoryFactory factFact = Teuchos::null,
Teuchos::RCP< FacadeClassFactory facadeFact = Teuchos::null 
)

Constructor that accepts a user-provided ParameterList.

Constructor for parameter list interpreter which directly interprets Teuchos::ParameterLists

The parameter list can be either in the easy parameter list format or in the factory driven parameter list format.

Parameters
[in]paramList(Teuchos::ParameterList): ParameterList containing the MueLu parameters
[in]comm(RCP<Teuchos::Comm<int> >): Optional RCP of a Teuchos communicator (default: Teuchos::null)
[in]factFact(RCP<FactoryFactory>): Optional parameter allowing to define user-specific factory interpreters for user-specific extensions of the XML interface. (default: Teuchos::null)
[in]facadeFact(RCP<FacadeFactory>): Optional parameter containing a FacadeFactory class. The user can register its own facade classes in the FacadeFactory and provide it to the ParameterListInterpreter. (default: Teuchos::null, means, only standard FacadeClass that come with MueLu are available)

Definition at line 131 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::ParameterListInterpreter ( const std::string &  xmlFileName,
const Teuchos::Comm< int > &  comm,
Teuchos::RCP< FactoryFactory factFact = Teuchos::null,
Teuchos::RCP< FacadeClassFactory facadeFact = Teuchos::null 
)

Constructor that reads parameters from an XML file.

XML options are converted to ParameterList entries by Teuchos.

Parameters
[in]xmlFileName(std::string): XML file to read
[in]comm(Teuchos::Comm<int>): Teuchos communicator
[in]factFact(RCP<FactoryFactory>): Optional parameter allowing to define user-specific factory interpreters for user-specific extensions of the XML interface. (default: Teuchos::null)
[in]facadeFact(RCP<FacadeFactory>): Optional parameter containing a FacadeFactory class. The user can register its own facade classes in the FacadeFactory and provide it to the ParameterListInterpreter. (default: Teuchos::null, means, only standard FacadeClass that come with MueLu are available)

Definition at line 157 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
virtual MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::~ParameterListInterpreter ( )
inlinevirtual

Destructor.

Definition at line 166 of file MueLu_ParameterListInterpreter_decl.hpp.

Member Function Documentation

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetParameterList ( const Teuchos::ParameterList paramList)

Set parameter list for Parameter list interpreter.

The routine checks whether it is a parameter list in the easy parameter format or the more advanced factory-based parameter format and calls the corresponding interpreter routine.

When finished, the parameter list is set that will used by the hierarchy build phase.

This method includes validation and some pre-parsing of the list for:

  • verbosity level
  • data to export
  • cycle type
  • max coarse size
  • max levels
  • number of equations
Parameters
[in]paramList,:ParameterList containing the MueLu parameters.

Definition at line 170 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetupHierarchy ( Hierarchy H) const
virtual

Call the SetupHierarchy routine from the HiearchyManager object.

Reimplemented from MueLu::HierarchyManager< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 2322 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetupOperator ( Operator &  A) const
privatevirtual

Setup Operator object.

Reimplemented from MueLu::HierarchyManager< Scalar, LocalOrdinal, GlobalOrdinal, Node >.

Definition at line 2302 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetEasyParameterList ( const Teuchos::ParameterList paramList)
private

Definition at line 248 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Validate ( const Teuchos::ParameterList paramList) const
private

Definition at line 1826 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 519 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Smoothers ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 657 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_CoarseSolvers ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 867 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Aggregation_TentativeP ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 926 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Restriction ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1255 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_RAP ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1103 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Coordinates ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1222 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Repartition ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps,
RCP< Factory > &  nullSpaceFactory 
) const
private

Definition at line 1318 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Nullspace ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps,
RCP< Factory > &  nullSpaceFactory 
) const
private

Definition at line 1524 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_SemiCoarsen ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1553 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_PCoarsen ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1613 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_SA ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1665 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Emin ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1715 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_PG ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1783 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::UpdateFactoryManager_Matlab ( Teuchos::ParameterList paramList,
const Teuchos::ParameterList defaultList,
FactoryManager manager,
int  levelID,
std::vector< keep_pair > &  keeps 
) const
private

Definition at line 1803 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::SetFactoryParameterList ( const Teuchos::ParameterList paramList)
private

Factory interpreter stuff.

Definition at line 1905 of file MueLu_ParameterListInterpreter_def.hpp.

template<class Scalar , class LocalOrdinal , class GlobalOrdinal , class Node >
void MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::BuildFactoryMap ( const Teuchos::ParameterList paramList,
const FactoryMap factoryMapIn,
FactoryMap factoryMapOut,
FactoryManagerMap factoryManagers 
) const
private

Interpret "Factories" sublist.

Parameters
paramList[in]: "Factories" ParameterList
factoryMapIn[in]: FactoryMap maps variable names to factories. This factory map is used to resolve data dependencies of previously defined factories.
factoryMapOut[out]: FactoryMap maps variable names to factories. New factory entries are added to that FactoryMap. Usually, factoryMapIn and factoryMapOut should use the same object, such that new factories are added. We have to distinguish input and output if we build sub-factory managers, though.
factoryManagers[in/out]: FacotryManagerMap maps group names to a FactoryManager object.

Interpret "Factories" parameter list. For each "factory" entry, add a new entry in the factoryMapOut map or create a new FacotryManager

Parameter List Parsing:

Create an entry in factoryMapOut for each parameter of the list paramList

<ParameterList name="..."> <Parameter name="smootherFact0" type="string" value="TrilinosSmoother">

<ParameterList name="smootherFact1"> <Parameter name="type" type="string" value="TrilinosSmoother"> ... </ParameterList> </ParameterList>


Group factories We can group factories using parameter sublists with the "group" parameter

<ParameterList name="myFirstGroup"> <Parameter name="group" type="string" value="FactoryManager"> <Parameter name="A" type="string" value="mySubBlockAFactory1"> <Parameter name="P" type="string" value="myTentativePFact1"> <Parameter name="Aggregates" type="string" value="myAggFact1"> <Parameter name="Nullspace" type="string" value="myNspFact1"> <Parameter name="CoarseMap" type="string" value="myCoarseMap1"> </ParameterList> <ParameterList name="mySecondGroup"> <Parameter name="group" type="string" value="FactoryManager"> <Parameter name="A" type="string" value="mySubBlockAFactory2"> <Parameter name="P" type="string" value="myTentativePFact2"> <Parameter name="Aggregates" type="string" value="myAggFact1"> <Parameter name="Nullspace" type="string" value="myNspFact2"> <Parameter name="CoarseMap" type="string" value="myCoarseMap2"> </ParameterList>

These factory groups can be used with factories for blocked operators (such as the BlockedPFactory) to easily define the operations on the sub-blocks.

<ParameterList name="myBlockedPFact"> <Parameter name="factory" type="string" value="BlockedPFactory">

<ParameterList name="block1"> <Parameter name="group" type="string" value="myFirstGroup"> </ParameterList>

<ParameterList name="block2"> <Parameter name="group" type="string" value="mySecondGroup"> </ParameterList> </ParameterList>

As an alternative one can also directly specify the factories in the sublists "block1", "block2", etc..., of course. But using blocks has the advantage that one can reuse them in all blocked factories.

<ParameterList name="myBlockedPFact"> <Parameter name="factory" type="string" value="BlockedPFactory">

<ParameterList name="block1"> <Parameter name="A" type="string" value="mySubBlockAFactory1"> <Parameter name="P" type="string" value="myTentativePFact1"> <Parameter name="Aggregates" type="string" value="myAggFact1"> <Parameter name="Nullspace" type="string" value="myNspFact1"> <Parameter name="CoarseMap" type="string" value="myCoarseMap1"> </ParameterList>

<ParameterList name="block2"> <Parameter name="A" type="string" value="mySubBlockAFactory2"> <Parameter name="P" type="string" value="myTentativePFact2"> <Parameter name="Aggregates" type="string" value="myAggFact1"> <Parameter name="Nullspace" type="string" value="myNspFact2"> <Parameter name="CoarseMap" type="string" value="myCoarseMap2"> </ParameterList> </ParameterList>

As an alternative one can also directly specify the factories in the sublists "block1", "block2", etc..., of course.------— add more dependencies (circular dependencies)

The NullspaceFactory needs to know which factory generates the null space on the coarse level (e.g., the TentativePFactory or the RebalancedPFactory). However, we cannot set the information in this place in the xml file, since the tentative prolongator facotry is typically defined later. We have to add that dependency later to the NullspaceFactory:

<ParameterList name="myNspFact"> <Parameter name="factory" type="string" value="NullspaceFactory">

</ParameterList>

<ParameterList name="myTentativePFact"> <Parameter name="factory" type="string" value="TentativePFactory"> <...> <Parameter name="Nullspace" type="string" value="myNspFact"> <Parameter name="CoarseMap" type="string" value="myCoarseMap"> </ParameterList>

<ParameterList name="myRebalanceProlongatorFact"> <Parameter name="factory" type="string" value="RebalanceTransferFactory"> <...> <Parameter name="Nullspace" type="string" value="myTentativePFact"> </ParameterList>

After the definition of the generating factory for the nullspace (in this case myRebalanceProlongatorFact) we add that dependency to the NullspaceFactory instance myNspFact

<ParameterList name="myNspFactDeps"> <Parameter name="dependency for" type="string" value="myNspFact"> <Parameter name="Nullspace" type="string" value="myRebalanceProlongatorFact"> </ParameterList>

We have to create a new block (with a different name than myNspFact). In the example we use "myNspFactDeps". It should contain a parameter "dependency for" with the name of the factory that we want the dependencies to be addded to. With above block we do not need the entry for the Nullspace in the global FactoryManager any more.

Definition at line 2199 of file MueLu_ParameterListInterpreter_def.hpp.

Member Data Documentation

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
int MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::blockSize_
private

block size of matrix (fixed block size)

Definition at line 195 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
CycleType MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::Cycle_
private

multigrid cycle type (V-cycle or W-cycle)

Definition at line 196 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
int MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::WCycleStartLevel_
private

in case of W-cycle, level on which cycle should start

Definition at line 197 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
double MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::scalingFactor_
private

prolongator scaling factor

Definition at line 198 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
GlobalOrdinal MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::dofOffset_
private

global offset variable describing offset of DOFs in operator

Definition at line 199 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
bool MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::changedPRrebalance_
private

Easy interpreter stuff.

Definition at line 204 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
bool MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::changedImplicitTranspose_
private

Definition at line 205 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
bool MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::useCoordinates_
private

Definition at line 245 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
bool MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::useKokkos_
private

Definition at line 246 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
Teuchos::RCP<FactoryFactory> MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::factFact_
private

Internal factory for factories.

Definition at line 266 of file MueLu_ParameterListInterpreter_decl.hpp.

template<class Scalar = DefaultScalar, class LocalOrdinal = DefaultLocalOrdinal, class GlobalOrdinal = DefaultGlobalOrdinal, class Node = DefaultNode>
Teuchos::RCP<MueLu::FacadeClassFactory<Scalar, LocalOrdinal, GlobalOrdinal, Node> > MueLu::ParameterListInterpreter< Scalar, LocalOrdinal, GlobalOrdinal, Node >::facadeFact_
private

FacadeClass factory.

Definition at line 269 of file MueLu_ParameterListInterpreter_decl.hpp.


The documentation for this class was generated from the following files: