10#ifndef IFPACK2_TRIDICONTAINER_DEF_HPP
11#define IFPACK2_TRIDICONTAINER_DEF_HPP
13#include "Teuchos_LAPACK.hpp"
18#include "Teuchos_DefaultMpiComm.hpp"
20#include "Teuchos_DefaultSerialComm.hpp"
25template <
class MatrixType,
class LocalScalarType>
28 const Teuchos::Array<Teuchos::Array<LO> >& partitions,
29 const Teuchos::RCP<const import_type>&,
31 :
ContainerImpl<MatrixType, LocalScalarType>(matrix, partitions, pointIndexed)
32 , ipiv_(this->blockRows_.size() * this->scalarsPerRow_)
33 , scalarOffsets_(this->numBlocks_) {
34 TEUCHOS_TEST_FOR_EXCEPTION(
35 !matrix->hasColMap(), std::invalid_argument,
36 "Ifpack2::TriDiContainer: "
37 "The constructor's input matrix must have a column Map.");
43 scalarOffsets_[i] = scalarTotal;
46 if (actualBlockRows == 1) {
51 scalarTotal += 4 * (actualBlockRows - 1);
55 scalars_.resize(scalarTotal);
56 diagBlocks_.reserve(this->numBlocks_);
58 diagBlocks_.emplace_back(Teuchos::View, scalars_.data() + scalarOffsets_[i], this->
blockSizes_[i] * this->
scalarsPerRow_);
62template <
class MatrixType,
class LocalScalarType>
66template <
class MatrixType,
class LocalScalarType>
68 for (
int i = 0; i < this->numBlocks_; i++)
69 diagBlocks_[i].putScalar(Teuchos::ScalarTraits<LSC>::zero());
70 this->IsInitialized_ =
true;
73 this->IsComputed_ =
false;
76template <
class MatrixType,
class LocalScalarType>
79 TEUCHOS_TEST_FOR_EXCEPTION(
80 ipiv_.size() != this->blockRows_.size() * this->scalarsPerRow_, std::logic_error,
81 "Ifpack2::TriDiContainer::compute: ipiv_ array has the wrong size. "
82 "Please report this bug to the Ifpack2 developers.");
85 this->IsComputed_ =
false;
86 if (!this->isInitialized()) {
91 this->IsComputed_ =
true;
94template <
class MatrixType,
class LocalScalarType>
97 using Teuchos::ArrayView;
98 SC zero = Teuchos::ScalarTraits<SC>::zero();
99 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
106 if (this->scalarsPerRow_ > 1) {
107 Array<LO> colToBlockOffset(this->inputBlockMatrix_->getLocalNumCols(), INVALID);
108 for (
int i = 0; i < this->numBlocks_; i++) {
110 LO blockStart = this->blockOffsets_[i];
111 LO blockEnd = blockStart + this->blockSizes_[i];
112 ArrayView<const LO> blockRows = this->getBlockRows(i);
116 for (
size_t j = 0; j < size_t(blockRows.size()); j++) {
117 LO localCol = this->translateRowToCol(blockRows[j]);
118 colToBlockOffset[localCol] = blockStart + j;
120 using h_inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
121 using h_vals_type =
typename block_crs_matrix_type::values_host_view_type;
122 for (LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++) {
126 LO inputRow = this->blockRows_[blockStart + blockRow];
127 this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values);
128 LO numEntries = (LO)indices.size();
129 for (LO k = 0; k < numEntries; k++) {
130 LO colOffset = colToBlockOffset[indices[k]];
131 if (blockStart <= colOffset && colOffset < blockEnd) {
135 LO blockCol = colOffset - blockStart;
136 for (LO bc = 0; bc < this->bcrsBlockSize_; bc++) {
137 for (LO br = 0; br < this->bcrsBlockSize_; br++) {
138 LO r = this->bcrsBlockSize_ * blockRow + br;
139 LO c = this->bcrsBlockSize_ * blockCol + bc;
140 auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
142 diagBlocks_[i](r, c) = val;
152 Array<LO> colToBlockOffset(this->inputMatrix_->getLocalNumCols() * this->bcrsBlockSize_, INVALID);
153 for (
int i = 0; i < this->numBlocks_; i++) {
155 LO blockStart = this->blockOffsets_[i];
156 LO blockEnd = blockStart + this->blockSizes_[i];
157 ArrayView<const LO> blockRows = this->getBlockRows(i);
161 for (
size_t j = 0; j < size_t(blockRows.size()); j++) {
163 LO localCol = this->translateRowToCol(blockRows[j]);
164 colToBlockOffset[localCol] = blockStart + j;
166 for (
size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++) {
168 LO inputPointRow = this->blockRows_[blockStart + blockRow];
169 auto rowView = this->getInputRowView(inputPointRow);
170 for (
size_t k = 0; k < rowView.size(); k++) {
171 LO colOffset = colToBlockOffset[rowView.ind(k)];
172 if (blockStart <= colOffset && colOffset < blockEnd) {
173 LO blockCol = colOffset - blockStart;
174 auto val = rowView.val(k);
176 diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
184template <
class MatrixType,
class LocalScalarType>
185void TriDiContainer<MatrixType, LocalScalarType>::clearBlocks() {
188 ContainerImpl<MatrixType, LocalScalarType>::clearBlocks();
191template <
class MatrixType,
class LocalScalarType>
192void TriDiContainer<MatrixType, LocalScalarType>::factor() {
193 for (
int i = 0; i < this->numBlocks_; i++) {
194 Teuchos::LAPACK<int, LSC> lapack;
196 int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
197 lapack.GTTRF(diagBlocks_[i].numRowsCols(),
201 diagBlocks_[i].DU2(),
204 TEUCHOS_TEST_FOR_EXCEPTION(
205 INFO < 0, std::logic_error,
206 "Ifpack2::TriDiContainer::factor: "
207 "LAPACK's _GTTRF (LU factorization with partial pivoting) was called "
208 "incorrectly. INFO = "
210 "Please report this bug to the Ifpack2 developers.");
214 TEUCHOS_TEST_FOR_EXCEPTION(
215 INFO > 0, std::runtime_error,
216 "Ifpack2::TriDiContainer::factor: "
217 "LAPACK's _GTTRF (LU factorization with partial pivoting) reports that the "
218 "computed U factor is exactly singular. U("
219 << INFO <<
"," << INFO <<
") "
220 "(one-based index i) is exactly zero. This probably means that the input "
221 "matrix has a singular diagonal block.");
225template <
class MatrixType,
class LocalScalarType>
230 Teuchos::ETransp mode,
233 LSC zero = Teuchos::ScalarTraits<LSC>::zero();
234 size_t numVecs = X.extent(1);
235 size_t numRows = X.extent(0);
238 TEUCHOS_TEST_FOR_EXCEPTION(
239 X.extent(0) != Y.extent(0),
241 "Ifpack2::TriDiContainer::solveBlock: X and Y have "
242 "incompatible dimensions ("
243 << X.extent(0) <<
" resp. "
244 << Y.extent(0) <<
"). Please report this bug to "
245 "the Ifpack2 developers.");
246 TEUCHOS_TEST_FOR_EXCEPTION(
247 X.extent(0) !=
static_cast<size_t>(diagBlocks_[blockIndex].numRowsCols()),
249 "Ifpack2::TriDiContainer::solveBlock: The input "
250 "multivector X has incompatible dimensions from those of the "
252 << X.extent(0) <<
" vs. "
253 << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRowsCols() : diagBlocks_[blockIndex].numRowsCols())
254 <<
"). Please report this bug to the Ifpack2 developers.");
255 TEUCHOS_TEST_FOR_EXCEPTION(
256 Y.extent(0) !=
static_cast<size_t>(diagBlocks_[blockIndex].numRowsCols()),
258 "Ifpack2::TriDiContainer::solveBlock: The output "
259 "multivector Y has incompatible dimensions from those of the "
261 << Y.extent(0) <<
" vs. "
262 << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRowsCols() : diagBlocks_[blockIndex].numRowsCols())
263 <<
"). Please report this bug to the Ifpack2 developers.");
271 for (
size_t j = 0; j < numVecs; j++)
272 for (
size_t i = 0; i < numRows; i++)
275 for (
size_t j = 0; j < numVecs; j++)
276 for (
size_t i = 0; i < numRows; i++)
277 Y(i, j) *=
ISC(beta);
280 Teuchos::LAPACK<int, LSC> lapack;
285 std::vector<LSC> yTemp(numVecs * numRows);
286 for (
size_t j = 0; j < numVecs; j++) {
287 for (
size_t i = 0; i < numRows; i++)
288 yTemp[j * numRows + i] = X(i, j);
293 (mode == Teuchos::CONJ_TRANS ?
'C' : (mode == Teuchos::TRANS ?
'T' :
'N'));
294 int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->scalarsPerRow_];
296 diagBlocks_[blockIndex].numRowsCols(),
298 diagBlocks_[blockIndex].DL(),
299 diagBlocks_[blockIndex].D(),
300 diagBlocks_[blockIndex].DU(),
301 diagBlocks_[blockIndex].DU2(),
308 for (
size_t j = 0; j < numVecs; j++) {
309 for (
size_t i = 0; i < numRows; i++) {
310 Y(i, j) *=
ISC(beta);
311 Y(i, j) +=
ISC(alpha * yTemp[j * numRows + i]);
315 for (
size_t j = 0; j < numVecs; j++) {
316 for (
size_t i = 0; i < numRows; i++)
317 Y(i, j) =
ISC(yTemp[j * numRows + i]);
321 TEUCHOS_TEST_FOR_EXCEPTION(
322 INFO != 0, std::runtime_error,
323 "Ifpack2::TriDiContainer::solveBlock: "
324 "LAPACK's _GETRS (solve using LU factorization with partial pivoting) "
325 "failed with INFO = "
326 << INFO <<
" != 0.");
330template <
class MatrixType,
class LocalScalarType>
332 Teuchos::FancyOStream fos(Teuchos::rcp(&os,
false));
333 fos.setOutputToRootOnly(0);
338template <
class MatrixType,
class LocalScalarType>
340 std::ostringstream oss;
341 oss << Teuchos::Describable::description();
342 if (this->isInitialized()) {
343 if (this->isComputed()) {
344 oss <<
"{status = initialized, computed";
346 oss <<
"{status = initialized, not computed";
349 oss <<
"{status = not initialized, not computed";
356template <
class MatrixType,
class LocalScalarType>
359 const Teuchos::EVerbosityLevel verbLevel)
const {
361 if (verbLevel == Teuchos::VERB_NONE)
return;
362 os <<
"================================================================================" << endl;
363 os <<
"Ifpack2::TriDiContainer" << endl;
364 os <<
"Number of blocks = " << this->numBlocks_ << endl;
365 os <<
"isInitialized() = " << this->IsInitialized_ << endl;
366 os <<
"isComputed() = " << this->IsComputed_ << endl;
367 os <<
"================================================================================" << endl;
371template <
class MatrixType,
class LocalScalarType>
376#define IFPACK2_TRIDICONTAINER_INSTANT(S, LO, GO, N) \
377 template class Ifpack2::TriDiContainer<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
typename KokkosKernels::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View.
Definition Ifpack2_Container_decl.hpp:101
The implementation of the numerical features of Container (Jacobi, Gauss-Seidel, SGS)....
Definition Ifpack2_Container_decl.hpp:307
static bool debug()
Whether Ifpack2 is in debug mode.
Definition Ifpack2_Details_Behavior.cpp:31
Store and solve a local TriDi linear problem.
Definition Ifpack2_TriDiContainer_decl.hpp:73
virtual std::string description() const
A one-line description of this object.
Definition Ifpack2_TriDiContainer_def.hpp:339
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition Ifpack2_TriDiContainer_def.hpp:67
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition Ifpack2_TriDiContainer_def.hpp:77
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_TriDiContainer_def.hpp:358
virtual ~TriDiContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition Ifpack2_TriDiContainer_def.hpp:64
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition Ifpack2_TriDiContainer_def.hpp:331
static std::string getName()
Get the name of this container type for Details::constructContainer()
Definition Ifpack2_TriDiContainer_def.hpp:372
TriDiContainer(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_TriDiContainer_def.hpp:27
void solveBlock(ConstHostSubviewLocal X, HostSubviewLocal Y, int blockIndex, Teuchos::ETransp mode, LSC alpha, LSC beta) const
Definition Ifpack2_TriDiContainer_def.hpp:227
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40