10#ifndef IFPACK2_DIAGONAL_DEF_HPP
11#define IFPACK2_DIAGONAL_DEF_HPP
13#include "Ifpack2_Diagonal_decl.hpp"
14#include "Tpetra_CrsMatrix.hpp"
18template <
class MatrixType>
21 , initializeTime_(0.0)
27 , isInitialized_(false)
28 , isComputed_(false) {}
30template <
class MatrixType>
33 , initializeTime_(0.0)
39 , isInitialized_(false)
40 , isComputed_(false) {}
42template <
class MatrixType>
44 : userInverseDiag_(diag)
46 , initializeTime_(0.0)
52 , isInitialized_(false)
53 , isComputed_(false) {}
55template <
class MatrixType>
58template <
class MatrixType>
59Teuchos::RCP<const typename Diagonal<MatrixType>::map_type>
61 if (matrix_.is_null()) {
62 if (userInverseDiag_.is_null()) {
63 TEUCHOS_TEST_FOR_EXCEPTION(
64 true, std::runtime_error,
65 "Ifpack2::Diagonal::getDomainMap: "
66 "The input matrix A is null, and you did not provide a vector of "
67 "inverse diagonal entries. Please call setMatrix() with a nonnull "
68 "input matrix before calling this method.");
70 return userInverseDiag_->getMap();
73 return matrix_->getDomainMap();
77template <
class MatrixType>
78Teuchos::RCP<const typename Diagonal<MatrixType>::map_type>
80 if (matrix_.is_null()) {
81 if (userInverseDiag_.is_null()) {
82 TEUCHOS_TEST_FOR_EXCEPTION(
83 true, std::runtime_error,
84 "Ifpack2::Diagonal::getRangeMap: "
85 "The input matrix A is null, and you did not provide a vector of "
86 "inverse diagonal entries. Please call setMatrix() with a nonnull "
87 "input matrix before calling this method.");
89 return userInverseDiag_->getMap();
92 return matrix_->getRangeMap();
96template <
class MatrixType>
100template <
class MatrixType>
102 inverseDiag_ = Teuchos::null;
103 offsets_ = offsets_type();
104 isInitialized_ =
false;
108template <
class MatrixType>
110 setMatrix(
const Teuchos::RCP<const row_matrix_type>& A) {
111 if (A.getRawPtr() != matrix_.getRawPtr()) {
117template <
class MatrixType>
121 TEUCHOS_TEST_FOR_EXCEPTION(
122 matrix_.is_null() && userInverseDiag_.is_null(), std::runtime_error,
123 "Ifpack2::Diagonal::initialize: The matrix to precondition is null, "
124 "and you did not provide a Tpetra::Vector of diagonal entries. "
125 "Please call setMatrix() with a nonnull input before calling this method.");
133 if (!matrix_.is_null()) {
138 Teuchos::RCP<const crs_matrix_type> A_crs =
139 Teuchos::rcp_dynamic_cast<const crs_matrix_type>(matrix_);
141 if (A_crs.is_null()) {
142 offsets_ = offsets_type();
144 const size_t lclNumRows = A_crs->getLocalNumRows();
145 if (offsets_.extent(0) < lclNumRows) {
146 offsets_ = offsets_type();
147 offsets_ = offsets_type(
"offsets", lclNumRows);
149 A_crs->getCrsGraph()->getLocalDiagOffsets(offsets_);
153 isInitialized_ =
true;
157template <
class MatrixType>
161 TEUCHOS_TEST_FOR_EXCEPTION(
162 matrix_.is_null() && userInverseDiag_.is_null(), std::runtime_error,
163 "Ifpack2::Diagonal::compute: The matrix to precondition is null, "
164 "and you did not provide a Tpetra::Vector of diagonal entries. "
165 "Please call setMatrix() with a nonnull input before calling this method.");
167 if (!isInitialized_) {
177 if (matrix_.is_null()) {
178 inverseDiag_ = userInverseDiag_;
180 Teuchos::RCP<vector_type> tmpVec(
new vector_type(matrix_->getRowMap()));
181 Teuchos::RCP<const crs_matrix_type> A_crs =
182 Teuchos::rcp_dynamic_cast<const crs_matrix_type>(matrix_);
183 if (A_crs.is_null()) {
185 matrix_->getLocalDiagCopy(*tmpVec);
189 A_crs->getLocalDiagCopy(*tmpVec, offsets_);
191 tmpVec->reciprocal(*tmpVec);
192 inverseDiag_ = tmpVec;
199template <
class MatrixType>
201 apply(
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
202 Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
205 scalar_type beta)
const {
206 TEUCHOS_TEST_FOR_EXCEPTION(
207 !isComputed(), std::runtime_error,
208 "Ifpack2::Diagonal::apply: You "
209 "must first call compute() before you may call apply(). Once you have "
210 "called compute(), you need not call it again unless the values in the "
211 "matrix have changed, or unless you have called setMatrix().");
218 Y.elementWiseMultiply(alpha, *inverseDiag_, X, beta);
222template <
class MatrixType>
224 return numInitialize_;
227template <
class MatrixType>
232template <
class MatrixType>
237template <
class MatrixType>
239 return initializeTime_;
242template <
class MatrixType>
247template <
class MatrixType>
252template <
class MatrixType>
254 std::ostringstream out;
259 out <<
"\"Ifpack2::Diagonal\": "
261 if (this->getObjectLabel() !=
"") {
262 out <<
"Label: \"" << this->getObjectLabel() <<
"\", ";
264 if (matrix_.is_null()) {
265 out <<
"Matrix: null";
267 out <<
"Matrix: not null"
268 <<
", Global matrix dimensions: ["
269 << matrix_->getGlobalNumRows() <<
", "
270 << matrix_->getGlobalNumCols() <<
"]";
277template <
class MatrixType>
279 describe(Teuchos::FancyOStream& out,
280 const Teuchos::EVerbosityLevel verbLevel)
const {
283 const Teuchos::EVerbosityLevel vl =
284 (verbLevel == Teuchos::VERB_DEFAULT) ? Teuchos::VERB_LOW : verbLevel;
285 if (vl != Teuchos::VERB_NONE) {
286 Teuchos::OSTab tab0(out);
287 out <<
"\"Ifpack2::Diagonal\":";
288 Teuchos::OSTab tab1(out);
289 out <<
"Template parameter: "
290 << Teuchos::TypeNameTraits<MatrixType>::name() << endl;
291 if (this->getObjectLabel() !=
"") {
292 out <<
"Label: \"" << this->getObjectLabel() <<
"\", ";
294 out <<
"Number of initialize calls: " << numInitialize_ << endl
295 <<
"Number of compute calls: " << numCompute_ << endl
296 <<
"Number of apply calls: " << numApply_ << endl;
302#define IFPACK2_DIAGONAL_INSTANT(S, LO, GO, N) \
303 template class Ifpack2::Diagonal<Tpetra::RowMatrix<S, LO, GO, N> >;
Diagonal preconditioner.
Definition Ifpack2_Diagonal_decl.hpp:45
virtual ~Diagonal()
Destructor.
Definition Ifpack2_Diagonal_def.hpp:56
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
Teuchos::RCP< const map_type > getDomainMap() const
The Tpetra::Map representing this operator's domain.
Definition Ifpack2_Diagonal_def.hpp:60
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's range.
Definition Ifpack2_Diagonal_def.hpp:79
int getNumCompute() const
Return the number of calls to compute().
Definition Ifpack2_Diagonal_def.hpp:228
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_Diagonal_def.hpp:110
void setParameters(const Teuchos::ParameterList ¶ms)
Sets parameters on this object.
Definition Ifpack2_Diagonal_def.hpp:98
double getComputeTime() const
Return the time spent in compute().
Definition Ifpack2_Diagonal_def.hpp:243
int getNumApply() const
Return the number of calls to apply().
Definition Ifpack2_Diagonal_def.hpp:233
Diagonal(const Teuchos::RCP< const row_matrix_type > &A)
Constructor that takes a Tpetra::RowMatrix.
Definition Ifpack2_Diagonal_def.hpp:19
int getNumInitialize() const
Return the number of calls to initialize().
Definition Ifpack2_Diagonal_def.hpp:223
double getApplyTime() const
Return the time spent in apply().
Definition Ifpack2_Diagonal_def.hpp:248
void initialize()
Initialize.
Definition Ifpack2_Diagonal_def.hpp:118
void compute()
Compute the preconditioner.
Definition Ifpack2_Diagonal_def.hpp:158
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
double getInitializeTime() const
Return the time spent in initialize().
Definition Ifpack2_Diagonal_def.hpp:238
std::string description() const
Return a one-line description of this object.
Definition Ifpack2_Diagonal_def.hpp:253
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40