MueLu  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
MueLu_PgPFactory_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // MueLu: A package for multigrid based preconditioning
4 //
5 // Copyright 2012 NTESS and the MueLu contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef MUELU_PGPFACTORY_DECL_HPP_
11 #define MUELU_PGPFACTORY_DECL_HPP_
12 
13 #include <Xpetra_Vector_fwd.hpp>
15 #include <Xpetra_Matrix_fwd.hpp>
16 #include <Xpetra_Import_fwd.hpp>
18 #include <Xpetra_Export_fwd.hpp>
20 
21 #include "MueLu_ConfigDefs.hpp"
22 
23 #include "MueLu_PerfUtils_fwd.hpp"
24 #include "MueLu_PFactory.hpp"
25 #include "MueLu_Utilities_fwd.hpp"
26 
27 namespace MueLu {
28 
29 /* Options defining how to pick-up the next root node in the local aggregation procedure */
31  ANORM = 0, /* A norm */
32  L2NORM = 1, /* L2 norm */
33  DINVANORM = 2 /* Dinv A norm */
34 };
35 
42 template <class Scalar = DefaultScalar,
45  class Node = DefaultNode>
46 class PgPFactory : public PFactory {
47 #undef MUELU_PGPFACTORY_SHORT
48 #include "MueLu_UseShortNames.hpp"
49 
50  public:
52 
53 
58 
60  virtual ~PgPFactory() {}
61 
63 
65 
66 
68 
71 
75 
77 
78 
79  void DeclareInput(Level& fineLevel, Level& coarseLevel) const;
81 
83 
84 
90  void Build(Level& fineLevel, Level& coarseLevel) const;
91 
92  void BuildP(Level& fineLevel, Level& coarseLevel) const;
93 
95 
96  void ReUseDampingParameters(bool bReuse);
97 
98  private:
100 
102 
103  void ComputeRowBasedOmega(Level& fineLevel, Level& coarseLevel, const RCP<Matrix>& A, const RCP<Matrix>& P0, const RCP<Matrix>& DinvAP0, RCP<Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> >& RowBasedOmega) const;
104 
105  private:
107  std::string diagonalView_; // TODO do we need this?
108 };
109 
110 } // namespace MueLu
111 
112 #define MUELU_PGPFACTORY_SHORT
113 #endif /* MUELU_PGPFACTORY_DECL_HPP_ */
MueLu::DefaultLocalOrdinal LocalOrdinal
Tpetra::KokkosClassic::DefaultNode::DefaultNodeType DefaultNode
void DeclareInput(Level &fineLevel, Level &coarseLevel) const
Input.
void BuildP(Level &fineLevel, Level &coarseLevel) const
Abstract Build method.
void ComputeRowBasedOmega(Level &fineLevel, Level &coarseLevel, const RCP< Matrix > &A, const RCP< Matrix > &P0, const RCP< Matrix > &DinvAP0, RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &RowBasedOmega) const
MueLu::DefaultNode Node
void SetMinimizationMode(MinimizationNorm minnorm)
Set minimization mode (L2NORM for cheapest, ANORM more expensive, DINVANORM = default) ...
MueLu::DefaultScalar Scalar
Tpetra::Details::DefaultTypes::scalar_type DefaultScalar
void Build(Level &fineLevel, Level &coarseLevel) const
Build method.
MueLu::DefaultGlobalOrdinal GlobalOrdinal
Class that holds all level-specific information.
Definition: MueLu_Level.hpp:63
void MultiplySelfAll(const RCP< Matrix > &Op, Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &InnerProdVec) const
std::string diagonalView_
Factory parameters.
MinimizationNorm GetMinimizationMode()
return minimization mode
void ReUseDampingParameters(bool bReuse)
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
PgPFactory()
Constructor. User can supply a factory for generating the tentative prolongator.
virtual ~PgPFactory()
Destructor.
Factory that provides an interface for a concrete implementation of a prolongation operator...
Factory for building Petrov-Galerkin Smoothed Aggregation prolongators.
void MultiplyAll(const RCP< Matrix > &left, const RCP< Matrix > &right, Teuchos::RCP< Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &InnerProdVec) const