Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_Diagonal_decl.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4 //
5 // Copyright 2009 NTESS and the Ifpack2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef IFPACK2_DIAGONAL_DECL_HPP
11 #define IFPACK2_DIAGONAL_DECL_HPP
12 
15 #include "Tpetra_CrsMatrix_decl.hpp"
16 #include <type_traits>
17 
18 namespace Ifpack2 {
19 
37 template<class MatrixType>
38 class Diagonal :
39  virtual public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
40  typename MatrixType::local_ordinal_type,
41  typename MatrixType::global_ordinal_type,
42  typename MatrixType::node_type>,
43  virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
44  typename MatrixType::local_ordinal_type,
45  typename MatrixType::global_ordinal_type,
46  typename MatrixType::node_type> >
47 {
48 public:
49  typedef typename MatrixType::scalar_type scalar_type;
50  typedef typename MatrixType::local_ordinal_type local_ordinal_type;
51  typedef typename MatrixType::global_ordinal_type global_ordinal_type;
52  typedef typename MatrixType::node_type node_type;
53  typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
54 
56  typedef Tpetra::RowMatrix<scalar_type,
57  local_ordinal_type,
58  global_ordinal_type,
59  node_type> row_matrix_type;
60 
61  static_assert(std::is_same<MatrixType, row_matrix_type>::value, "Ifpack2::Diagonal: The template parameter MatrixType must be a Tpetra::RowMatrix specialization. Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore. The constructor can take either a RowMatrix or a CrsMatrix just fine.");
62 
64  typedef Tpetra::CrsMatrix<scalar_type,
65  local_ordinal_type,
66  global_ordinal_type,
67  node_type> crs_matrix_type;
69  typedef Tpetra::Vector<scalar_type,
70  local_ordinal_type,
71  global_ordinal_type,
72  node_type> vector_type;
74  typedef Tpetra::Map<local_ordinal_type,
75  global_ordinal_type,
76  node_type> map_type;
77 
82 
91 
104 
106  virtual ~Diagonal ();
107 
109 
112  void setParameters (const Teuchos::ParameterList& params);
113 
115  void initialize ();
116 
118  bool isInitialized () const {
119  return isInitialized_;
120  }
121 
123  void compute ();
124 
126  bool isComputed () const {
127  return isComputed_;
128  }
129 
131 
133 
156  virtual void
158 
160 
162 
168  void
169  apply (const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
170  Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
172  scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
173  scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
174 
177 
180 
182 
184 
186  //Teuchos::RCP<const Teuchos::Comm<int> > getComm () const;
187 
194  return matrix_;
195  }
196 
198  double getComputeFlops() const;
199 
201  double getApplyFlops() const;
202 
204  int getNumInitialize() const;
205 
207  int getNumCompute() const;
208 
210  int getNumApply() const;
211 
213  double getInitializeTime() const;
214 
216  double getComputeTime() const;
217 
219  double getApplyTime() const;
220 
222 
224 
226  std::string description() const;
227 
229  void
231  const Teuchos::EVerbosityLevel verbLevel =
234 
235 private:
237  void reset ();
238 
241 
246  Teuchos::RCP<const vector_type> userInverseDiag_;
247 
249  Teuchos::RCP<const vector_type> inverseDiag_;
250 
251  typedef Kokkos::View<size_t*, typename node_type::device_type> offsets_type;
252  offsets_type offsets_;
253 
254  double initializeTime_;
255  double computeTime_;
256  mutable double applyTime_;
257 
258  int numInitialize_;
259  int numCompute_;
260  mutable int numApply_;
261 
262  bool isInitialized_;
263  bool isComputed_;
264 };
265 
282 template<class MatrixType, class VectorType>
283 Teuchos::RCP<Ifpack2::Diagonal<Tpetra::RowMatrix<typename MatrixType::scalar_type,
284  typename MatrixType::local_ordinal_type,
285  typename MatrixType::global_ordinal_type,
286  typename MatrixType::node_type> > >
288 {
289  typedef Tpetra::RowMatrix<typename MatrixType::scalar_type,
290  typename MatrixType::local_ordinal_type,
291  typename MatrixType::global_ordinal_type,
292  typename MatrixType::node_type> row_matrix_type;
293 
294  return Teuchos::rcp (new Ifpack2::Diagonal<row_matrix_type> (invdiag));
295 }
296 
297 }//namespace Ifpack2
298 
299 #endif
std::string description() const
Return a one-line description of this object.
Definition: Ifpack2_Diagonal_def.hpp:265
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_Diagonal_decl.hpp:59
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:60
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the preconditioner to X, putting the result in Y.
Definition: Ifpack2_Diagonal_def.hpp:213
Diagonal preconditioner.
Definition: Ifpack2_Diagonal_decl.hpp:38
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_Diagonal_decl.hpp:118
double getApplyFlops() const
Return the number of flops for the application of the preconditioner.
Teuchos::RCP< const map_type > getDomainMap() const
The Tpetra::Map representing this operator&#39;s domain.
Definition: Ifpack2_Diagonal_def.hpp:64
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
double getComputeFlops() const
Return the number of flops in the computation phase.
Diagonal(const Teuchos::RCP< const row_matrix_type > &A)
Constructor that takes a Tpetra::RowMatrix.
Definition: Ifpack2_Diagonal_def.hpp:19
void compute()
Compute the preconditioner.
Definition: Ifpack2_Diagonal_def.hpp:167
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
virtual ~Diagonal()
Destructor.
Definition: Ifpack2_Diagonal_def.hpp:59
Teuchos::RCP< Ifpack2::Diagonal< Tpetra::RowMatrix< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > > createDiagonalPreconditioner(const Teuchos::RCP< const VectorType > &invdiag)
Definition: Ifpack2_Diagonal_decl.hpp:287
double getComputeTime() const
Return the time spent in compute().
Definition: Ifpack2_Diagonal_def.hpp:255
int getNumCompute() const
Return the number of calls to compute().
Definition: Ifpack2_Diagonal_def.hpp:240
Declaration of interface for preconditioners that can change their matrix after construction.
double getInitializeTime() const
Return the time spent in initialize().
Definition: Ifpack2_Diagonal_def.hpp:250
static const EVerbosityLevel verbLevel_default
void setParameters(const Teuchos::ParameterList &params)
Sets parameters on this object.
Definition: Ifpack2_Diagonal_def.hpp:102
int getNumInitialize() const
Return the number of calls to initialize().
Definition: Ifpack2_Diagonal_def.hpp:235
Teuchos::RCP< const row_matrix_type > getMatrix() const
Return the communicator associated with this matrix operator.
Definition: Ifpack2_Diagonal_decl.hpp:193
double getApplyTime() const
Return the time spent in apply().
Definition: Ifpack2_Diagonal_def.hpp:260
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Diagonal_def.hpp:116
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.
Definition: Ifpack2_Diagonal_def.hpp:293
Teuchos::RCP< const map_type > getRangeMap() const
The Tpetra::Map representing this operator&#39;s range.
Definition: Ifpack2_Diagonal_def.hpp:83
void initialize()
Initialize.
Definition: Ifpack2_Diagonal_def.hpp:125
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Tpetra::Map specialization used by this class.
Definition: Ifpack2_Diagonal_decl.hpp:76
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class.
Definition: Ifpack2_Diagonal_decl.hpp:61
int getNumApply() const
Return the number of calls to apply().
Definition: Ifpack2_Diagonal_def.hpp:245
Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > vector_type
Tpetra::Vector specialization used by this class.
Definition: Ifpack2_Diagonal_decl.hpp:72
bool isComputed() const
Return true if compute() has been called.
Definition: Ifpack2_Diagonal_decl.hpp:126