10#ifndef IFPACK2_BANDEDCONTAINER_DEF_HPP
11#define IFPACK2_BANDEDCONTAINER_DEF_HPP
13#include "Teuchos_LAPACK.hpp"
14#include "Tpetra_CrsMatrix.hpp"
20# include "Teuchos_DefaultMpiComm.hpp"
22# include "Teuchos_DefaultSerialComm.hpp"
27template<
class MatrixType,
class LocalScalarType>
30 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
31 const Teuchos::RCP<const import_type>&,
33 ContainerImpl<MatrixType, LocalScalarType>(matrix, partitions, pointIndexed),
39 TEUCHOS_TEST_FOR_EXCEPTION(
40 ! matrix->hasColMap (), std::invalid_argument,
"Ifpack2::BandedContainer: "
41 "The constructor's input matrix must have a column Map.");
44template<
class MatrixType,
class LocalScalarType>
48template<
class MatrixType,
class LocalScalarType>
52 if(List.isParameter(
"relaxation: banded container superdiagonals"))
54 int ku = List.get<
int>(
"relaxation: banded container superdiagonals");
58 if(List.isParameter(
"relaxation: banded container subdiagonals"))
60 int kl = List.get<
int>(
"relaxation: banded container subdiagonals");
66template<
class MatrixType,
class LocalScalarType>
71 using Teuchos::ArrayView;
73 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
74 size_t colToOffsetSize = this->
inputMatrix_->getLocalNumCols();
77 Array<LO> colToBlockOffset(colToOffsetSize, INVALID);
94 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
97 colToBlockOffset[localCol] = blockStart + j;
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++)
107 LO inputRow = this->
blockRows_[blockStart + blockRow];
109 LO numEntries = (LO) indices.size();
110 for(LO k = 0; k < numEntries; k++)
112 LO colOffset = colToBlockOffset[indices[k]];
113 if(blockStart <= colOffset && colOffset < blockEnd)
118 LO blockCol = colOffset - blockStart;
121 for(LO br = 0; br < this->bcrsBlockSize_; br++)
123 LO r = this->bcrsBlockSize_ * blockRow + br;
124 LO c = this->bcrsBlockSize_ * blockCol + bc;
144 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
148 colToBlockOffset[localCol] = blockStart + j;
150 for(LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++)
153 LO inputSplitRow = this->
blockRows_[blockStart + blockRow];
155 for(
size_t k = 0; k < rowView.size(); k++)
157 LO colOffset = colToBlockOffset[rowView.ind(k)];
158 if(blockStart <= colOffset && colOffset < blockEnd)
160 LO blockCol = colOffset - blockStart;
161 maxSub = std::max(maxSub, blockRow - blockCol);
162 maxSuper = std::max(maxSuper, blockCol - blockRow);
172template<
class MatrixType,
class LocalScalarType>
186 LO stride = 2 * kl_[b] + ku_[b] + 1;
187 scalarOffsets_[b] = totalScalars;
190 scalars_.resize(totalScalars);
191 for(
int b = 0; b < this->numBlocks_; b++)
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());
199 std::fill (ipiv_.begin (), ipiv_.end (), 0);
206template<
class MatrixType,
class LocalScalarType>
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.");
227template<
class MatrixType,
class LocalScalarType>
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();
241 if(this->scalarsPerRow_ > 1)
243 Array<LO> colToBlockOffset(this->inputBlockMatrix_->getLocalNumCols(), INVALID);
244 for(
int i = 0; i < this->numBlocks_; i++)
247 LO blockStart = this->blockOffsets_[i];
248 LO blockEnd = blockStart + this->blockSizes_[i];
249 ArrayView<const LO> blockRows = this->getBlockRows(i);
253 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
255 LO localCol = this->translateRowToCol(blockRows[j]);
256 colToBlockOffset[localCol] = blockStart + j;
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++)
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++)
270 LO colOffset = colToBlockOffset[indices[k]];
271 if(blockStart <= colOffset && colOffset < blockEnd)
276 LO blockCol = colOffset - blockStart;
277 for(LO bc = 0; bc < this->bcrsBlockSize_; bc++)
279 for(LO br = 0; br < this->bcrsBlockSize_; br++)
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)];
285 diagBlocks_[i](r, c) = val;
297 Array<LO> colToBlockOffset(this->inputMatrix_->getLocalNumCols() * this->bcrsBlockSize_, INVALID);
298 for(
int i = 0; i < this->numBlocks_; i++)
301 LO blockStart = this->blockOffsets_[i];
302 LO blockEnd = blockStart + this->blockSizes_[i];
303 ArrayView<const LO> blockRows = this->getBlockRows(i);
307 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
310 LO localCol = this->translateRowToCol(blockRows[j]);
311 colToBlockOffset[localCol] = blockStart + j;
313 for(
size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++)
316 LO inputPointRow = this->blockRows_[blockStart + blockRow];
317 auto rowView = this->getInputRowView(inputPointRow);
318 for(
size_t k = 0; k < rowView.size(); k++)
320 LO colOffset = colToBlockOffset[rowView.ind(k)];
321 if(blockStart <= colOffset && colOffset < blockEnd)
323 LO blockCol = colOffset - blockStart;
324 auto val = rowView.val(k);
326 diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
334template<
class MatrixType,
class LocalScalarType>
341 ContainerImpl<MatrixType, LocalScalarType>::clearBlocks();
344template<
class MatrixType,
class LocalScalarType>
349 Teuchos::LAPACK<int, LSC> lapack;
352 for(
int i = 0; i < this->numBlocks_; i++)
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(),
364 diagBlocks_[i].values(),
365 diagBlocks_[i].stride(),
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.");
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.");
387template<
class MatrixType,
class LocalScalarType>
393 Teuchos::ETransp mode,
395 const LSC beta)
const
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.");
420 size_t numRows = X.extent (0);
421 size_t numVecs = X.extent (1);
423 LSC zero = Teuchos::ScalarTraits<LSC>::zero ();
429 for(
size_t j = 0; j < numVecs; j++)
430 for(
size_t i = 0; i < numRows; i++)
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);
440 Teuchos::LAPACK<int, LSC> lapack;
444 std::vector<LSC> yTemp(numVecs * numRows);
445 for(
size_t j = 0; j < numVecs; j++)
447 for(
size_t i = 0; i < numRows; i++)
448 yTemp[j * numRows + i] = X(i, j);
452 const char trans = (mode == Teuchos::CONJ_TRANS ?
'C' : (mode == Teuchos::TRANS ?
'T' :
'N'));
454 const int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
457 diagBlocks_[blockIndex].numCols(),
458 diagBlocks_[blockIndex].lowerBandwidth(),
460 diagBlocks_[blockIndex].upperBandwidth() - diagBlocks_[blockIndex].lowerBandwidth(),
462 diagBlocks_[blockIndex].values(),
463 diagBlocks_[blockIndex].stride(),
470 for(
size_t j = 0; j < numVecs; j++)
472 for(
size_t i = 0; i < numRows; i++)
474 Y(i, j) *= ISC(beta);
475 Y(i, j) += ISC(alpha * yTemp[j * numRows + i]);
480 for(
size_t j = 0; j < numVecs; j++)
482 for(
size_t i = 0; i < numRows; i++)
483 Y(i, j) = ISC(yTemp[j * numRows + i]);
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.");
494template<
class MatrixType,
class LocalScalarType>
497print (std::ostream& os)
const
499 Teuchos::FancyOStream fos (Teuchos::rcpFromRef (os));
500 fos.setOutputToRootOnly (0);
505template<
class MatrixType,
class LocalScalarType>
510 std::ostringstream oss;
511 oss << Teuchos::Describable::description();
514 oss <<
"{status = initialized, computed";
517 oss <<
"{status = initialized, not computed";
521 oss <<
"{status = not initialized, not computed";
527template<
class MatrixType,
class LocalScalarType>
531 const Teuchos::EVerbosityLevel verbLevel)
const
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++)
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;
543 os <<
"isComputed() = " << this->
IsComputed_ << std::endl;
544 os <<
"================================================================================" << std::endl;
548template<
class MatrixType,
class LocalScalarType>
556#define IFPACK2_BANDEDCONTAINER_INSTANT(S,LO,GO,N) \
557 template class Ifpack2::BandedContainer< Tpetra::RowMatrix<S, LO, GO, N>, S >;
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