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 : virtual public Ifpack2::Preconditioner<typename MatrixType::scalar_type,
39  typename MatrixType::local_ordinal_type,
40  typename MatrixType::global_ordinal_type,
41  typename MatrixType::node_type>,
42  virtual public Ifpack2::Details::CanChangeMatrix<Tpetra::RowMatrix<typename MatrixType::scalar_type,
43  typename MatrixType::local_ordinal_type,
44  typename MatrixType::global_ordinal_type,
45  typename MatrixType::node_type> > {
46  public:
47  typedef typename MatrixType::scalar_type scalar_type;
48  typedef typename MatrixType::local_ordinal_type local_ordinal_type;
49  typedef typename MatrixType::global_ordinal_type global_ordinal_type;
50  typedef typename MatrixType::node_type node_type;
51  typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType magnitude_type;
52 
54  typedef Tpetra::RowMatrix<scalar_type,
55  local_ordinal_type,
56  global_ordinal_type,
57  node_type>
59 
60  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.");
61 
63  typedef Tpetra::CrsMatrix<scalar_type,
64  local_ordinal_type,
65  global_ordinal_type,
66  node_type>
69  typedef Tpetra::Vector<scalar_type,
70  local_ordinal_type,
71  global_ordinal_type,
72  node_type>
75  typedef Tpetra::Map<local_ordinal_type,
76  global_ordinal_type,
77  node_type>
79 
84 
93 
106 
108  virtual ~Diagonal();
109 
111 
114  void setParameters(const Teuchos::ParameterList& params);
115 
117  void initialize();
118 
120  bool isInitialized() const {
121  return isInitialized_;
122  }
123 
125  void compute();
126 
128  bool isComputed() const {
129  return isComputed_;
130  }
131 
133 
135 
158  virtual void
160 
162 
164 
170  void
171  apply(const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
172  Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
174  scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
175  scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero()) const;
176 
179 
182 
184 
186 
188  // Teuchos::RCP<const Teuchos::Comm<int> > getComm () const;
189 
196  return matrix_;
197  }
198 
200  double getComputeFlops() const;
201 
203  double getApplyFlops() const;
204 
206  int getNumInitialize() const;
207 
209  int getNumCompute() const;
210 
212  int getNumApply() const;
213 
215  double getInitializeTime() const;
216 
218  double getComputeTime() const;
219 
221  double getApplyTime() const;
222 
224 
226 
228  std::string description() const;
229 
231  void
233  const Teuchos::EVerbosityLevel verbLevel =
236 
237  private:
239  void reset();
240 
243 
248  Teuchos::RCP<const vector_type> userInverseDiag_;
249 
251  Teuchos::RCP<const vector_type> inverseDiag_;
252 
253  typedef Kokkos::View<size_t*, typename node_type::device_type> offsets_type;
254  offsets_type offsets_;
255 
256  double initializeTime_;
257  double computeTime_;
258  mutable double applyTime_;
259 
260  int numInitialize_;
261  int numCompute_;
262  mutable int numApply_;
263 
264  bool isInitialized_;
265  bool isComputed_;
266 };
267 
284 template <class MatrixType, class VectorType>
285 Teuchos::RCP<Ifpack2::Diagonal<Tpetra::RowMatrix<typename MatrixType::scalar_type,
286  typename MatrixType::local_ordinal_type,
287  typename MatrixType::global_ordinal_type,
288  typename MatrixType::node_type> > >
290  typedef Tpetra::RowMatrix<typename MatrixType::scalar_type,
291  typename MatrixType::local_ordinal_type,
292  typename MatrixType::global_ordinal_type,
293  typename MatrixType::node_type>
294  row_matrix_type;
295 
297 }
298 
299 } // namespace Ifpack2
300 
301 #endif
std::string description() const
Return a one-line description of this object.
Definition: Ifpack2_Diagonal_def.hpp:253
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:58
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:201
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:120
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:60
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:158
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
virtual ~Diagonal()
Destructor.
Definition: Ifpack2_Diagonal_def.hpp:56
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:289
double getComputeTime() const
Return the time spent in compute().
Definition: Ifpack2_Diagonal_def.hpp:243
int getNumCompute() const
Return the number of calls to compute().
Definition: Ifpack2_Diagonal_def.hpp:228
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:238
static const EVerbosityLevel verbLevel_default
void setParameters(const Teuchos::ParameterList &params)
Sets parameters on this object.
Definition: Ifpack2_Diagonal_def.hpp:98
int getNumInitialize() const
Return the number of calls to initialize().
Definition: Ifpack2_Diagonal_def.hpp:223
Teuchos::RCP< const row_matrix_type > getMatrix() const
Return the communicator associated with this matrix operator.
Definition: Ifpack2_Diagonal_decl.hpp:195
double getApplyTime() const
Return the time spent in apply().
Definition: Ifpack2_Diagonal_def.hpp:248
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Diagonal_def.hpp:110
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:279
Teuchos::RCP< const map_type > getRangeMap() const
The Tpetra::Map representing this operator&#39;s range.
Definition: Ifpack2_Diagonal_def.hpp:79
void initialize()
Initialize.
Definition: Ifpack2_Diagonal_def.hpp:118
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Tpetra::Map specialization used by this class.
Definition: Ifpack2_Diagonal_decl.hpp:78
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:60
int getNumApply() const
Return the number of calls to apply().
Definition: Ifpack2_Diagonal_def.hpp:233
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:73
bool isComputed() const
Return true if compute() has been called.
Definition: Ifpack2_Diagonal_decl.hpp:128