10 #ifndef IFPACK2_DENSECONTAINER_DEF_HPP
11 #define IFPACK2_DENSECONTAINER_DEF_HPP
13 #include "Tpetra_CrsMatrix.hpp"
15 #include "Tpetra_BlockMultiVector.hpp"
21 # include "Teuchos_DefaultSerialComm.hpp"
26 template<
class MatrixType,
class LocalScalarType>
32 ContainerImpl<MatrixType, LocalScalarType> (matrix, partitions, pointIndexed),
33 scalarOffsets_ (this->numBlocks_)
36 !matrix->hasColMap(), std::invalid_argument,
"Ifpack2::DenseContainer: "
37 "The constructor's input matrix must have a column Map.");
43 scalarOffsets_[i] = totalScalars;
47 scalars_.
resize(totalScalars);
52 diagBlocks_.emplace_back(
Teuchos::View, scalars_.
data() + scalarOffsets_[i], denseRows, denseRows, denseRows);
58 template<
class MatrixType,
class LocalScalarType>
62 template<
class MatrixType,
class LocalScalarType>
68 for(
int i = 0; i < this->numBlocks_; i++)
70 std::fill (ipiv_.begin (), ipiv_.end (), 0);
72 this->IsInitialized_ =
true;
75 this->IsComputed_ =
false;
78 template<
class MatrixType,
class LocalScalarType>
89 this->IsComputed_ =
false;
90 if (!this->isInitialized ()) {
96 this->IsComputed_ =
true;
99 template<
class MatrixType,
class LocalScalarType>
105 SC zero = STS::zero();
113 if(this->scalarsPerRow_ > 1)
115 Array<LO> colToBlockOffset(this->inputBlockMatrix_->getLocalNumCols(), INVALID);
116 for(
int i = 0; i < this->numBlocks_; i++)
119 LO blockStart = this->blockOffsets_[i];
120 LO blockEnd = blockStart + this->blockSizes_[i];
121 ArrayView<const LO> blockRows = this->getBlockRows(i);
125 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
127 LO localCol = this->translateRowToCol(blockRows[j]);
128 colToBlockOffset[localCol] = blockStart + j;
130 using h_inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
131 using h_vals_type =
typename block_crs_matrix_type::values_host_view_type;
132 for(LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++)
137 LO inputRow = this->blockRows_[blockStart + blockRow];
138 this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values);
139 LO numEntries = (LO) indices.size();
140 for(LO k = 0; k < numEntries; k++)
142 LO colOffset = colToBlockOffset[indices[k]];
143 if(blockStart <= colOffset && colOffset < blockEnd)
148 LO blockCol = colOffset - blockStart;
149 for(LO bc = 0; bc < this->bcrsBlockSize_; bc++)
151 for(LO br = 0; br < this->bcrsBlockSize_; br++)
153 LO r = this->bcrsBlockSize_ * blockRow + br;
154 LO c = this->bcrsBlockSize_ * blockCol + bc;
155 auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
157 diagBlocks_[i](r, c) = val;
169 Array<LO> colToBlockOffset(this->inputMatrix_->getLocalNumCols() * this->bcrsBlockSize_, INVALID);
170 for(
int i = 0; i < this->numBlocks_; i++)
173 LO blockStart = this->blockOffsets_[i];
174 LO blockEnd = blockStart + this->blockSizes_[i];
175 ArrayView<const LO> blockRows = this->getBlockRows(i);
179 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
182 LO localCol = this->translateRowToCol(blockRows[j]);
183 colToBlockOffset[localCol] = blockStart + j;
185 for(
size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++)
188 LO inputPointRow = this->blockRows_[blockStart + blockRow];
189 auto rowView = this->getInputRowView(inputPointRow);
190 for(
size_t k = 0; k < rowView.size(); k++)
192 LO colOffset = colToBlockOffset[rowView.ind(k)];
193 if(blockStart <= colOffset && colOffset < blockEnd)
195 LO blockCol = colOffset - blockStart;
196 auto val = rowView.val(k);
198 diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
207 template<
class MatrixType,
class LocalScalarType>
209 DenseContainer<MatrixType, LocalScalarType>::
213 for(
int i = 0; i < this->numBlocks_; i++)
216 int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
217 lapack.
GETRF(diagBlocks_[i].numRows(),
218 diagBlocks_[i].numCols(),
219 diagBlocks_[i].values(),
220 diagBlocks_[i].stride(),
224 INFO < 0, std::logic_error,
"Ifpack2::DenseContainer::factor: "
225 "LAPACK's _GETRF (LU factorization with partial pivoting) was called "
226 "incorrectly. INFO = " << INFO <<
" < 0. "
227 "Please report this bug to the Ifpack2 developers.");
232 INFO > 0, std::runtime_error,
"Ifpack2::DenseContainer::factor: "
233 "LAPACK's _GETRF (LU factorization with partial pivoting) reports that the "
234 "computed U factor is exactly singular. U(" << INFO <<
"," << INFO <<
") "
235 "(one-based index i) is exactly zero. This probably means that the input "
236 "matrix has a singular diagonal block.");
240 template<
class MatrixType,
class LocalScalarType>
242 DenseContainer<MatrixType, LocalScalarType>::
243 solveBlock(ConstHostSubviewLocal X,
250 #ifdef HAVE_IFPACK2_DEBUG
252 X.extent (0) != Y.extent (0),
253 std::logic_error,
"Ifpack2::DenseContainer::solveBlock: X and Y have "
254 "incompatible dimensions (" << X.extent (0) <<
" resp. "
255 << Y.extent (0) <<
"). Please report this bug to "
256 "the Ifpack2 developers.");
259 X.extent (1) != Y.extent(1),
260 std::logic_error,
"Ifpack2::DenseContainer::solveBlock: X and Y have "
261 "incompatible numbers of vectors (" << X.extent (1) <<
" resp. "
262 << Y.extent (1) <<
"). Please report this bug to "
263 "the Ifpack2 developers.");
267 size_t numRows = X.extent(0);
268 size_t numVecs = X.extent(1);
269 if(alpha == STS::zero()) {
270 if(beta == STS::zero()) {
274 for(
size_t j = 0; j < numVecs; j++)
276 for(
size_t i = 0; i < numRows; i++)
277 Y(i, j) = STS::zero();
281 for(
size_t j = 0; j < numVecs; j++)
283 for(
size_t i = 0; i < numRows; i++)
284 Y(i, j) = beta * Y(i, j);
292 std::vector<LSC> yTemp(numVecs * numRows);
293 for(
size_t j = 0; j < numVecs; j++)
295 for(
size_t i = 0; i < numRows; i++)
296 yTemp[j * numRows + i] = X(i, j);
300 int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
304 diagBlocks_[blockIndex].numRows (),
306 diagBlocks_[blockIndex].values (),
307 diagBlocks_[blockIndex].stride (),
313 if (beta != STS::zero ()) {
314 for(
size_t j = 0; j < numVecs; j++)
316 for(
size_t i = 0; i < numRows; i++)
318 Y(i, j) *= ISC(beta);
319 Y(i, j) += ISC(alpha * yTemp[j * numRows + i]);
324 for(
size_t j = 0; j < numVecs; j++)
326 for(
size_t i = 0; i < numRows; i++)
327 Y(i, j) = ISC(yTemp[j * numRows + i]);
332 INFO != 0, std::runtime_error,
"Ifpack2::DenseContainer::solveBlock: "
333 "LAPACK's _GETRS (solve using LU factorization with partial pivoting) "
334 "failed with INFO = " << INFO <<
" != 0.");
338 template<
class MatrixType,
class LocalScalarType>
345 this->describe (fos);
349 template<
class MatrixType,
class LocalScalarType>
354 std::ostringstream oss;
355 oss <<
"Ifpack::DenseContainer: ";
356 if (this->isInitialized()) {
357 if (this->isComputed()) {
358 oss <<
"{status = initialized, computed";
361 oss <<
"{status = initialized, not computed";
365 oss <<
"{status = not initialized, not computed";
372 template<
class MatrixType,
class LocalScalarType>
380 os <<
"================================================================================" << endl;
381 os <<
"Ifpack2::DenseContainer" << endl;
382 for(
int i = 0; i < this->numBlocks_; i++)
384 os <<
"Block " << i <<
" number of rows = " << this->blockSizes_[i] << endl;
386 os <<
"isInitialized() = " << this->IsInitialized_ << endl;
387 os <<
"isComputed() = " << this->IsComputed_ << endl;
388 os <<
"================================================================================" << endl;
392 template<
class MatrixType,
class LocalScalarType>
400 template<
class MatrixType,
class LocalScalarType>
412 #define IFPACK2_DENSECONTAINER_INSTANT(S,LO,GO,N) \
413 template class Ifpack2::DenseContainer< Tpetra::RowMatrix<S, LO, GO, N>, S >;
415 #endif // IFPACK2_DENSECONTAINER_DEF_HPP
Teuchos::Array< LO > blockRows_
Local indices of the rows of the input matrix that belong to this block.
Definition: Ifpack2_Container_decl.hpp:261
Teuchos::Array< LO > blockSizes_
Number of rows in each block.
Definition: Ifpack2_Container_decl.hpp:263
int numBlocks_
The number of blocks (partitions) in the container.
Definition: Ifpack2_Container_decl.hpp:259
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:310
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition: Ifpack2_DenseContainer_def.hpp:81
LO scalarsPerRow_
Definition: Ifpack2_Container_decl.hpp:284
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition: Ifpack2_DenseContainer_def.hpp:401
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:28
void GETRF(const OrdinalType &m, const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *IPIV, OrdinalType *info) const
Store and solve a local dense linear problem.
Definition: Ifpack2_DenseContainer_decl.hpp:70
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition: Ifpack2_DenseContainer_def.hpp:341
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:375
basic_FancyOStream & setOutputToRootOnly(const int rootRank)
void resize(size_type new_size, const value_type &x=value_type())
virtual ~DenseContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition: Ifpack2_DenseContainer_def.hpp:60
void GETRS(const char &TRANS, const OrdinalType &n, const OrdinalType &nrhs, const ScalarType *A, const OrdinalType &lda, const OrdinalType *IPIV, ScalarType *B, const OrdinalType &ldb, OrdinalType *info) const
virtual std::string description() const
A one-line description of this object.
Definition: Ifpack2_DenseContainer_def.hpp:352
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition: Ifpack2_DenseContainer_def.hpp:65