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),
31template<
class MatrixType>
34 initializeTime_ (0.0),
40 isInitialized_ (false),
44template<
class MatrixType>
46 userInverseDiag_ (diag),
48 initializeTime_ (0.0),
54 isInitialized_ (false),
58template<
class MatrixType>
62template<
class MatrixType>
63Teuchos::RCP<const typename Diagonal<MatrixType>::map_type>
66 if (matrix_.is_null ()) {
67 if (userInverseDiag_.is_null ()) {
68 TEUCHOS_TEST_FOR_EXCEPTION(
69 true, std::runtime_error,
"Ifpack2::Diagonal::getDomainMap: "
70 "The input matrix A is null, and you did not provide a vector of "
71 "inverse diagonal entries. Please call setMatrix() with a nonnull "
72 "input matrix before calling this method.");
74 return userInverseDiag_->getMap ();
77 return matrix_->getDomainMap ();
81template<
class MatrixType>
82Teuchos::RCP<const typename Diagonal<MatrixType>::map_type>
85 if (matrix_.is_null ()) {
86 if (userInverseDiag_.is_null ()) {
87 TEUCHOS_TEST_FOR_EXCEPTION(
88 true, std::runtime_error,
"Ifpack2::Diagonal::getRangeMap: "
89 "The input matrix A is null, and you did not provide a vector of "
90 "inverse diagonal entries. Please call setMatrix() with a nonnull "
91 "input matrix before calling this method.");
93 return userInverseDiag_->getMap ();
96 return matrix_->getRangeMap ();
100template<
class MatrixType>
105template<
class MatrixType>
106void Diagonal<MatrixType>::reset ()
108 inverseDiag_ = Teuchos::null;
109 offsets_ = offsets_type ();
110 isInitialized_ =
false;
114template<
class MatrixType>
116setMatrix (
const Teuchos::RCP<const row_matrix_type>& A)
118 if (A.getRawPtr () != matrix_.getRawPtr ()) {
124template<
class MatrixType>
129 TEUCHOS_TEST_FOR_EXCEPTION(
130 matrix_.is_null () && userInverseDiag_.is_null (), std::runtime_error,
131 "Ifpack2::Diagonal::initialize: The matrix to precondition is null, "
132 "and you did not provide a Tpetra::Vector of diagonal entries. "
133 "Please call setMatrix() with a nonnull input before calling this method.");
141 if (! matrix_.is_null ()) {
146 Teuchos::RCP<const crs_matrix_type> A_crs =
147 Teuchos::rcp_dynamic_cast<const crs_matrix_type> (matrix_);
149 if (A_crs.is_null ()) {
150 offsets_ = offsets_type ();
153 const size_t lclNumRows = A_crs->getLocalNumRows ();
154 if (offsets_.extent (0) < lclNumRows) {
155 offsets_ = offsets_type ();
156 offsets_ = offsets_type (
"offsets", lclNumRows);
158 A_crs->getCrsGraph ()->getLocalDiagOffsets (offsets_);
162 isInitialized_ =
true;
166template<
class MatrixType>
171 TEUCHOS_TEST_FOR_EXCEPTION(
172 matrix_.is_null () && userInverseDiag_.is_null (), std::runtime_error,
173 "Ifpack2::Diagonal::compute: The matrix to precondition is null, "
174 "and you did not provide a Tpetra::Vector of diagonal entries. "
175 "Please call setMatrix() with a nonnull input before calling this method.");
177 if (! isInitialized_) {
187 if (matrix_.is_null ()) {
188 inverseDiag_ = userInverseDiag_;
191 Teuchos::RCP<vector_type> tmpVec (
new vector_type (matrix_->getRowMap ()));
192 Teuchos::RCP<const crs_matrix_type> A_crs =
193 Teuchos::rcp_dynamic_cast<const crs_matrix_type> (matrix_);
194 if (A_crs.is_null ()) {
196 matrix_->getLocalDiagCopy (*tmpVec);
201 A_crs->getLocalDiagCopy (*tmpVec, offsets_);
203 tmpVec->reciprocal (*tmpVec);
204 inverseDiag_ = tmpVec;
211template<
class MatrixType>
213apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
214 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
217 scalar_type beta)
const
219 TEUCHOS_TEST_FOR_EXCEPTION(
220 !
isComputed (), std::runtime_error,
"Ifpack2::Diagonal::apply: You "
221 "must first call compute() before you may call apply(). Once you have "
222 "called compute(), you need not call it again unless the values in the "
223 "matrix have changed, or unless you have called setMatrix().");
230 Y.elementWiseMultiply (alpha, *inverseDiag_, X, beta);
234template <
class MatrixType>
236 return numInitialize_;
239template <
class MatrixType>
244template <
class MatrixType>
249template <
class MatrixType>
251 return initializeTime_;
254template<
class MatrixType>
259template<
class MatrixType>
264template <
class MatrixType>
267 std::ostringstream out;
272 out <<
"\"Ifpack2::Diagonal\": "
274 if (this->getObjectLabel () !=
"") {
275 out <<
"Label: \"" << this->getObjectLabel () <<
"\", ";
277 if (matrix_.is_null ()) {
278 out <<
"Matrix: null";
281 out <<
"Matrix: not null"
282 <<
", Global matrix dimensions: ["
283 << matrix_->getGlobalNumRows () <<
", "
284 << matrix_->getGlobalNumCols () <<
"]";
291template <
class MatrixType>
293describe (Teuchos::FancyOStream &out,
294 const Teuchos::EVerbosityLevel verbLevel)
const
298 const Teuchos::EVerbosityLevel vl =
299 (verbLevel == Teuchos::VERB_DEFAULT) ? Teuchos::VERB_LOW : verbLevel;
300 if (vl != Teuchos::VERB_NONE) {
301 Teuchos::OSTab tab0 (out);
302 out <<
"\"Ifpack2::Diagonal\":";
303 Teuchos::OSTab tab1 (out);
304 out <<
"Template parameter: "
305 << Teuchos::TypeNameTraits<MatrixType>::name () << endl;
306 if (this->getObjectLabel () !=
"") {
307 out <<
"Label: \"" << this->getObjectLabel () <<
"\", ";
309 out <<
"Number of initialize calls: " << numInitialize_ << endl
310 <<
"Number of compute calls: " << numCompute_ << endl
311 <<
"Number of apply calls: " << numApply_ << endl;
317#define IFPACK2_DIAGONAL_INSTANT(S,LO,GO,N) \
318 template class Ifpack2::Diagonal< Tpetra::RowMatrix<S, LO, GO, N> >;
virtual ~Diagonal()
Destructor.
Definition Ifpack2_Diagonal_def.hpp:59
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
Teuchos::RCP< const map_type > getDomainMap() const
The Tpetra::Map representing this operator's domain.
Definition Ifpack2_Diagonal_def.hpp:64
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's range.
Definition Ifpack2_Diagonal_def.hpp:83
int getNumCompute() const
Return the number of calls to compute().
Definition Ifpack2_Diagonal_def.hpp:240
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition Ifpack2_Diagonal_def.hpp:116
void setParameters(const Teuchos::ParameterList ¶ms)
Sets parameters on this object.
Definition Ifpack2_Diagonal_def.hpp:102
double getComputeTime() const
Return the time spent in compute().
Definition Ifpack2_Diagonal_def.hpp:255
int getNumApply() const
Return the number of calls to apply().
Definition Ifpack2_Diagonal_def.hpp:245
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:235
double getApplyTime() const
Return the time spent in apply().
Definition Ifpack2_Diagonal_def.hpp:260
void initialize()
Initialize.
Definition Ifpack2_Diagonal_def.hpp:125
void compute()
Compute the preconditioner.
Definition Ifpack2_Diagonal_def.hpp:167
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
double getInitializeTime() const
Return the time spent in initialize().
Definition Ifpack2_Diagonal_def.hpp:250
std::string description() const
Return a one-line description of this object.
Definition Ifpack2_Diagonal_def.hpp:265
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41