Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_BandedContainer_def.hpp
1 /*@HEADER
2 // ***********************************************************************
3 //
4 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
5 // Copyright (2009) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ***********************************************************************
40 //@HEADER
41 */
42 
43 #ifndef IFPACK2_BANDEDCONTAINER_DEF_HPP
44 #define IFPACK2_BANDEDCONTAINER_DEF_HPP
45 
46 #include "Teuchos_LAPACK.hpp"
47 #include "Tpetra_CrsMatrix.hpp"
48 #include <iostream>
49 #include <sstream>
50 
51 #ifdef HAVE_MPI
52 # include <mpi.h>
54 #else
55 # include "Teuchos_DefaultSerialComm.hpp"
56 #endif // HAVE_MPI
57 
58 namespace Ifpack2 {
59 
60 template<class MatrixType, class LocalScalarType>
63  const Teuchos::Array<Teuchos::Array<LO> >& partitions,
65  bool pointIndexed) :
66  ContainerImpl<MatrixType, LocalScalarType>(matrix, partitions, pointIndexed),
67  ipiv_(this->blockRows_.size() * this->scalarsPerRow_),
68  kl_(this->numBlocks_, -1),
69  ku_(this->numBlocks_, -1),
70  scalarOffsets_(this->numBlocks_)
71 {
73  ! matrix->hasColMap (), std::invalid_argument, "Ifpack2::BandedContainer: "
74  "The constructor's input matrix must have a column Map.");
75 }
76 
77 template<class MatrixType, class LocalScalarType>
80 
81 template<class MatrixType, class LocalScalarType>
84 {
85  if(List.isParameter("relaxation: banded container superdiagonals"))
86  {
87  int ku = List.get<int>("relaxation: banded container superdiagonals");
88  for(LO b = 0; b < this->numBlocks_; b++)
89  ku_[b] = ku;
90  }
91  if(List.isParameter("relaxation: banded container subdiagonals"))
92  {
93  int kl = List.get<int>("relaxation: banded container subdiagonals");
94  for(LO b = 0; b < this->numBlocks_; b++)
95  kl_[b] = kl;
96  }
97 }
98 
99 template<class MatrixType, class LocalScalarType>
102 {
103  using Teuchos::Array;
104  using Teuchos::ArrayView;
105  //now, for any block where kl_ or ku_ has not already been set, compute the actual bandwidth
106  const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
107  size_t colToOffsetSize = this->inputMatrix_->getNodeNumCols();
108  if(this->pointIndexed_)
109  colToOffsetSize *= this->bcrsBlockSize_;
110  Array<LO> colToBlockOffset(colToOffsetSize, INVALID);
111  //Same logic as extract() to find entries in blocks efficiently
112  //(but here, just to find bandwidth, no scalars used)
113  for(int i = 0; i < this->numBlocks_; i++)
114  {
115  //maxSub, maxSuper are the maximum lower and upper bandwidth
116  LO maxSub = 0;
117  LO maxSuper = 0;
118  if(this->scalarsPerRow_ > 1)
119  {
120  //Get the interval where block i is defined in blockRows_
121  LO blockStart = this->blockOffsets_[i];
122  LO blockEnd = (i == this->numBlocks_ - 1) ? this->blockRows_.size() : this->blockOffsets_[i + 1];
123  ArrayView<const LO> blockRows = this->getBlockRows(i);
124  //Set the lookup table entries for the columns appearing in block i.
125  //If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
126  //this is OK. The values updated here are only needed to process block i's entries.
127  for(size_t j = 0; j < size_t(blockRows.size()); j++)
128  {
129  LO localCol = this->translateRowToCol(blockRows[j]);
130  colToBlockOffset[localCol] = blockStart + j;
131  }
132  for(LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++)
133  {
134  //get a raw view of the whole block row
135  const LO* indices;
136  SC* values;
137  LO numEntries;
138  LO inputRow = this->blockRows_[blockStart + blockRow];
139  this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values, numEntries);
140  for(LO k = 0; k < numEntries; k++)
141  {
142  LO colOffset = colToBlockOffset[indices[k]];
143  if(blockStart <= colOffset && colOffset < blockEnd)
144  {
145  //This entry does appear in the diagonal block.
146  //(br, bc) identifies the scalar's position in the BlockCrs block.
147  //Convert this to (r, c) which is its position in the container block.
148  LO blockCol = colOffset - blockStart;
149  for(LO bc = 0; bc < this->bcrsBlockSize_; bc++)
150  {
151  for(LO br = 0; br < this->bcrsBlockSize_; br++)
152  {
153  LO r = this->bcrsBlockSize_ * blockRow + br;
154  LO c = this->bcrsBlockSize_ * blockCol + bc;
155  if(r - c > maxSub)
156  maxSub = r - c;
157  if(c - r > maxSuper)
158  maxSuper = c - r;
159  }
160  }
161  }
162  }
163  }
164  }
165  else
166  {
167  //Get the interval where block i is defined in blockRows_
168  LO blockStart = this->blockOffsets_[i];
169  LO blockEnd = (i == this->numBlocks_ - 1) ? this->blockRows_.size() : this->blockOffsets_[i + 1];
170  ArrayView<const LO> blockRows = this->getBlockRows(i);
171  //Set the lookup table entries for the columns appearing in block i.
172  //If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
173  //this is OK. The values updated here are only needed to process block i's entries.
174  for(size_t j = 0; j < size_t(blockRows.size()); j++)
175  {
176  //translateRowToCol will return the corresponding split column
177  LO localCol = this->translateRowToCol(blockRows[j]);
178  colToBlockOffset[localCol] = blockStart + j;
179  }
180  for(LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++)
181  {
182  //get a view of the general row
183  LO inputSplitRow = this->blockRows_[blockStart + blockRow];
184  auto rowView = this->getInputRowView(inputSplitRow);
185  for(size_t k = 0; k < rowView.size(); k++)
186  {
187  LO colOffset = colToBlockOffset[rowView.ind(k)];
188  if(blockStart <= colOffset && colOffset < blockEnd)
189  {
190  LO blockCol = colOffset - blockStart;
191  maxSub = std::max(maxSub, blockRow - blockCol);
192  maxSuper = std::max(maxSuper, blockCol - blockRow);
193  }
194  }
195  }
196  }
197  kl_[i] = maxSub;
198  ku_[i] = maxSuper;
199  }
200 }
201 
202 template<class MatrixType, class LocalScalarType>
203 void
206 {
207  //note: in BlockRelaxation, Container_->setParameters() immediately
208  //precedes Container_->initialize(). So no matter what, either all
209  //the block bandwidths were set (in setParameters()) or none of them were.
210  //If none were they must be computed individually.
211  if(kl_[0] == -1)
212  computeBandwidth();
213  GO totalScalars = 0;
214  for(LO b = 0; b < this->numBlocks_; b++)
215  {
216  LO stride = 2 * kl_[b] + ku_[b] + 1;
217  scalarOffsets_[b] = totalScalars;
218  totalScalars += stride * this->blockSizes_[b] * this->scalarsPerRow_;
219  }
220  scalars_.resize(totalScalars);
221  for(int b = 0; b < this->numBlocks_; b++)
222  {
223  //NOTE: the stride and upper bandwidth used to construct the SerialBandDenseMatrix looks
224  //too large, but the extra kl_ in upper band space is needed by the LAPACK factorization routine.
225  LO nrows = this->blockSizes_[b] * this->scalarsPerRow_;
226  diagBlocks_.emplace_back(Teuchos::View, scalars_.data() + scalarOffsets_[b], 2 * kl_[b] + ku_[b] + 1, nrows, nrows, kl_[b], kl_[b] + ku_[b]);
227  diagBlocks_[b].putScalar(Teuchos::ScalarTraits<LSC>::zero());
228  }
229  std::fill (ipiv_.begin (), ipiv_.end (), 0);
230  // We assume that if you called this method, you intend to recompute
231  // everything.
232  this->IsComputed_ = false;
233  this->IsInitialized_ = true;
234 }
235 
236 template<class MatrixType, class LocalScalarType>
237 void
240 {
242  ipiv_.size () != this->blockRows_.size() * this->scalarsPerRow_, std::logic_error,
243  "Ifpack2::BandedContainer::compute: ipiv_ array has the wrong size. "
244  "Please report this bug to the Ifpack2 developers.");
245 
246  this->IsComputed_ = false;
247  if (!this->isInitialized ()) {
248  this->initialize ();
249  }
250 
251  extract (); // Extract the submatrices
252  factor (); // Factor them
253 
254  this->IsComputed_ = true;
255 }
256 
257 template<class MatrixType, class LocalScalarType>
259 {
260  using Teuchos::Array;
261  using Teuchos::ArrayView;
262  using STS = Teuchos::ScalarTraits<SC>;
263  SC zero = STS::zero();
264  const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
265  //To extract diagonal blocks, need to translate local rows to local columns.
266  //Strategy: make a lookup table that translates local cols in the matrix to offsets in blockRows_:
267  //blockOffsets_[b] <= offset < blockOffsets_[b+1]: tests whether the column is in block b.
268  //offset - blockOffsets_[b]: gives the column within block b.
269  //
270  //This provides the block and col within a block in O(1).
271  if(this->scalarsPerRow_ > 1)
272  {
273  Array<LO> colToBlockOffset(this->inputBlockMatrix_->getNodeNumCols(), INVALID);
274  for(int i = 0; i < this->numBlocks_; i++)
275  {
276  //Get the interval where block i is defined in blockRows_
277  LO blockStart = this->blockOffsets_[i];
278  LO blockEnd = blockStart + this->blockSizes_[i];
279  ArrayView<const LO> blockRows = this->getBlockRows(i);
280  //Set the lookup table entries for the columns appearing in block i.
281  //If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
282  //this is OK. The values updated here are only needed to process block i's entries.
283  for(size_t j = 0; j < size_t(blockRows.size()); j++)
284  {
285  LO localCol = this->translateRowToCol(blockRows[j]);
286  colToBlockOffset[localCol] = blockStart + j;
287  }
288  for(LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++)
289  {
290  //get a raw view of the whole block row
291  const LO* indices;
292  SC* values;
293  LO numEntries;
294  LO inputRow = this->blockRows_[blockStart + blockRow];
295  this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values, numEntries);
296  for(LO k = 0; k < numEntries; k++)
297  {
298  LO colOffset = colToBlockOffset[indices[k]];
299  if(blockStart <= colOffset && colOffset < blockEnd)
300  {
301  //This entry does appear in the diagonal block.
302  //(br, bc) identifies the scalar's position in the BlockCrs block.
303  //Convert this to (r, c) which is its position in the container block.
304  LO blockCol = colOffset - blockStart;
305  for(LO bc = 0; bc < this->bcrsBlockSize_; bc++)
306  {
307  for(LO br = 0; br < this->bcrsBlockSize_; br++)
308  {
309  LO r = this->bcrsBlockSize_ * blockRow + br;
310  LO c = this->bcrsBlockSize_ * blockCol + bc;
311  auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
312  if(val != zero)
313  diagBlocks_[i](r, c) = val;
314  }
315  }
316  }
317  }
318  }
319  }
320  }
321  else
322  {
323  //get the mapping from point-indexed matrix columns to offsets in blockRows_
324  //(this includes regular CrsMatrix columns, in which case bcrsBlockSize_ == 1)
325  Array<LO> colToBlockOffset(this->inputMatrix_->getNodeNumCols() * this->bcrsBlockSize_, INVALID);
326  for(int i = 0; i < this->numBlocks_; i++)
327  {
328  //Get the interval where block i is defined in blockRows_
329  LO blockStart = this->blockOffsets_[i];
330  LO blockEnd = blockStart + this->blockSizes_[i];
331  ArrayView<const LO> blockRows = this->getBlockRows(i);
332  //Set the lookup table entries for the columns appearing in block i.
333  //If OverlapLevel_ > 0, then this may overwrite values for previous blocks, but
334  //this is OK. The values updated here are only needed to process block i's entries.
335  for(size_t j = 0; j < size_t(blockRows.size()); j++)
336  {
337  //translateRowToCol will return the corresponding split column
338  LO localCol = this->translateRowToCol(blockRows[j]);
339  colToBlockOffset[localCol] = blockStart + j;
340  }
341  for(size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++)
342  {
343  //get a view of the split row
344  LO inputPointRow = this->blockRows_[blockStart + blockRow];
345  auto rowView = this->getInputRowView(inputPointRow);
346  for(size_t k = 0; k < rowView.size(); k++)
347  {
348  LO colOffset = colToBlockOffset[rowView.ind(k)];
349  if(blockStart <= colOffset && colOffset < blockEnd)
350  {
351  LO blockCol = colOffset - blockStart;
352  auto val = rowView.val(k);
353  if(val != zero)
354  diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
355  }
356  }
357  }
358  }
359  }
360 }
361 
362 template<class MatrixType, class LocalScalarType>
363 void
364 BandedContainer<MatrixType, LocalScalarType>::
365 clearBlocks ()
366 {
367  diagBlocks_.clear();
368  scalars_.clear();
369  ContainerImpl<MatrixType, LocalScalarType>::clearBlocks();
370 }
371 
372 template<class MatrixType, class LocalScalarType>
373 void
374 BandedContainer<MatrixType, LocalScalarType>::
375 factor ()
376 {
378  int INFO = 0;
379 
380  for(int i = 0; i < this->numBlocks_; i++)
381  {
382  // Plausibility checks for matrix
383  TEUCHOS_TEST_FOR_EXCEPTION(diagBlocks_[i].values()==0, std::invalid_argument,
384  "BandedContainer<T>::factor: Diagonal block is an empty SerialBandDenseMatrix<T>!");
385  TEUCHOS_TEST_FOR_EXCEPTION(diagBlocks_[i].upperBandwidth() < diagBlocks_[i].lowerBandwidth(), std::invalid_argument,
386  "BandedContainer<T>::factor: Diagonal block needs kl additional superdiagonals for factorization! However, the number of superdiagonals is smaller than the number of subdiagonals!");
387  int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
388  lapack.GBTRF (diagBlocks_[i].numRows(),
389  diagBlocks_[i].numCols(),
390  diagBlocks_[i].lowerBandwidth(),
391  diagBlocks_[i].upperBandwidth() - diagBlocks_[i].lowerBandwidth(), /* enter the real number of superdiagonals (see Teuchos_SerialBandDenseSolver)*/
392  diagBlocks_[i].values(),
393  diagBlocks_[i].stride(),
394  blockIpiv,
395  &INFO);
396 
397  // INFO < 0 is a bug.
399  INFO < 0, std::logic_error, "Ifpack2::BandedContainer::factor: "
400  "LAPACK's _GBTRF (LU factorization with partial pivoting) was called "
401  "incorrectly. INFO = " << INFO << " < 0. "
402  "Please report this bug to the Ifpack2 developers.");
403  // INFO > 0 means the matrix is singular. This is probably an issue
404  // either with the choice of rows the rows we extracted, or with the
405  // input matrix itself.
407  INFO > 0, std::runtime_error, "Ifpack2::BandedContainer::factor: "
408  "LAPACK's _GBTRF (LU factorization with partial pivoting) reports that the "
409  "computed U factor is exactly singular. U(" << INFO << "," << INFO << ") "
410  "(one-based index i) is exactly zero. This probably means that the input "
411  "matrix has a singular diagonal block.");
412  }
413 }
414 
415 template<class MatrixType, class LocalScalarType>
416 void
417 BandedContainer<MatrixType, LocalScalarType>::
418 solveBlock(HostSubviewLocal X,
419  HostSubviewLocal Y,
420  int blockIndex,
421  Teuchos::ETransp mode,
422  const LSC alpha,
423  const LSC beta) const
424 {
425  #ifdef HAVE_IFPACK2_DEBUG
427  X.extent (0) != Y.extent (0),
428  std::logic_error, "Ifpack2::BandedContainer::solveBlock: X and Y have "
429  "incompatible dimensions (" << X.extent (0) << " resp. "
430  << Y.extent (0) << "). Please report this bug to "
431  "the Ifpack2 developers.");
433  X.extent (0) != static_cast<size_t> (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numCols() : diagBlocks_[blockIndex].numRows()),
434  std::logic_error, "Ifpack2::BandedContainer::solveBlock: The input "
435  "multivector X has incompatible dimensions from those of the "
436  "inverse operator (" << X.extent (0) << " vs. "
437  << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numCols() : diagBlocks_[blockIndex].numRows())
438  << "). Please report this bug to the Ifpack2 developers.");
440  Y.extent (0) != static_cast<size_t> (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRows() : diagBlocks_[blockIndex].numCols()),
441  std::logic_error, "Ifpack2::BandedContainer::solveBlock: The output "
442  "multivector Y has incompatible dimensions from those of the "
443  "inverse operator (" << Y.extent (0) << " vs. "
444  << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRows() : diagBlocks_[blockIndex].numCols())
445  << "). Please report this bug to the Ifpack2 developers.");
446  #endif
447 
448  size_t numRows = X.extent (0);
449  size_t numVecs = X.extent (1);
450 
451  LSC zero = Teuchos::ScalarTraits<LSC>::zero ();
452  if (alpha == zero) { // don't need to solve the linear system
453  if (beta == zero) {
454  // Use BLAS AXPY semantics for beta == 0: overwrite, clobbering
455  // any Inf or NaN values in Y (rather than multiplying them by
456  // zero, resulting in NaN values).
457  for(size_t j = 0; j < numVecs; j++)
458  for(size_t i = 0; i < numRows; i++)
459  Y(i, j) = zero;
460  }
461  else { // beta != 0
462  for(size_t j = 0; j < numVecs; j++)
463  for(size_t i = 0; i < numRows; i++)
464  Y(i, j) = beta * Y(i, j);
465  }
466  }
467  else { // alpha != 0; must solve the linear system
469  // If beta is nonzero or Y is not constant stride, we have to use
470  // a temporary output multivector. It gets a copy of X, since
471  // GBTRS overwrites its (multi)vector input with its output.
472  std::vector<LSC> yTemp(numVecs * numRows);
473  for(size_t j = 0; j < numVecs; j++)
474  {
475  for(size_t i = 0; i < numRows; i++)
476  yTemp[j * numRows + i] = X(i, j);
477  }
478 
479  int INFO = 0;
480  const char trans = (mode == Teuchos::CONJ_TRANS ? 'C' : (mode == Teuchos::TRANS ? 'T' : 'N'));
481 
482  const int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
483 
484  lapack.GBTRS(trans,
485  diagBlocks_[blockIndex].numCols(),
486  diagBlocks_[blockIndex].lowerBandwidth(),
487  /* enter the real number of superdiagonals (see Teuchos_SerialBandDenseSolver)*/
488  diagBlocks_[blockIndex].upperBandwidth() - diagBlocks_[blockIndex].lowerBandwidth(),
489  numVecs,
490  diagBlocks_[blockIndex].values(),
491  diagBlocks_[blockIndex].stride(),
492  blockIpiv,
493  yTemp.data(),
494  numRows,
495  &INFO);
496 
497  if (beta != zero) {
498  for(size_t j = 0; j < numVecs; j++)
499  {
500  for(size_t i = 0; i < numRows; i++)
501  {
502  Y(i, j) *= ISC(beta);
503  Y(i, j) += ISC(alpha * yTemp[j * numRows + i]);
504  }
505  }
506  }
507  else {
508  for(size_t j = 0; j < numVecs; j++)
509  {
510  for(size_t i = 0; i < numRows; i++)
511  Y(i, j) = ISC(yTemp[j * numRows + i]);
512  }
513  }
514 
516  INFO != 0, std::runtime_error, "Ifpack2::BandedContainer::solveBlock: "
517  "LAPACK's _GBTRS (solve using LU factorization with partial pivoting) "
518  "failed with INFO = " << INFO << " != 0.");
519  }
520 }
521 
522 template<class MatrixType, class LocalScalarType>
523 std::ostream&
525 print (std::ostream& os) const
526 {
527  Teuchos::FancyOStream fos (Teuchos::rcpFromRef (os));
528  fos.setOutputToRootOnly (0);
529  describe (fos);
530  return os;
531 }
532 
533 template<class MatrixType, class LocalScalarType>
534 std::string
536 description () const
537 {
538  std::ostringstream oss;
540  if (this->isInitialized()) {
541  if (this->isComputed()) {
542  oss << "{status = initialized, computed";
543  }
544  else {
545  oss << "{status = initialized, not computed";
546  }
547  }
548  else {
549  oss << "{status = not initialized, not computed";
550  }
551  oss << "}";
552  return oss.str();
553 }
554 
555 template<class MatrixType, class LocalScalarType>
556 void
559  const Teuchos::EVerbosityLevel verbLevel) const
560 {
561  if(verbLevel==Teuchos::VERB_NONE) return;
562  os << "================================================================================" << std::endl;
563  os << "Ifpack2::BandedContainer" << std::endl;
564  for(int i = 0; i < this->numBlocks_; i++)
565  {
566  os << "Block " << i << ": Number of rows = " << this->blockSizes_[i] << std::endl;
567  os << "Block " << i << ": Number of subdiagonals = " << diagBlocks_[i].lowerBandwidth() << std::endl;
568  os << "Block " << i << ": Number of superdiagonals = " << diagBlocks_[i].upperBandwidth() << std::endl;
569  }
570  os << "isInitialized() = " << this->IsInitialized_ << std::endl;
571  os << "isComputed() = " << this->IsComputed_ << std::endl;
572  os << "================================================================================" << std::endl;
573  os << std::endl;
574 }
575 
576 template<class MatrixType, class LocalScalarType>
578 {
579  return "Banded";
580 }
581 
582 } // namespace Ifpack2
583 
584 #define IFPACK2_BANDEDCONTAINER_INSTANT(S,LO,GO,N) \
585  template class Ifpack2::BandedContainer< Tpetra::RowMatrix<S, LO, GO, N>, S >;
586 
587 #endif // IFPACK2_BANDEDCONTAINER_HPP
virtual void initialize() override
Do all set-up operations that only require matrix structure.
Definition: Ifpack2_BandedContainer_def.hpp:205
T & get(const std::string &name, T def_value)
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS). This class allows a custom scalar type (LocalScalarType) to be used for storing blocks and solving the block systems. Hiding this template parameter from the Container interface simplifies the BlockRelaxation and ContainerFactory classes.
Definition: Ifpack2_Container_decl.hpp:342
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void GBTRF(const OrdinalType &m, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const
virtual void compute() override
Initialize and compute each block.
Definition: Ifpack2_BandedContainer_def.hpp:239
bool isParameter(const std::string &name) const
virtual std::string description() const
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition: Ifpack2_BandedContainer_def.hpp:577
virtual std::ostream & print(std::ostream &os) const override
Print information about this object to the given output stream.
Definition: Ifpack2_BandedContainer_def.hpp:525
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
virtual ~BandedContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition: Ifpack2_BandedContainer_def.hpp:79
void computeBandwidth()
Definition: Ifpack2_BandedContainer_def.hpp:101
void GBTRS(const char &TRANS, const OrdinalType &n, const OrdinalType &kl, const OrdinalType &ku, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
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:62
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:558
virtual void setParameters(const Teuchos::ParameterList &List) override
Set all necessary parameters (super- and sub-diagonal bandwidth)
Definition: Ifpack2_BandedContainer_def.hpp:83
Store and solve a local Banded linear problem.
Definition: Ifpack2_BandedContainer_decl.hpp:102
virtual std::string description() const override
A one-line description of this object.
Definition: Ifpack2_BandedContainer_def.hpp:536