Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Diagonal_def.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_DEF_HPP
11#define IFPACK2_DIAGONAL_DEF_HPP
12
13#include "Ifpack2_Diagonal_decl.hpp"
14#include "Tpetra_CrsMatrix.hpp"
15
16namespace Ifpack2 {
17
18template<class MatrixType>
19Diagonal<MatrixType>::Diagonal (const Teuchos::RCP<const row_matrix_type>& A) :
20 matrix_ (A),
21 initializeTime_ (0.0),
22 computeTime_ (0.0),
23 applyTime_ (0.0),
24 numInitialize_ (0),
25 numCompute_ (0),
26 numApply_ (0),
27 isInitialized_ (false),
28 isComputed_ (false)
29{}
30
31template<class MatrixType>
32Diagonal<MatrixType>::Diagonal (const Teuchos::RCP<const crs_matrix_type>& A) :
33 matrix_ (A),
34 initializeTime_ (0.0),
35 computeTime_ (0.0),
36 applyTime_ (0.0),
37 numInitialize_ (0),
38 numCompute_ (0),
39 numApply_ (0),
40 isInitialized_ (false),
41 isComputed_ (false)
42{}
43
44template<class MatrixType>
45Diagonal<MatrixType>::Diagonal (const Teuchos::RCP<const vector_type>& diag) :
46 userInverseDiag_ (diag),
47 inverseDiag_ (diag),
48 initializeTime_ (0.0),
49 computeTime_ (0.0),
50 applyTime_ (0.0),
51 numInitialize_ (0),
52 numCompute_ (0),
53 numApply_ (0),
54 isInitialized_ (false),
55 isComputed_ (false)
56{}
57
58template<class MatrixType>
61
62template<class MatrixType>
63Teuchos::RCP<const typename Diagonal<MatrixType>::map_type>
65{
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.");
73 } else {
74 return userInverseDiag_->getMap ();
75 }
76 } else {
77 return matrix_->getDomainMap ();
78 }
79}
80
81template<class MatrixType>
82Teuchos::RCP<const typename Diagonal<MatrixType>::map_type>
84{
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.");
92 } else {
93 return userInverseDiag_->getMap ();
94 }
95 } else {
96 return matrix_->getRangeMap ();
97 }
98}
99
100template<class MatrixType>
102setParameters (const Teuchos::ParameterList& /*params*/)
103{}
104
105template<class MatrixType>
106void Diagonal<MatrixType>::reset ()
107{
108 inverseDiag_ = Teuchos::null;
109 offsets_ = offsets_type ();
110 isInitialized_ = false;
111 isComputed_ = false;
112}
113
114template<class MatrixType>
116setMatrix (const Teuchos::RCP<const row_matrix_type>& A)
117{
118 if (A.getRawPtr () != matrix_.getRawPtr ()) { // it's a different matrix
119 reset ();
120 matrix_ = A;
121 }
122}
123
124template<class MatrixType>
126{
127 // Either the matrix to precondition must be nonnull, or the user
128 // must have provided a Vector of diagonal entries.
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.");
134
135 // If the user did provide an input matrix, then that takes
136 // precedence over the vector of inverse diagonal entries, if they
137 // provided one earlier. This is only possible if they created this
138 // Diagonal instance using the constructor that takes a
139 // Tpetra::Vector pointer, and then called setMatrix() with a
140 // nonnull input matrix.
141 if (! matrix_.is_null ()) {
142 // If you call initialize(), it means that you are asserting that
143 // the structure of the input sparse matrix may have changed.
144 // This means we should always recompute the diagonal offsets, if
145 // the input matrix is a Tpetra::CrsMatrix.
146 Teuchos::RCP<const crs_matrix_type> A_crs =
147 Teuchos::rcp_dynamic_cast<const crs_matrix_type> (matrix_);
148
149 if (A_crs.is_null ()) {
150 offsets_ = offsets_type (); // offsets are no longer valid
151 }
152 else {
153 const size_t lclNumRows = A_crs->getLocalNumRows ();
154 if (offsets_.extent (0) < lclNumRows) {
155 offsets_ = offsets_type (); // clear first to save memory
156 offsets_ = offsets_type ("offsets", lclNumRows);
157 }
158 A_crs->getCrsGraph ()->getLocalDiagOffsets (offsets_);
159 }
160 }
161
162 isInitialized_ = true;
163 ++numInitialize_;
164}
165
166template<class MatrixType>
168{
169 // Either the matrix to precondition must be nonnull, or the user
170 // must have provided a Vector of diagonal entries.
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.");
176
177 if (! isInitialized_) {
178 initialize ();
179 }
180
181 // If the user did provide an input matrix, then that takes
182 // precedence over the vector of inverse diagonal entries, if they
183 // provided one earlier. This is only possible if they created this
184 // Diagonal instance using the constructor that takes a
185 // Tpetra::Vector pointer, and then called setMatrix() with a
186 // nonnull input matrix.
187 if (matrix_.is_null ()) { // accept the user's diagonal
188 inverseDiag_ = userInverseDiag_;
189 }
190 else {
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 ()) {
195 // Get the diagonal entries from the Tpetra::RowMatrix.
196 matrix_->getLocalDiagCopy (*tmpVec);
197 }
198 else {
199 // Get the diagonal entries from the Tpetra::CrsMatrix using the
200 // precomputed offsets.
201 A_crs->getLocalDiagCopy (*tmpVec, offsets_);
202 }
203 tmpVec->reciprocal (*tmpVec); // invert the diagonal entries
204 inverseDiag_ = tmpVec;
205 }
206
207 isComputed_ = true;
208 ++numCompute_;
209}
210
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,
215 Teuchos::ETransp /*mode*/,
216 scalar_type alpha,
217 scalar_type beta) const
218{
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().");
224
225 // FIXME (mfh 12 Sep 2014) This assumes that row Map == range Map ==
226 // domain Map. If the preconditioner has a matrix, we should ask
227 // the matrix whether we need to do an Import before and/or an
228 // Export after.
229
230 Y.elementWiseMultiply (alpha, *inverseDiag_, X, beta);
231 ++numApply_;
232}
233
234template <class MatrixType>
236 return numInitialize_;
237}
238
239template <class MatrixType>
241 return numCompute_;
242}
243
244template <class MatrixType>
246 return numApply_;
247}
248
249template <class MatrixType>
251 return initializeTime_;
252}
253
254template<class MatrixType>
256 return computeTime_;
257}
258
259template<class MatrixType>
261 return applyTime_;
262}
263
264template <class MatrixType>
266{
267 std::ostringstream out;
268
269 // Output is a valid YAML dictionary in flow style. If you don't
270 // like everything on a single line, you should call describe()
271 // instead.
272 out << "\"Ifpack2::Diagonal\": "
273 << "{";
274 if (this->getObjectLabel () != "") {
275 out << "Label: \"" << this->getObjectLabel () << "\", ";
276 }
277 if (matrix_.is_null ()) {
278 out << "Matrix: null";
279 }
280 else {
281 out << "Matrix: not null"
282 << ", Global matrix dimensions: ["
283 << matrix_->getGlobalNumRows () << ", "
284 << matrix_->getGlobalNumCols () << "]";
285 }
286
287 out << "}";
288 return out.str ();
289}
290
291template <class MatrixType>
293describe (Teuchos::FancyOStream &out,
294 const Teuchos::EVerbosityLevel verbLevel) const
295{
296 using std::endl;
297
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 () << "\", ";
308 }
309 out << "Number of initialize calls: " << numInitialize_ << endl
310 << "Number of compute calls: " << numCompute_ << endl
311 << "Number of apply calls: " << numApply_ << endl;
312 }
313}
314
315} // namespace Ifpack2
316
317#define IFPACK2_DIAGONAL_INSTANT(S,LO,GO,N) \
318 template class Ifpack2::Diagonal< Tpetra::RowMatrix<S, LO, GO, N> >;
319
320#endif
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 &params)
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