Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_BandedContainer_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_BANDEDCONTAINER_DEF_HPP
11#define IFPACK2_BANDEDCONTAINER_DEF_HPP
12
13#include "Teuchos_LAPACK.hpp"
14#include "Tpetra_CrsMatrix.hpp"
16#include <iostream>
17#include <sstream>
18
19#ifdef HAVE_MPI
20#include <mpi.h>
21#include "Teuchos_DefaultMpiComm.hpp"
22#else
23#include "Teuchos_DefaultSerialComm.hpp"
24#endif // HAVE_MPI
25
26namespace Ifpack2 {
27
28template <class MatrixType, class LocalScalarType>
30 BandedContainer(const Teuchos::RCP<const row_matrix_type>& matrix,
31 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
32 const Teuchos::RCP<const import_type>&,
33 bool pointIndexed)
34 : ContainerImpl<MatrixType, LocalScalarType>(matrix, partitions, pointIndexed)
35 , ipiv_(this->blockRows_.size() * this->scalarsPerRow_)
36 , kl_(this->numBlocks_, -1)
37 , ku_(this->numBlocks_, -1)
38 , scalarOffsets_(this->numBlocks_) {
39 TEUCHOS_TEST_FOR_EXCEPTION(
40 !matrix->hasColMap(), std::invalid_argument,
41 "Ifpack2::BandedContainer: "
42 "The constructor's input matrix must have a column Map.");
43}
44
45template <class MatrixType, class LocalScalarType>
48
49template <class MatrixType, class LocalScalarType>
51 setParameters(const Teuchos::ParameterList& List) {
52 if (List.isParameter("relaxation: banded container superdiagonals")) {
53 int ku = List.get<int>("relaxation: banded container superdiagonals");
54 for (LO b = 0; b < this->numBlocks_; b++)
55 ku_[b] = ku;
56 }
57 if (List.isParameter("relaxation: banded container subdiagonals")) {
58 int kl = List.get<int>("relaxation: banded container subdiagonals");
59 for (LO b = 0; b < this->numBlocks_; b++)
60 kl_[b] = kl;
61 }
62}
63
64template <class MatrixType, class LocalScalarType>
67 using Teuchos::Array;
68 using Teuchos::ArrayView;
69 // now, for any block where kl_ or ku_ has not already been set, compute the actual bandwidth
70 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
71 size_t colToOffsetSize = this->inputMatrix_->getLocalNumCols();
72 if (this->pointIndexed_)
73 colToOffsetSize *= this->bcrsBlockSize_;
74 Array<LO> colToBlockOffset(colToOffsetSize, INVALID);
75 // Same logic as extract() to find entries in blocks efficiently
76 //(but here, just to find bandwidth, no scalars used)
77 for (int i = 0; i < this->numBlocks_; i++) {
78 // maxSub, maxSuper are the maximum lower and upper bandwidth
79 LO maxSub = 0;
80 LO maxSuper = 0;
81 if (this->scalarsPerRow_ > 1) {
82 // Get the interval where block i is defined in blockRows_
83 LO blockStart = this->blockOffsets_[i];
84 LO blockEnd = (i == this->numBlocks_ - 1) ? this->blockRows_.size() : this->blockOffsets_[i + 1];
85 ArrayView<const LO> blockRows = this->getBlockRows(i);
86 // Set the lookup table entries for the columns appearing in block i.
87 // If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
88 // this is OK. The values updated here are only needed to process block i's entries.
89 for (size_t j = 0; j < size_t(blockRows.size()); j++) {
90 LO localCol = this->translateRowToCol(blockRows[j]);
91 colToBlockOffset[localCol] = blockStart + j;
92 }
93
94 using h_inds_type = typename block_crs_matrix_type::local_inds_host_view_type;
95 using h_vals_type = typename block_crs_matrix_type::values_host_view_type;
96 for (LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++) {
97 // get a raw view of the whole block row
98 h_inds_type indices;
99 h_vals_type values;
100 LO inputRow = this->blockRows_[blockStart + blockRow];
101 this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values);
102 LO numEntries = (LO)indices.size();
103 for (LO k = 0; k < numEntries; k++) {
104 LO colOffset = colToBlockOffset[indices[k]];
105 if (blockStart <= colOffset && colOffset < blockEnd) {
106 // This entry does appear in the diagonal block.
107 //(br, bc) identifies the scalar's position in the BlockCrs block.
108 // Convert this to (r, c) which is its position in the container block.
109 LO blockCol = colOffset - blockStart;
110 for (LO bc = 0; bc < this->bcrsBlockSize_; bc++) {
111 for (LO br = 0; br < this->bcrsBlockSize_; br++) {
112 LO r = this->bcrsBlockSize_ * blockRow + br;
113 LO c = this->bcrsBlockSize_ * blockCol + bc;
114 if (r - c > maxSub)
115 maxSub = r - c;
116 if (c - r > maxSuper)
117 maxSuper = c - r;
118 }
119 }
120 }
121 }
122 }
123 } else {
124 // Get the interval where block i is defined in blockRows_
125 LO blockStart = this->blockOffsets_[i];
126 LO blockEnd = (i == this->numBlocks_ - 1) ? this->blockRows_.size() : this->blockOffsets_[i + 1];
127 ArrayView<const LO> blockRows = this->getBlockRows(i);
128 // Set the lookup table entries for the columns appearing in block i.
129 // If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
130 // this is OK. The values updated here are only needed to process block i's entries.
131 for (size_t j = 0; j < size_t(blockRows.size()); j++) {
132 // translateRowToCol will return the corresponding split column
133 LO localCol = this->translateRowToCol(blockRows[j]);
134 colToBlockOffset[localCol] = blockStart + j;
135 }
136 for (LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++) {
137 // get a view of the general row
138 LO inputSplitRow = this->blockRows_[blockStart + blockRow];
139 auto rowView = this->getInputRowView(inputSplitRow);
140 for (size_t k = 0; k < rowView.size(); k++) {
141 LO colOffset = colToBlockOffset[rowView.ind(k)];
142 if (blockStart <= colOffset && colOffset < blockEnd) {
143 LO blockCol = colOffset - blockStart;
144 maxSub = std::max(maxSub, blockRow - blockCol);
145 maxSuper = std::max(maxSuper, blockCol - blockRow);
146 }
147 }
148 }
149 }
150 kl_[i] = maxSub;
151 ku_[i] = maxSuper;
152 }
153}
154
155template <class MatrixType, class LocalScalarType>
157 initialize() {
158 // note: in BlockRelaxation, Container_->setParameters() immediately
159 // precedes Container_->initialize(). So no matter what, either all
160 // the block bandwidths were set (in setParameters()) or none of them were.
161 // If none were they must be computed individually.
162 if (kl_[0] == -1)
163 computeBandwidth();
164 GO totalScalars = 0;
165 for (LO b = 0; b < this->numBlocks_; b++) {
166 LO stride = 2 * kl_[b] + ku_[b] + 1;
167 scalarOffsets_[b] = totalScalars;
168 totalScalars += stride * this->blockSizes_[b] * this->scalarsPerRow_;
169 }
170 scalars_.resize(totalScalars);
171 for (int b = 0; b < this->numBlocks_; b++) {
172 // NOTE: the stride and upper bandwidth used to construct the SerialBandDenseMatrix looks
173 // too large, but the extra kl_ in upper band space is needed by the LAPACK factorization routine.
174 LO nrows = this->blockSizes_[b] * this->scalarsPerRow_;
175 diagBlocks_.emplace_back(Teuchos::View, scalars_.data() + scalarOffsets_[b], 2 * kl_[b] + ku_[b] + 1, nrows, nrows, kl_[b], kl_[b] + ku_[b]);
176 diagBlocks_[b].putScalar(Teuchos::ScalarTraits<LSC>::zero());
177 }
178 std::fill(ipiv_.begin(), ipiv_.end(), 0);
179 // We assume that if you called this method, you intend to recompute
180 // everything.
181 this->IsComputed_ = false;
182 this->IsInitialized_ = true;
183}
184
185template <class MatrixType, class LocalScalarType>
187 compute() {
188 TEUCHOS_TEST_FOR_EXCEPTION(
189 ipiv_.size() != this->blockRows_.size() * this->scalarsPerRow_, std::logic_error,
190 "Ifpack2::BandedContainer::compute: ipiv_ array has the wrong size. "
191 "Please report this bug to the Ifpack2 developers.");
192
193 this->IsComputed_ = false;
194 if (!this->isInitialized()) {
195 this->initialize();
196 }
197
198 extract(); // Extract the submatrices
199 factor(); // Factor them
200
201 this->IsComputed_ = true;
202}
203
204template <class MatrixType, class LocalScalarType>
206 using Teuchos::Array;
207 using Teuchos::ArrayView;
208 using STS = Teuchos::ScalarTraits<SC>;
209 SC zero = STS::zero();
210 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
211 // To extract diagonal blocks, need to translate local rows to local columns.
212 // Strategy: make a lookup table that translates local cols in the matrix to offsets in blockRows_:
213 // blockOffsets_[b] <= offset < blockOffsets_[b+1]: tests whether the column is in block b.
214 // offset - blockOffsets_[b]: gives the column within block b.
215 //
216 // This provides the block and col within a block in O(1).
217 if (this->scalarsPerRow_ > 1) {
218 Array<LO> colToBlockOffset(this->inputBlockMatrix_->getLocalNumCols(), INVALID);
219 for (int i = 0; i < this->numBlocks_; i++) {
220 // Get the interval where block i is defined in blockRows_
221 LO blockStart = this->blockOffsets_[i];
222 LO blockEnd = blockStart + this->blockSizes_[i];
223 ArrayView<const LO> blockRows = this->getBlockRows(i);
224 // Set the lookup table entries for the columns appearing in block i.
225 // If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
226 // this is OK. The values updated here are only needed to process block i's entries.
227 for (size_t j = 0; j < size_t(blockRows.size()); j++) {
228 LO localCol = this->translateRowToCol(blockRows[j]);
229 colToBlockOffset[localCol] = blockStart + j;
230 }
231 using h_inds_type = typename block_crs_matrix_type::local_inds_host_view_type;
232 using h_vals_type = typename block_crs_matrix_type::values_host_view_type;
233 for (LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++) {
234 // get a raw view of the whole block row
235 h_inds_type indices;
236 h_vals_type values;
237 LO inputRow = this->blockRows_[blockStart + blockRow];
238 this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values);
239 LO numEntries = (LO)indices.size();
240 for (LO k = 0; k < numEntries; k++) {
241 LO colOffset = colToBlockOffset[indices[k]];
242 if (blockStart <= colOffset && colOffset < blockEnd) {
243 // This entry does appear in the diagonal block.
244 //(br, bc) identifies the scalar's position in the BlockCrs block.
245 // Convert this to (r, c) which is its position in the container block.
246 LO blockCol = colOffset - blockStart;
247 for (LO bc = 0; bc < this->bcrsBlockSize_; bc++) {
248 for (LO br = 0; br < this->bcrsBlockSize_; br++) {
249 LO r = this->bcrsBlockSize_ * blockRow + br;
250 LO c = this->bcrsBlockSize_ * blockCol + bc;
251 auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
252 if (val != zero)
253 diagBlocks_[i](r, c) = val;
254 }
255 }
256 }
257 }
258 }
259 }
260 } else {
261 // get the mapping from point-indexed matrix columns to offsets in blockRows_
262 //(this includes regular CrsMatrix columns, in which case bcrsBlockSize_ == 1)
263 Array<LO> colToBlockOffset(this->inputMatrix_->getLocalNumCols() * this->bcrsBlockSize_, INVALID);
264 for (int i = 0; i < this->numBlocks_; i++) {
265 // Get the interval where block i is defined in blockRows_
266 LO blockStart = this->blockOffsets_[i];
267 LO blockEnd = blockStart + this->blockSizes_[i];
268 ArrayView<const LO> blockRows = this->getBlockRows(i);
269 // Set the lookup table entries for the columns appearing in block i.
270 // If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
271 // this is OK. The values updated here are only needed to process block i's entries.
272 for (size_t j = 0; j < size_t(blockRows.size()); j++) {
273 // translateRowToCol will return the corresponding split column
274 LO localCol = this->translateRowToCol(blockRows[j]);
275 colToBlockOffset[localCol] = blockStart + j;
276 }
277 for (size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++) {
278 // get a view of the split row
279 LO inputPointRow = this->blockRows_[blockStart + blockRow];
280 auto rowView = this->getInputRowView(inputPointRow);
281 for (size_t k = 0; k < rowView.size(); k++) {
282 LO colOffset = colToBlockOffset[rowView.ind(k)];
283 if (blockStart <= colOffset && colOffset < blockEnd) {
284 LO blockCol = colOffset - blockStart;
285 auto val = rowView.val(k);
286 if (val != zero)
287 diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
288 }
289 }
290 }
291 }
292 }
293}
294
295template <class MatrixType, class LocalScalarType>
296void BandedContainer<MatrixType, LocalScalarType>::
297 clearBlocks() {
298 diagBlocks_.clear();
299 scalars_.clear();
300 ContainerImpl<MatrixType, LocalScalarType>::clearBlocks();
301}
302
303template <class MatrixType, class LocalScalarType>
304void BandedContainer<MatrixType, LocalScalarType>::
305 factor() {
306 Teuchos::LAPACK<int, LSC> lapack;
307 int INFO = 0;
308
309 for (int i = 0; i < this->numBlocks_; i++) {
310 // Plausibility checks for matrix
311 TEUCHOS_TEST_FOR_EXCEPTION(diagBlocks_[i].values() == 0, std::invalid_argument,
312 "BandedContainer<T>::factor: Diagonal block is an empty SerialBandDenseMatrix<T>!");
313 TEUCHOS_TEST_FOR_EXCEPTION(diagBlocks_[i].upperBandwidth() < diagBlocks_[i].lowerBandwidth(), std::invalid_argument,
314 "BandedContainer<T>::factor: Diagonal block needs kl additional superdiagonals for factorization! However, the number of superdiagonals is smaller than the number of subdiagonals!");
315 int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
316 lapack.GBTRF(diagBlocks_[i].numRows(),
317 diagBlocks_[i].numCols(),
318 diagBlocks_[i].lowerBandwidth(),
319 diagBlocks_[i].upperBandwidth() - diagBlocks_[i].lowerBandwidth(), /* enter the real number of superdiagonals (see Teuchos_SerialBandDenseSolver)*/
320 diagBlocks_[i].values(),
321 diagBlocks_[i].stride(),
322 blockIpiv,
323 &INFO);
324
325 // INFO < 0 is a bug.
326 TEUCHOS_TEST_FOR_EXCEPTION(
327 INFO < 0, std::logic_error,
328 "Ifpack2::BandedContainer::factor: "
329 "LAPACK's _GBTRF (LU factorization with partial pivoting) was called "
330 "incorrectly. INFO = "
331 << INFO << " < 0. "
332 "Please report this bug to the Ifpack2 developers.");
333 // INFO > 0 means the matrix is singular. This is probably an issue
334 // either with the choice of rows the rows we extracted, or with the
335 // input matrix itself.
336 TEUCHOS_TEST_FOR_EXCEPTION(
337 INFO > 0, std::runtime_error,
338 "Ifpack2::BandedContainer::factor: "
339 "LAPACK's _GBTRF (LU factorization with partial pivoting) reports that the "
340 "computed U factor is exactly singular. U("
341 << INFO << "," << INFO << ") "
342 "(one-based index i) is exactly zero. This probably means that the input "
343 "matrix has a singular diagonal block.");
344 }
345}
346
347template <class MatrixType, class LocalScalarType>
348void BandedContainer<MatrixType, LocalScalarType>::
349 solveBlock(ConstHostSubviewLocal X,
350 HostSubviewLocal Y,
351 int blockIndex,
352 Teuchos::ETransp mode,
353 const LSC alpha,
354 const LSC beta) const {
356 TEUCHOS_TEST_FOR_EXCEPTION(
357 X.extent(0) != Y.extent(0),
358 std::logic_error,
359 "Ifpack2::BandedContainer::solveBlock: X and Y have "
360 "incompatible dimensions ("
361 << X.extent(0) << " resp. "
362 << Y.extent(0) << "). Please report this bug to "
363 "the Ifpack2 developers.");
364 TEUCHOS_TEST_FOR_EXCEPTION(
365 X.extent(0) != static_cast<size_t>(mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numCols() : diagBlocks_[blockIndex].numRows()),
366 std::logic_error,
367 "Ifpack2::BandedContainer::solveBlock: The input "
368 "multivector X has incompatible dimensions from those of the "
369 "inverse operator ("
370 << X.extent(0) << " vs. "
371 << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numCols() : diagBlocks_[blockIndex].numRows())
372 << "). Please report this bug to the Ifpack2 developers.");
373 TEUCHOS_TEST_FOR_EXCEPTION(
374 Y.extent(0) != static_cast<size_t>(mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRows() : diagBlocks_[blockIndex].numCols()),
375 std::logic_error,
376 "Ifpack2::BandedContainer::solveBlock: The output "
377 "multivector Y has incompatible dimensions from those of the "
378 "inverse operator ("
379 << Y.extent(0) << " vs. "
380 << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRows() : diagBlocks_[blockIndex].numCols())
381 << "). Please report this bug to the Ifpack2 developers.");
382 }
383
384 size_t numRows = X.extent(0);
385 size_t numVecs = X.extent(1);
386
387 LSC zero = Teuchos::ScalarTraits<LSC>::zero();
388 if (alpha == zero) { // don't need to solve the linear system
389 if (beta == zero) {
390 // Use BLAS AXPY semantics for beta == 0: overwrite, clobbering
391 // any Inf or NaN values in Y (rather than multiplying them by
392 // zero, resulting in NaN values).
393 for (size_t j = 0; j < numVecs; j++)
394 for (size_t i = 0; i < numRows; i++)
395 Y(i, j) = zero;
396 } else { // beta != 0
397 for (size_t j = 0; j < numVecs; j++)
398 for (size_t i = 0; i < numRows; i++)
399 Y(i, j) = beta * Y(i, j);
400 }
401 } else { // alpha != 0; must solve the linear system
402 Teuchos::LAPACK<int, LSC> lapack;
403 // If beta is nonzero or Y is not constant stride, we have to use
404 // a temporary output multivector. It gets a copy of X, since
405 // GBTRS overwrites its (multi)vector input with its output.
406 std::vector<LSC> yTemp(numVecs * numRows);
407 for (size_t j = 0; j < numVecs; j++) {
408 for (size_t i = 0; i < numRows; i++)
409 yTemp[j * numRows + i] = X(i, j);
410 }
411
412 int INFO = 0;
413 const char trans = (mode == Teuchos::CONJ_TRANS ? 'C' : (mode == Teuchos::TRANS ? 'T' : 'N'));
414
415 const int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
416
417 lapack.GBTRS(trans,
418 diagBlocks_[blockIndex].numCols(),
419 diagBlocks_[blockIndex].lowerBandwidth(),
420 /* enter the real number of superdiagonals (see Teuchos_SerialBandDenseSolver)*/
421 diagBlocks_[blockIndex].upperBandwidth() - diagBlocks_[blockIndex].lowerBandwidth(),
422 numVecs,
423 diagBlocks_[blockIndex].values(),
424 diagBlocks_[blockIndex].stride(),
425 blockIpiv,
426 yTemp.data(),
427 numRows,
428 &INFO);
429
430 if (beta != zero) {
431 for (size_t j = 0; j < numVecs; j++) {
432 for (size_t i = 0; i < numRows; i++) {
433 Y(i, j) *= ISC(beta);
434 Y(i, j) += ISC(alpha * yTemp[j * numRows + i]);
435 }
436 }
437 } else {
438 for (size_t j = 0; j < numVecs; j++) {
439 for (size_t i = 0; i < numRows; i++)
440 Y(i, j) = ISC(yTemp[j * numRows + i]);
441 }
442 }
443
444 TEUCHOS_TEST_FOR_EXCEPTION(
445 INFO != 0, std::runtime_error,
446 "Ifpack2::BandedContainer::solveBlock: "
447 "LAPACK's _GBTRS (solve using LU factorization with partial pivoting) "
448 "failed with INFO = "
449 << INFO << " != 0.");
450 }
451}
452
453template <class MatrixType, class LocalScalarType>
454std::ostream&
456 print(std::ostream& os) const {
457 Teuchos::FancyOStream fos(Teuchos::rcpFromRef(os));
458 fos.setOutputToRootOnly(0);
459 describe(fos);
460 return os;
461}
462
463template <class MatrixType, class LocalScalarType>
464std::string
466 description() const {
467 std::ostringstream oss;
468 oss << Teuchos::Describable::description();
469 if (this->isInitialized()) {
470 if (this->isComputed()) {
471 oss << "{status = initialized, computed";
472 } else {
473 oss << "{status = initialized, not computed";
474 }
475 } else {
476 oss << "{status = not initialized, not computed";
477 }
478 oss << "}";
479 return oss.str();
480}
481
482template <class MatrixType, class LocalScalarType>
484 describe(Teuchos::FancyOStream& os,
485 const Teuchos::EVerbosityLevel verbLevel) const {
486 if (verbLevel == Teuchos::VERB_NONE) return;
487 os << "================================================================================" << std::endl;
488 os << "Ifpack2::BandedContainer" << std::endl;
489 for (int i = 0; i < this->numBlocks_; i++) {
490 os << "Block " << i << ": Number of rows = " << this->blockSizes_[i] << std::endl;
491 os << "Block " << i << ": Number of subdiagonals = " << diagBlocks_[i].lowerBandwidth() << std::endl;
492 os << "Block " << i << ": Number of superdiagonals = " << diagBlocks_[i].upperBandwidth() << std::endl;
493 }
494 os << "isInitialized() = " << this->IsInitialized_ << std::endl;
495 os << "isComputed() = " << this->IsComputed_ << std::endl;
496 os << "================================================================================" << std::endl;
497 os << std::endl;
498}
499
500template <class MatrixType, class LocalScalarType>
502 return "Banded";
503}
504
505} // namespace Ifpack2
506
507#define IFPACK2_BANDEDCONTAINER_INSTANT(S, LO, GO, N) \
508 template class Ifpack2::BandedContainer<Tpetra::RowMatrix<S, LO, GO, N>, S>;
509
510#endif // IFPACK2_BANDEDCONTAINER_HPP
Declaration of Ifpack2::Details::Behavior, a class that describes Ifpack2's run-time behavior.
Store and solve a local Banded linear problem.
Definition Ifpack2_BandedContainer_decl.hpp:70
BandedContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, const Teuchos::RCP< const import_type > &, bool pointIndexed)
Constructor.
Definition Ifpack2_BandedContainer_def.hpp:30
virtual void compute() override
Initialize and compute each block.
Definition Ifpack2_BandedContainer_def.hpp:187
virtual void initialize() override
Do all set-up operations that only require matrix structure.
Definition Ifpack2_BandedContainer_def.hpp:157
virtual std::ostream & print(std::ostream &os) const override
Print information about this object to the given output stream.
Definition Ifpack2_BandedContainer_def.hpp:456
virtual ~BandedContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition Ifpack2_BandedContainer_def.hpp:47
virtual void setParameters(const Teuchos::ParameterList &List) override
Set all necessary parameters (super- and sub-diagonal bandwidth)
Definition Ifpack2_BandedContainer_def.hpp:51
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition Ifpack2_BandedContainer_def.hpp:501
virtual std::string description() const override
A one-line description of this object.
Definition Ifpack2_BandedContainer_def.hpp:466
void computeBandwidth()
Definition Ifpack2_BandedContainer_def.hpp:66
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print the object with some verbosity level to the given FancyOStream.
Definition Ifpack2_BandedContainer_def.hpp:484
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS)....
Definition Ifpack2_Container_decl.hpp:307
static bool debug()
Whether Ifpack2 is in debug mode.
Definition Ifpack2_Details_Behavior.cpp:31
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40