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