10#ifndef IFPACK2_DENSECONTAINER_DEF_HPP
11#define IFPACK2_DENSECONTAINER_DEF_HPP
13#include "Tpetra_CrsMatrix.hpp"
14#include "Teuchos_LAPACK.hpp"
15#include "Tpetra_BlockMultiVector.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)
34 , scalarOffsets_(this->numBlocks_) {
35 TEUCHOS_TEST_FOR_EXCEPTION(
36 !matrix->hasColMap(), std::invalid_argument,
37 "Ifpack2::DenseContainer: "
38 "The constructor's input matrix must have a column Map.");
43 scalarOffsets_[i] = totalScalars;
47 scalars_.resize(totalScalars);
51 diagBlocks_.emplace_back(Teuchos::View, scalars_.data() + scalarOffsets_[i], denseRows, denseRows, denseRows);
54 ipiv_.resize(this->
blockRows_.size() * this->scalarsPerRow_);
57template <
class MatrixType,
class LocalScalarType>
61template <
class MatrixType,
class LocalScalarType>
65 for (
int i = 0; i < this->numBlocks_; i++)
66 diagBlocks_[i].putScalar(Teuchos::ScalarTraits<LSC>::zero());
67 std::fill(ipiv_.begin(), ipiv_.end(), 0);
69 this->IsInitialized_ =
true;
72 this->IsComputed_ =
false;
75template <
class MatrixType,
class LocalScalarType>
84 this->IsComputed_ =
false;
85 if (!this->isInitialized()) {
91 this->IsComputed_ =
true;
94template <
class MatrixType,
class LocalScalarType>
97 using Teuchos::ArrayView;
98 using STS = Teuchos::ScalarTraits<SC>;
99 SC zero = STS::zero();
100 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
107 if (this->scalarsPerRow_ > 1) {
108 Array<LO> colToBlockOffset(this->inputBlockMatrix_->getLocalNumCols(), INVALID);
109 for (
int i = 0; i < this->numBlocks_; i++) {
111 LO blockStart = this->blockOffsets_[i];
112 LO blockEnd = blockStart + this->blockSizes_[i];
113 ArrayView<const LO> blockRows = this->getBlockRows(i);
117 for (
size_t j = 0; j < size_t(blockRows.size()); j++) {
118 LO localCol = this->translateRowToCol(blockRows[j]);
119 colToBlockOffset[localCol] = blockStart + j;
121 using h_inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
122 using h_vals_type =
typename block_crs_matrix_type::values_host_view_type;
123 for (LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++) {
127 LO inputRow = this->blockRows_[blockStart + blockRow];
128 this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values);
129 LO numEntries = (LO)indices.size();
130 for (LO k = 0; k < numEntries; k++) {
131 LO colOffset = colToBlockOffset[indices[k]];
132 if (blockStart <= colOffset && colOffset < blockEnd) {
136 LO blockCol = colOffset - blockStart;
137 for (LO bc = 0; bc < this->bcrsBlockSize_; bc++) {
138 for (LO br = 0; br < this->bcrsBlockSize_; br++) {
139 LO r = this->bcrsBlockSize_ * blockRow + br;
140 LO c = this->bcrsBlockSize_ * blockCol + bc;
141 auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
143 diagBlocks_[i](r, c) = val;
153 Array<LO> colToBlockOffset(this->inputMatrix_->getLocalNumCols() * this->bcrsBlockSize_, INVALID);
154 for (
int i = 0; i < this->numBlocks_; i++) {
156 LO blockStart = this->blockOffsets_[i];
157 LO blockEnd = blockStart + this->blockSizes_[i];
158 ArrayView<const LO> blockRows = this->getBlockRows(i);
162 for (
size_t j = 0; j < size_t(blockRows.size()); j++) {
164 LO localCol = this->translateRowToCol(blockRows[j]);
165 colToBlockOffset[localCol] = blockStart + j;
167 for (
size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++) {
169 LO inputPointRow = this->blockRows_[blockStart + blockRow];
170 auto rowView = this->getInputRowView(inputPointRow);
171 for (
size_t k = 0; k < rowView.size(); k++) {
172 LO colOffset = colToBlockOffset[rowView.ind(k)];
173 if (blockStart <= colOffset && colOffset < blockEnd) {
174 LO blockCol = colOffset - blockStart;
175 auto val = rowView.val(k);
177 diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
185template <
class MatrixType,
class LocalScalarType>
186void DenseContainer<MatrixType, LocalScalarType>::
188 Teuchos::LAPACK<int, LSC> lapack;
189 for (
int i = 0; i < this->numBlocks_; i++) {
191 int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
192 lapack.GETRF(diagBlocks_[i].numRows(),
193 diagBlocks_[i].numCols(),
194 diagBlocks_[i].values(),
195 diagBlocks_[i].stride(),
198 TEUCHOS_TEST_FOR_EXCEPTION(
199 INFO < 0, std::logic_error,
200 "Ifpack2::DenseContainer::factor: "
201 "LAPACK's _GETRF (LU factorization with partial pivoting) was called "
202 "incorrectly. INFO = "
204 "Please report this bug to the Ifpack2 developers.");
208 TEUCHOS_TEST_FOR_EXCEPTION(
209 INFO > 0, std::runtime_error,
210 "Ifpack2::DenseContainer::factor: "
211 "LAPACK's _GETRF (LU factorization with partial pivoting) reports that the "
212 "computed U factor is exactly singular. U("
213 << INFO <<
"," << INFO <<
") "
214 "(one-based index i) is exactly zero. This probably means that the input "
215 "matrix has a singular diagonal block.");
219template <
class MatrixType,
class LocalScalarType>
220void DenseContainer<MatrixType, LocalScalarType>::
221 solveBlock(ConstHostSubviewLocal X,
224 Teuchos::ETransp mode,
228 TEUCHOS_TEST_FOR_EXCEPTION(
229 X.extent(0) != Y.extent(0),
231 "Ifpack2::DenseContainer::solveBlock: X and Y have "
232 "incompatible dimensions ("
233 << X.extent(0) <<
" resp. "
234 << Y.extent(0) <<
"). Please report this bug to "
235 "the Ifpack2 developers.");
237 TEUCHOS_TEST_FOR_EXCEPTION(
238 X.extent(1) != Y.extent(1),
240 "Ifpack2::DenseContainer::solveBlock: X and Y have "
241 "incompatible numbers of vectors ("
242 << X.extent(1) <<
" resp. "
243 << Y.extent(1) <<
"). Please report this bug to "
244 "the Ifpack2 developers.");
247 typedef Teuchos::ScalarTraits<LSC> STS;
248 size_t numRows = X.extent(0);
249 size_t numVecs = X.extent(1);
250 if (alpha == STS::zero()) {
251 if (beta == STS::zero()) {
255 for (
size_t j = 0; j < numVecs; j++) {
256 for (
size_t i = 0; i < numRows; i++)
257 Y(i, j) = STS::zero();
260 for (
size_t j = 0; j < numVecs; j++) {
261 for (
size_t i = 0; i < numRows; i++)
262 Y(i, j) = beta * Y(i, j);
265 Teuchos::LAPACK<int, LSC> lapack;
269 std::vector<LSC> yTemp(numVecs * numRows);
270 for (
size_t j = 0; j < numVecs; j++) {
271 for (
size_t i = 0; i < numRows; i++)
272 yTemp[j * numRows + i] = X(i, j);
276 int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
278 (mode == Teuchos::CONJ_TRANS ?
'C' : (mode == Teuchos::TRANS ?
'T' :
'N'));
280 diagBlocks_[blockIndex].numRows(),
282 diagBlocks_[blockIndex].values(),
283 diagBlocks_[blockIndex].stride(),
289 if (beta != STS::zero()) {
290 for (
size_t j = 0; j < numVecs; j++) {
291 for (
size_t i = 0; i < numRows; i++) {
292 Y(i, j) *= ISC(beta);
293 Y(i, j) += ISC(alpha * yTemp[j * numRows + i]);
297 for (
size_t j = 0; j < numVecs; j++) {
298 for (
size_t i = 0; i < numRows; i++)
299 Y(i, j) = ISC(yTemp[j * numRows + i]);
303 TEUCHOS_TEST_FOR_EXCEPTION(
304 INFO != 0, std::runtime_error,
305 "Ifpack2::DenseContainer::solveBlock: "
306 "LAPACK's _GETRS (solve using LU factorization with partial pivoting) "
307 "failed with INFO = "
308 << INFO <<
" != 0.");
312template <
class MatrixType,
class LocalScalarType>
315 print(std::ostream& os)
const {
316 Teuchos::FancyOStream fos(Teuchos::rcpFromRef(os));
317 fos.setOutputToRootOnly(0);
322template <
class MatrixType,
class LocalScalarType>
326 std::ostringstream oss;
327 oss <<
"Ifpack::DenseContainer: ";
328 if (this->isInitialized()) {
329 if (this->isComputed()) {
330 oss <<
"{status = initialized, computed";
332 oss <<
"{status = initialized, not computed";
335 oss <<
"{status = not initialized, not computed";
342template <
class MatrixType,
class LocalScalarType>
345 const Teuchos::EVerbosityLevel verbLevel)
const {
347 if (verbLevel == Teuchos::VERB_NONE)
return;
348 os <<
"================================================================================" << endl;
349 os <<
"Ifpack2::DenseContainer" << endl;
350 for (
int i = 0; i < this->numBlocks_; i++) {
351 os <<
"Block " << i <<
" number of rows = " << this->blockSizes_[i] << endl;
353 os <<
"isInitialized() = " << this->IsInitialized_ << endl;
354 os <<
"isComputed() = " << this->IsComputed_ << endl;
355 os <<
"================================================================================" << endl;
359template <
class MatrixType,
class LocalScalarType>
366template <
class MatrixType,
class LocalScalarType>
377#define IFPACK2_DENSECONTAINER_INSTANT(S, LO, GO, N) \
378 template class Ifpack2::DenseContainer<Tpetra::RowMatrix<S, LO, GO, N>, S>;
Declaration of Ifpack2::Details::Behavior, a class that describes Ifpack2's run-time behavior.
int numBlocks_
The number of blocks (partitions) in the container.
Definition Ifpack2_Container_decl.hpp:257
LO scalarsPerRow_
Definition Ifpack2_Container_decl.hpp:282
Teuchos::Array< LO > blockSizes_
Number of rows in each block.
Definition Ifpack2_Container_decl.hpp:261
Teuchos::Array< LO > blockRows_
Local indices of the rows of the input matrix that belong to this block.
Definition Ifpack2_Container_decl.hpp:259
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS)....
Definition Ifpack2_Container_decl.hpp:307
Store and solve a local dense linear problem.
Definition Ifpack2_DenseContainer_decl.hpp:71
virtual std::string description() const
A one-line description of this object.
Definition Ifpack2_DenseContainer_def.hpp:325
virtual void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to the given FancyOStream.
Definition Ifpack2_DenseContainer_def.hpp:344
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition Ifpack2_DenseContainer_def.hpp:367
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition Ifpack2_DenseContainer_def.hpp:315
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition Ifpack2_DenseContainer_def.hpp:77
DenseContainer(const Teuchos::RCP< const row_matrix_type > &matrix, const Teuchos::Array< Teuchos::Array< LO > > &partitions, const Teuchos::RCP< const import_type > &importer, bool pointIndexed)
Constructor.
Definition Ifpack2_DenseContainer_def.hpp:29
virtual ~DenseContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition Ifpack2_DenseContainer_def.hpp:59
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition Ifpack2_DenseContainer_def.hpp:63
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