MueLu
Version of the Day
|
Factory for creating a graph based on a given matrix. More...
#include <MueLu_CoalesceDropFactory_decl.hpp>
Public Member Functions | |
void | Build (Level ¤tLevel) const |
Build an object with this factory. More... | |
![]() | |
SingleLevelFactoryBase () | |
Constructor. More... | |
virtual | ~SingleLevelFactoryBase () |
Destructor. More... | |
virtual void | CallBuild (Level &requestedLevel) const |
virtual void | CallDeclareInput (Level &requestedLevel) const |
![]() | |
void | EnableMultipleCallCheck () const |
void | DisableMultipleCallCheck () const |
void | ResetDebugData () const |
Factory () | |
Constructor. More... | |
virtual | ~Factory () |
Destructor. More... | |
virtual void | SetFactory (const std::string &varName, const RCP< const FactoryBase > &factory) |
Configuration. More... | |
const RCP< const FactoryBase > | GetFactory (const std::string &varName) const |
Default implementation of FactoryAcceptor::GetFactory() More... | |
RCP< ParameterList > | RemoveFactoriesFromList (const ParameterList &list) const |
![]() | |
FactoryBase () | |
Constructor. More... | |
virtual | ~FactoryBase () |
Destructor. More... | |
int | GetID () const |
return unique factory id More... | |
![]() | |
virtual | ~BaseClass () |
Destructor. More... | |
![]() | |
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... | |
![]() | |
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) |
![]() | |
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 |
![]() | |
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... | |
![]() | |
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) |
![]() | |
virtual | ~FactoryAcceptor () |
![]() | |
ParameterListAcceptorImpl () | |
virtual | ~ParameterListAcceptorImpl ()=default |
virtual void | SetParameterList (const Teuchos::ParameterList ¶mList) |
Set parameters from a parameter list and return with default values. More... | |
virtual const Teuchos::ParameterList & | GetParameterList () const |
void | SetParameter (const std::string &name, const ParameterEntry &entry) |
Set a parameter directly as a ParameterEntry. More... | |
const ParameterEntry & | GetParameter (const std::string &name) const |
Retrieves a const entry with the name name. More... | |
virtual void | GetDocumentation (std::ostream &os) const |
![]() | |
ParameterListAcceptor () | |
virtual | ~ParameterListAcceptor ()=default |
Private Member Functions | |
void | MergeRows (const Matrix &A, const LO row, Array< LO > &cols, const Array< LO > &translation) const |
Method to merge rows of matrix for systems of PDEs. More... | |
void | MergeRowsWithDropping (const Matrix &A, const LO row, const ArrayRCP< const SC > &ghostedDiagVals, SC threshold, Array< LO > &cols, const Array< LO > &translation) const |
Teuchos::RCP< Xpetra::Matrix < Scalar, LocalOrdinal, GlobalOrdinal, Node > > | BlockDiagonalize (Level ¤tLevel, const RCP< Matrix > &A, bool generate_matrix) const |
void | BlockDiagonalizeGraph (const RCP< LWGraph > &inputGraph, const RCP< LocalOrdinalVector > &ghostedBlockNumber, RCP< LWGraph > &outputGraph, RCP< const Import > &importer) const |
Private Attributes | |
RCP< PreDropFunctionBaseClass > | predrop_ |
Constructors/Destructors. | |
CoalesceDropFactory () | |
Constructor. More... | |
virtual | ~CoalesceDropFactory () |
Destructor. More... | |
RCP< const ParameterList > | GetValidParameterList () const |
Return a const parameter list of valid parameters that setParameterList() will accept. More... | |
void | DeclareInput (Level ¤tLevel) const |
Input. More... | |
void | SetPreDropFunction (const RCP< MueLu::PreDropFunctionBaseClass< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &predrop) |
set predrop function More... | |
Additional Inherited Members | |
![]() | |
static void | EnableTimerSync () |
static void | DisableTimerSync () |
static void | EnableMultipleCheckGlobally () |
static void | DisableMultipleCheckGlobally () |
![]() | |
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 void | setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel) |
static EVerbosityLevel | getDefaultVerbLevel () |
![]() | |
static void | setDefaultOStream (const RCP< FancyOStream > &defaultOStream) |
static RCP< FancyOStream > | getDefaultOStream () |
![]() | |
static const EVerbosityLevel | verbLevel_default |
![]() | |
void | Input (Level &level, const std::string &varName) const |
void | Input (Level &level, const std::string &varName, const std::string &varParamName) const |
template<class T > | |
T | Get (Level &level, const std::string &varName) const |
template<class T > | |
T | Get (Level &level, const std::string &varName, const std::string &varParamName) const |
template<class T > | |
void | Set (Level &level, const std::string &varName, const T &data) const |
template<class T > | |
bool | IsType (Level &level, const std::string &varName) const |
bool | IsAvailable (Level &level, const std::string &varName) const |
![]() | |
void | initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null) |
![]() | |
void | initializeVerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null) |
virtual void | informUpdatedVerbosityState () const |
![]() | |
static bool | timerSync_ = false |
Factory for creating a graph based on a given matrix.
Factory for creating graphs from matrices with entries selectively dropped.
Both the classic dropping strategy as well as a coordinate-based distance laplacian method is implemented. For performance reasons there are four distinctive code paths for the classical method:
Additionally, there is a code path for the distance-laplacian mode.
Parameter | type | default | master.xml | validated | requested | description |
---|---|---|---|---|---|---|
A | Factory | null | * | * | Generating factory of the operator A | |
UnAmalgamationInfo | Factory | null | * | * | Generating factory of type AmalgamationFactory which generates the variable 'UnAmalgamationInfo'. Do not change the default unless you know what you are doing. | |
Coordinates | Factory | null | * | (*) | Generating factory for variable 'Coordinates'. The coordinates are only needed if "distance laplacian" is chosen for the parameter "aggregation: drop scheme" | |
"aggregation: drop scheme" | std::string | "classical" | * | * | Coalescing algorithm. You can choose either "classical" (=default) or "distance laplacian" | |
"aggregation: drop tol" | double | 0.0 | * | * | Threshold parameter for dropping small entries | |
"aggregation: Dirichlet threshold" | double | 0.0 | * | * | Threshold for determining whether entries are zero during Dirichlet row detection | |
"lightweight wrap" | bool | true | * | hidden switch between fast implementation based on MueLu::LWGraph and a failsafe slower implementation based on Xpetra::Graph (for comparison). The user should not change the default value (=true) |
The * in the master.xml
column denotes that the parameter is defined in the master.xml
file.
The * in the validated
column means that the parameter is declared in the list of valid input parameters (see CoalesceDropFactory::GetValidParameters).
The * in the requested
column states that the data is requested as input with all dependencies (see CoalesceDropFactory::DeclareInput).
After CoalesceDropFactory::Build the following data is available (if requested)
Parameter | generated by | description |
---|---|---|
Graph | CoalesceDropFactory | Graph of matrix A |
DofsPerNode | CoalesceDropFactory | number of DOFs per node. Note, that we assume a constant number of DOFs per node for all nodes associated with the operator A. |
The CoalesceDropFactory is internally using the AmalgamationFactory for amalgamating the dof-based maps to node-based maps. The AmalgamationFactory creates the "UnAmalgamationInfo" container which basically stores all the necessary information for translating dof based data to node based data and vice versa. The container is used, since this way the amalgamation is only done once and later reused by other factories.
Of course, often one does not need the information from the "UnAmalgamationInfo" container since the same information could be extracted of the "Graph" or the map from the "Coordinates" vector. However, there are also some situations (e.g. when doing rebalancing based on HyperGraph partitioning without coordinate information) where one has not access to a "Graph" or "Coordinates" variable.
Definition at line 98 of file MueLu_CoalesceDropFactory_decl.hpp.
MueLu::CoalesceDropFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node >::CoalesceDropFactory | ( | ) |
Constructor.
Definition at line 123 of file MueLu_CoalesceDropFactory_def.hpp.
|
inlinevirtual |
Destructor.
Definition at line 110 of file MueLu_CoalesceDropFactory_decl.hpp.
|
virtual |
Return a const parameter list of valid parameters that setParameterList() will accept.
Also define the default values of parameters according to the input parameter list.
Reimplemented from MueLu::Factory.
Definition at line 90 of file MueLu_CoalesceDropFactory_def.hpp.
|
virtual |
Input.
Implements MueLu::SingleLevelFactoryBase.
Definition at line 127 of file MueLu_CoalesceDropFactory_def.hpp.
|
inline |
set predrop function
Definition at line 122 of file MueLu_CoalesceDropFactory_decl.hpp.
|
virtual |
Build an object with this factory.
Implements MueLu::SingleLevelFactoryBase.
Definition at line 146 of file MueLu_CoalesceDropFactory_def.hpp.
|
private |
Method to merge rows of matrix for systems of PDEs.
Definition at line 1652 of file MueLu_CoalesceDropFactory_def.hpp.
|
private |
Definition at line 1711 of file MueLu_CoalesceDropFactory_def.hpp.
|
private |
Definition at line 1786 of file MueLu_CoalesceDropFactory_def.hpp.
|
private |
Definition at line 1909 of file MueLu_CoalesceDropFactory_def.hpp.
|
mutableprivate |
Definition at line 130 of file MueLu_CoalesceDropFactory_decl.hpp.