|
MueLu Version of the Day
|
Factory for building scalar Laplace operator (that is used as fake operator for variable dof size problems) More...
#include <MueLu_VariableDofLaplacianFactory_decl.hpp>
Constructors/Destructors. | |
| VariableDofLaplacianFactory () | |
| Constructor. | |
| virtual | ~VariableDofLaplacianFactory () |
| Destructor. | |
| RCP< const ParameterList > | GetValidParameterList () const |
| Return a const parameter list of valid parameters that setParameterList() will accept. | |
| void | DeclareInput (Level ¤tLevel) const |
| Input. | |
| void | Build (Level ¤tLevel) const |
| Build an object with this factory. | |
| void | buildPaddedMap (const Teuchos::ArrayRCP< const LocalOrdinal > &dofPresent, std::vector< LocalOrdinal > &map, size_t nDofs) const |
| void | assignGhostLocalNodeIds (const Teuchos::RCP< const Map > &rowDofMap, const Teuchos::RCP< const Map > &colDofMap, std::vector< LocalOrdinal > &myLocalNodeIds, const std::vector< LocalOrdinal > &dofMap, size_t maxDofPerNode, size_t &nLocalNodes, size_t &nLocalPlusGhostNodes, Teuchos::RCP< const Teuchos::Comm< int > > comm) const |
| void | squeezeOutNnzs (Teuchos::ArrayRCP< size_t > &rowPtr, Teuchos::ArrayRCP< LocalOrdinal > &cols, Teuchos::ArrayRCP< Scalar > &vals, const std::vector< bool > &keep) const |
| void | buildLaplacian (const Teuchos::ArrayRCP< size_t > &rowPtr, const Teuchos::ArrayRCP< LocalOrdinal > &cols, Teuchos::ArrayRCP< Scalar > &vals, const size_t &numdim, const RCP< Xpetra::MultiVector< typename Teuchos::ScalarTraits< Scalar >::magnitudeType, LocalOrdinal, GlobalOrdinal, Node > > &ghostedCoords) const |
| template<class listType > | |
| void | MueLu_az_sort (listType list[], size_t N, size_t list2[], Scalar list3[]) const |
Additional Inherited Members | |
Public Member Functions inherited from MueLu::SingleLevelFactoryBase | |
| SingleLevelFactoryBase () | |
| Constructor. | |
| virtual | ~SingleLevelFactoryBase () |
| Destructor. | |
| virtual void | CallBuild (Level &requestedLevel) const |
| virtual void | CallDeclareInput (Level &requestedLevel) const |
Public Member Functions inherited from MueLu::Factory | |
| Factory () | |
| Constructor. | |
| virtual | ~Factory () |
| Destructor. | |
| virtual void | SetFactory (const std::string &varName, const RCP< const FactoryBase > &factory) |
| Configuration. | |
| const RCP< const FactoryBase > | GetFactory (const std::string &varName) const |
| Default implementation of FactoryAcceptor::GetFactory() | |
| RCP< ParameterList > | RemoveFactoriesFromList (const ParameterList &list) const |
| void | EnableMultipleCallCheck () const |
| void | DisableMultipleCallCheck () const |
| void | ResetDebugData () const |
Public Member Functions inherited from MueLu::FactoryBase | |
| FactoryBase () | |
| Constructor. | |
| virtual | ~FactoryBase () |
| Destructor. | |
| int | GetID () const |
| return unique factory id | |
Public Member Functions inherited from MueLu::BaseClass | |
| virtual | ~BaseClass () |
| Destructor. | |
Public Member Functions inherited from MueLu::VerboseObject | |
| VerboseObject () | |
| virtual | ~VerboseObject () |
| Destructor. | |
| VerbLevel | GetVerbLevel () const |
| Get the verbosity level. | |
| void | SetVerbLevel (const VerbLevel verbLevel) |
| Set the verbosity level of this object. | |
| int | GetProcRankVerbose () const |
| Get proc rank used for printing. Do not use this information for any other purpose. | |
| int | SetProcRankVerbose (int procRank) const |
| Set proc rank used for printing. | |
| bool | IsPrint (MsgType type, int thisProcRankOnly=-1) const |
| Find out whether we need to print out information for a specific message type. | |
| Teuchos::FancyOStream & | GetOStream (MsgType type, int thisProcRankOnly=0) const |
| Get an output stream for outputting the input message type. | |
| Teuchos::FancyOStream & | GetBlackHole () const |
Public Member Functions inherited from MueLu::Describable | |
| virtual | ~Describable () |
| Destructor. | |
| 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. | |
| 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. | |
| virtual std::string | ShortClassName () const |
| Return the class name of the object, without template parameters and without namespace. | |
Public Member Functions inherited from MueLu::FactoryAcceptor | |
| virtual | ~FactoryAcceptor () |
Public Member Functions inherited from MueLu::ParameterListAcceptorImpl | |
| ParameterListAcceptorImpl () | |
| virtual | ~ParameterListAcceptorImpl ()=default |
| virtual void | SetParameterList (const Teuchos::ParameterList ¶mList) |
| Set parameters from a parameter list and return with default values. | |
| virtual const Teuchos::ParameterList & | GetParameterList () const |
| virtual const Teuchos::ParameterList & | GetParameterListWithoutValidation () const |
| void | SetParameter (const std::string &name, const ParameterEntry &entry) |
| Set a parameter directly as a ParameterEntry. | |
| const ParameterEntry & | GetParameter (const std::string &name) const |
| Retrieves a const entry with the name name. | |
| virtual void | GetDocumentation (std::ostream &os) const |
Public Member Functions inherited from MueLu::ParameterListAcceptor | |
| ParameterListAcceptor () | |
| virtual | ~ParameterListAcceptor ()=default |
Static Public Member Functions inherited from MueLu::Factory | |
| static void | EnableTimerSync () |
| static void | DisableTimerSync () |
| static void | EnableMultipleCheckGlobally () |
| static void | DisableMultipleCheckGlobally () |
Static Public Member Functions inherited from MueLu::VerboseObject | |
| static void | SetDefaultVerbLevel (const VerbLevel defaultVerbLevel) |
| Set the default (global) verbosity level. | |
| static VerbLevel | GetDefaultVerbLevel () |
| Get the default (global) verbosity level. | |
| static void | SetMueLuOStream (const Teuchos::RCP< Teuchos::FancyOStream > &mueluOStream) |
| static void | SetMueLuOFileStream (const std::string &filename) |
| static Teuchos::RCP< Teuchos::FancyOStream > | GetMueLuOStream () |
Protected Member Functions inherited from MueLu::Factory | |
| 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 |
Static Protected Attributes inherited from MueLu::Factory | |
| static bool | timerSync_ = false |
Factory for building scalar Laplace operator (that is used as fake operator for variable dof size problems)
Build distance Laplacian associated with input matrix A (which might have a variable number of DOFs per node). Coordinates are needed to calculate the distance laplacian values. The user-provided array "DofPresent" stores whether an array is present (=1) or not (=0) in the matrix. The length of the array is number of nodes * maxDofPerNode and therefore it might be larger or equal than the number of rows in the input matrix.
The factory produces the distance laplacian matrix A as output (with one dof per node) as well as the coarse version of the DofStatus (needed for the next coarser level), containing information about (artificial) Dirichlet rows in the matrix.
| Parameter | type | default | master.xml | validated | requested | description |
|---|---|---|---|---|---|---|
| A | Factory | null | * | * | Generating factory of the input matrix A with potentially variable number of DOFs. Might be padded or non-padded. Padded means, that the matrix has additional artificial rows and columns to have a constant number of DOFs per node. | |
| Coordinates | Factory | null | * | * | Generating factory for Coordinates needed for building distance laplacian. | |
| DofPresent | Teuchos::ArrayRCP<LocalOrdinal> | NoFactory | (*) | Optional array containing information whether DOF is actually present in matrix or not. | ||
| Advanced Dirichlet: threshold | double | 1e-5 | * | Drop tolerance for Dirichlet detection | ||
| Variable DOF amalgamation: threshold | double | 1.8e-9 | * | Drop tolerance for amalgamation process | ||
| maxDofPerNode | int | 1 | * | Maximum number of DOFs per node |
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 VariableDofLaplacianFactory::GetValidParameters).
The * in the requested column states that the data is requested as input with all dependencies (see VariableDofLaplacianFactory::DeclareInput).
After TentativePFactory::Build the following data is available (if requested)
| Parameter | generated by | description |
|---|---|---|
| A | VariableDofLaplacianFactory | Laplacian operator |
| DofStatus | VariableDofLaplacianFactory | Status array for next coarse level |
Definition at line 66 of file MueLu_VariableDofLaplacianFactory_decl.hpp.
| MueLu::VariableDofLaplacianFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node >::VariableDofLaplacianFactory | ( | ) |
Constructor.
Definition at line 34 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
inlinevirtual |
Destructor.
Definition at line 78 of file MueLu_VariableDofLaplacianFactory_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 20 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
virtual |
Input.
Implements MueLu::SingleLevelFactoryBase.
Definition at line 37 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
virtual |
Build an object with this factory.
Implements MueLu::SingleLevelFactoryBase.
Definition at line 48 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
private |
Definition at line 433 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
private |
Definition at line 441 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
private |
Definition at line 399 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
private |
Definition at line 344 of file MueLu_VariableDofLaplacianFactory_def.hpp.
|
inlineprivate |
Definition at line 100 of file MueLu_VariableDofLaplacianFactory_decl.hpp.