10#ifndef IFPACK2_TRIDICONTAINER_DEF_HPP
11#define IFPACK2_TRIDICONTAINER_DEF_HPP
13#include "Teuchos_LAPACK.hpp"
17# include "Teuchos_DefaultMpiComm.hpp"
19# 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),
35 TEUCHOS_TEST_FOR_EXCEPTION(
36 ! matrix->hasColMap (), std::invalid_argument,
"Ifpack2::TriDiContainer: "
37 "The constructor's input matrix must have a column Map.");
44 scalarOffsets_[i] = scalarTotal;
47 if(actualBlockRows == 1)
55 scalarTotal += 4 * (actualBlockRows - 1);
59 scalars_.resize(scalarTotal);
60 diagBlocks_.reserve(this->numBlocks_);
61 for(
int i = 0; i < this->numBlocks_; i++)
63 diagBlocks_.emplace_back(Teuchos::View, scalars_.data() + scalarOffsets_[i], this->blockSizes_[i] * this->scalarsPerRow_);
67template<
class MatrixType,
class LocalScalarType>
71template<
class MatrixType,
class LocalScalarType>
75 diagBlocks_[i].putScalar(Teuchos::ScalarTraits<LSC>::zero());
82template<
class MatrixType,
class LocalScalarType>
85 #ifdef HAVE_IFPACK2_DEBUG
86 TEUCHOS_TEST_FOR_EXCEPTION(
87 ipiv_.size () != this->blockRows_.size() * this->scalarsPerRow_, std::logic_error,
88 "Ifpack2::TriDiContainer::compute: ipiv_ array has the wrong size. "
89 "Please report this bug to the Ifpack2 developers.");
101template<
class MatrixType,
class LocalScalarType>
104 using Teuchos::Array;
105 using Teuchos::ArrayView;
106 SC zero = Teuchos::ScalarTraits<SC>::zero();
107 const LO INVALID = Teuchos::OrdinalTraits<LO>::invalid();
114 if(this->scalarsPerRow_ > 1)
116 Array<LO> colToBlockOffset(this->inputBlockMatrix_->getLocalNumCols(), INVALID);
117 for(
int i = 0; i < this->numBlocks_; i++)
120 LO blockStart = this->blockOffsets_[i];
121 LO blockEnd = blockStart + this->blockSizes_[i];
122 ArrayView<const LO> blockRows = this->getBlockRows(i);
126 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
128 LO localCol = this->translateRowToCol(blockRows[j]);
129 colToBlockOffset[localCol] = blockStart + j;
131 using h_inds_type =
typename block_crs_matrix_type::local_inds_host_view_type;
132 using h_vals_type =
typename block_crs_matrix_type::values_host_view_type;
133 for(LO blockRow = 0; blockRow < LO(blockRows.size()); blockRow++)
138 LO inputRow = this->blockRows_[blockStart + blockRow];
139 this->inputBlockMatrix_->getLocalRowView(inputRow, indices, values);
140 LO numEntries = (LO) indices.size();
141 for(LO k = 0; k < numEntries; k++)
143 LO colOffset = colToBlockOffset[indices[k]];
144 if(blockStart <= colOffset && colOffset < blockEnd)
149 LO blockCol = colOffset - blockStart;
150 for(LO bc = 0; bc < this->bcrsBlockSize_; bc++)
152 for(LO br = 0; br < this->bcrsBlockSize_; br++)
154 LO r = this->bcrsBlockSize_ * blockRow + br;
155 LO c = this->bcrsBlockSize_ * blockCol + bc;
156 auto val = values[k * (this->bcrsBlockSize_ * this->bcrsBlockSize_) + (br + this->bcrsBlockSize_ * bc)];
158 diagBlocks_[i](r, c) = val;
170 Array<LO> colToBlockOffset(this->inputMatrix_->getLocalNumCols() * this->bcrsBlockSize_, INVALID);
171 for(
int i = 0; i < this->numBlocks_; i++)
174 LO blockStart = this->blockOffsets_[i];
175 LO blockEnd = blockStart + this->blockSizes_[i];
176 ArrayView<const LO> blockRows = this->getBlockRows(i);
180 for(
size_t j = 0; j < size_t(blockRows.size()); j++)
183 LO localCol = this->translateRowToCol(blockRows[j]);
184 colToBlockOffset[localCol] = blockStart + j;
186 for(
size_t blockRow = 0; blockRow < size_t(blockRows.size()); blockRow++)
189 LO inputPointRow = this->blockRows_[blockStart + blockRow];
190 auto rowView = this->getInputRowView(inputPointRow);
191 for(
size_t k = 0; k < rowView.size(); k++)
193 LO colOffset = colToBlockOffset[rowView.ind(k)];
194 if(blockStart <= colOffset && colOffset < blockEnd)
196 LO blockCol = colOffset - blockStart;
197 auto val = rowView.val(k);
199 diagBlocks_[i](blockRow, blockCol) = rowView.val(k);
207template<
class MatrixType,
class LocalScalarType>
208void TriDiContainer<MatrixType, LocalScalarType>::clearBlocks ()
212 ContainerImpl<MatrixType, LocalScalarType>::clearBlocks();
215template<
class MatrixType,
class LocalScalarType>
216void TriDiContainer<MatrixType, LocalScalarType>::factor ()
218 for(
int i = 0; i < this->numBlocks_; i++)
220 Teuchos::LAPACK<int, LSC> lapack;
222 int* blockIpiv = &ipiv_[this->blockOffsets_[i] * this->scalarsPerRow_];
223 lapack.GTTRF (diagBlocks_[i].numRowsCols (),
227 diagBlocks_[i].DU2(),
230 TEUCHOS_TEST_FOR_EXCEPTION(
231 INFO < 0, std::logic_error,
"Ifpack2::TriDiContainer::factor: "
232 "LAPACK's _GTTRF (LU factorization with partial pivoting) was called "
233 "incorrectly. INFO = " << INFO <<
" < 0. "
234 "Please report this bug to the Ifpack2 developers.");
238 TEUCHOS_TEST_FOR_EXCEPTION(
239 INFO > 0, std::runtime_error,
"Ifpack2::TriDiContainer::factor: "
240 "LAPACK's _GTTRF (LU factorization with partial pivoting) reports that the "
241 "computed U factor is exactly singular. U(" << INFO <<
"," << INFO <<
") "
242 "(one-based index i) is exactly zero. This probably means that the input "
243 "matrix has a singular diagonal block.");
247template<
class MatrixType,
class LocalScalarType>
252 Teuchos::ETransp mode,
256 LSC zero = Teuchos::ScalarTraits<LSC>::zero();
257 size_t numVecs = X.extent(1);
258 size_t numRows = X.extent(0);
260 #ifdef HAVE_IFPACK2_DEBUG
261 TEUCHOS_TEST_FOR_EXCEPTION(
262 X.extent (0) != Y.extent (0),
263 std::logic_error,
"Ifpack2::TriDiContainer::solveBlock: X and Y have "
264 "incompatible dimensions (" << X.extent (0) <<
" resp. "
265 << Y.extent (0) <<
"). Please report this bug to "
266 "the Ifpack2 developers.");
267 TEUCHOS_TEST_FOR_EXCEPTION(
268 X.extent (0) !=
static_cast<size_t> (diagBlocks_[blockIndex].numRowsCols()),
269 std::logic_error,
"Ifpack2::TriDiContainer::solveBlock: The input "
270 "multivector X has incompatible dimensions from those of the "
271 "inverse operator (" << X.extent (0) <<
" vs. "
272 << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRowsCols () : diagBlocks_[blockIndex].numRowsCols())
273 <<
"). Please report this bug to the Ifpack2 developers.");
274 TEUCHOS_TEST_FOR_EXCEPTION(
275 Y.extent (0) !=
static_cast<size_t> (diagBlocks_[blockIndex].numRowsCols()),
276 std::logic_error,
"Ifpack2::TriDiContainer::solveBlock: The output "
277 "multivector Y has incompatible dimensions from those of the "
278 "inverse operator (" << Y.extent (0) <<
" vs. "
279 << (mode == Teuchos::NO_TRANS ? diagBlocks_[blockIndex].numRowsCols() : diagBlocks_[blockIndex].numRowsCols ())
280 <<
"). Please report this bug to the Ifpack2 developers.");
288 for(
size_t j = 0; j < numVecs; j++)
289 for(
size_t i = 0; i < numRows; i++)
293 for(
size_t j = 0; j < numVecs; j++)
294 for(
size_t i = 0; i < numRows; i++)
295 Y(i, j) *=
ISC(beta);
299 Teuchos::LAPACK<int, LSC> lapack;
304 std::vector<LSC> yTemp(numVecs * numRows);
305 for(
size_t j = 0; j < numVecs; j++)
307 for(
size_t i = 0; i < numRows; i++)
308 yTemp[j * numRows + i] = X(i, j);
313 (mode == Teuchos::CONJ_TRANS ?
'C' : (mode == Teuchos::TRANS ?
'T' :
'N'));
314 int* blockIpiv = &ipiv_[this->blockOffsets_[blockIndex] * this->
scalarsPerRow_];
316 diagBlocks_[blockIndex].numRowsCols(),
318 diagBlocks_[blockIndex].DL(),
319 diagBlocks_[blockIndex].D(),
320 diagBlocks_[blockIndex].DU(),
321 diagBlocks_[blockIndex].DU2(),
328 for(
size_t j = 0; j < numVecs; j++)
330 for(
size_t i = 0; i < numRows; i++)
332 Y(i, j) *=
ISC(beta);
333 Y(i, j) +=
ISC(alpha * yTemp[j * numRows + i]);
338 for(
size_t j = 0; j < numVecs; j++)
340 for(
size_t i = 0; i < numRows; i++)
341 Y(i, j) =
ISC(yTemp[j * numRows + i]);
345 TEUCHOS_TEST_FOR_EXCEPTION(
346 INFO != 0, std::runtime_error,
"Ifpack2::TriDiContainer::solveBlock: "
347 "LAPACK's _GETRS (solve using LU factorization with partial pivoting) "
348 "failed with INFO = " << INFO <<
" != 0.");
352template<
class MatrixType,
class LocalScalarType>
355 Teuchos::FancyOStream fos(Teuchos::rcp(&os,
false));
356 fos.setOutputToRootOnly(0);
361template<
class MatrixType,
class LocalScalarType>
364 std::ostringstream oss;
365 oss << Teuchos::Describable::description();
368 oss <<
"{status = initialized, computed";
371 oss <<
"{status = initialized, not computed";
375 oss <<
"{status = not initialized, not computed";
382template<
class MatrixType,
class LocalScalarType>
386 const Teuchos::EVerbosityLevel verbLevel)
const
389 if(verbLevel==Teuchos::VERB_NONE)
return;
390 os <<
"================================================================================" << endl;
391 os <<
"Ifpack2::TriDiContainer" << endl;
392 os <<
"Number of blocks = " << this->numBlocks_ << endl;
394 os <<
"isComputed() = " << this->
IsComputed_ << endl;
395 os <<
"================================================================================" << endl;
399template<
class MatrixType,
class LocalScalarType>
405#define IFPACK2_TRIDICONTAINER_INSTANT(S,LO,GO,N) \
406 template class Ifpack2::TriDiContainer< Tpetra::RowMatrix<S, LO, GO, N>, S >;
virtual void extract()=0
Extract the submatrices identified by the local indices set by the constructor. BlockMatrix may be an...
int numBlocks_
The number of blocks (partitions) in the container.
Definition Ifpack2_Container_decl.hpp:259
typename Kokkos::ArithTraits< SC >::val_type ISC
Internal representation of Scalar in Kokkos::View.
Definition Ifpack2_Container_decl.hpp:102
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
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
bool IsComputed_
If true, the container has been successfully computed.
Definition Ifpack2_Container_decl.hpp:289
virtual std::string description() const
A one-line description of this object.
Definition Ifpack2_TriDiContainer_def.hpp:362
virtual void initialize()
Do all set-up operations that only require matrix structure.
Definition Ifpack2_TriDiContainer_def.hpp:72
virtual void compute()
Extract the local diagonal block and prepare the solver.
Definition Ifpack2_TriDiContainer_def.hpp:83
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:385
virtual ~TriDiContainer()
Destructor (declared virtual for memory safety of derived classes).
Definition Ifpack2_TriDiContainer_def.hpp:69
virtual std::ostream & print(std::ostream &os) const
Print information about this object to the given output stream.
Definition Ifpack2_TriDiContainer_def.hpp:353
static std::string getName()
Get the name of this container type for Details::constructContainer().
Definition Ifpack2_TriDiContainer_def.hpp:400
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:249
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:41