43 #ifndef IFPACK2_EXPERIMENTAL_CRSRBILUK_DEF_HPP 44 #define IFPACK2_EXPERIMENTAL_CRSRBILUK_DEF_HPP 46 #include <Tpetra_Experimental_BlockMultiVector.hpp> 47 #include<Ifpack2_OverlappingRowMatrix.hpp> 48 #include<Ifpack2_LocalFilter.hpp> 49 #include <Ifpack2_Experimental_RBILUK.hpp> 51 #include <Ifpack2_RILUK.hpp> 54 #include <Kokkos_ArithTraits.hpp> 60 template<
class MatrixType>
64 A_block_(
Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(Matrix_in))
67 template<
class MatrixType>
74 template<
class MatrixType>
78 template<
class MatrixType>
86 if (A.getRawPtr () != A_block_.getRawPtr ())
88 this->isAllocated_ =
false;
89 this->isInitialized_ =
false;
90 this->isComputed_ =
false;
91 this->
Graph_ = Teuchos::null;
92 L_block_ = Teuchos::null;
93 U_block_ = Teuchos::null;
94 D_block_ = Teuchos::null;
101 template<
class MatrixType>
102 const typename RBILUK<MatrixType>::block_crs_matrix_type&
105 TEUCHOS_TEST_FOR_EXCEPTION(
106 L_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getL: The L factor " 107 "is null. Please call initialize() (and preferably also compute()) " 108 "before calling this method. If the input matrix has not yet been set, " 109 "you must first call setMatrix() with a nonnull input matrix before you " 110 "may call initialize() or compute().");
115 template<
class MatrixType>
116 const typename RBILUK<MatrixType>::block_crs_matrix_type&
119 TEUCHOS_TEST_FOR_EXCEPTION(
120 D_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getD: The D factor " 121 "(of diagonal entries) is null. Please call initialize() (and " 122 "preferably also compute()) before calling this method. If the input " 123 "matrix has not yet been set, you must first call setMatrix() with a " 124 "nonnull input matrix before you may call initialize() or compute().");
129 template<
class MatrixType>
130 const typename RBILUK<MatrixType>::block_crs_matrix_type&
133 TEUCHOS_TEST_FOR_EXCEPTION(
134 U_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getU: The U factor " 135 "is null. Please call initialize() (and preferably also compute()) " 136 "before calling this method. If the input matrix has not yet been set, " 137 "you must first call setMatrix() with a nonnull input matrix before you " 138 "may call initialize() or compute().");
142 template<
class MatrixType>
148 if (! this->isAllocated_) {
160 L_block_ = rcp(
new block_crs_matrix_type(*this->
Graph_->getL_Graph (), blockSize_) );
161 U_block_ = rcp(
new block_crs_matrix_type(*this->
Graph_->getU_Graph (), blockSize_) );
162 D_block_ = rcp(
new block_crs_matrix_type(*(Ifpack2::Details::computeDiagonalGraph(*(this->
Graph_->getOverlapGraph()))),
164 L_block_->setAllToScalar (STM::zero ());
165 U_block_->setAllToScalar (STM::zero ());
166 D_block_->setAllToScalar (STM::zero ());
169 this->isAllocated_ =
true;
172 template<
class MatrixType>
173 Teuchos::RCP<const typename RBILUK<MatrixType>::block_crs_matrix_type>
178 template<
class MatrixType>
183 using Teuchos::rcp_dynamic_cast;
185 if (A_block_.is_null()) {
186 TEUCHOS_TEST_FOR_EXCEPTION(A_.is_null(), std::runtime_error,
"Ifpack2::Experimental::RBILUK::initialize: " 187 "The matrix is null. Please call setMatrix() with a nonnull input before calling this method.");
189 TEUCHOS_TEST_FOR_EXCEPTION(filteredA.is_null(), std::runtime_error,
"Ifpack2::Experimental::RBILUK::initialize: " 190 "Cannot cast to filtered matrix.");
192 if (overlappedA != Teuchos::null) {
193 A_block_ = rcp_dynamic_cast<
const block_crs_matrix_type>(overlappedA->getUnderlyingMatrix());
196 A_block_ = rcp_dynamic_cast<
const block_crs_matrix_type>(filteredA->getUnderlyingMatrix());
200 TEUCHOS_TEST_FOR_EXCEPTION(
201 A_block_.is_null (), std::runtime_error,
"Ifpack2::Experimental::RBILUK::initialize: " 202 "The matrix is null. Please call setMatrix() with a nonnull input " 203 "before calling this method.");
205 blockSize_ = A_block_->getBlockSize();
207 Teuchos::Time timer (
"RBILUK::initialize");
209 Teuchos::TimeMonitor timeMon (timer);
218 this->isInitialized_ =
false;
219 this->isAllocated_ =
false;
220 this->isComputed_ =
false;
221 this->
Graph_ = Teuchos::null;
227 RCP<const crs_graph_type> matrixCrsGraph = Teuchos::rcpFromRef(A_block_->getCrsGraph() );
229 this->LevelOfFill_, 0));
231 this->Graph_->initialize ();
232 allocate_L_and_U_blocks ();
233 initAllValues (*A_block_);
236 this->isInitialized_ =
true;
237 this->numInitialize_ += 1;
238 this->initializeTime_ += timer.totalElapsedTime ();
242 template<
class MatrixType>
249 typedef Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> map_type;
252 bool DiagFound =
false;
253 size_t NumNonzeroDiags = 0;
254 size_t MaxNumEntries = A.getNodeMaxNumRowEntries();
260 Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getNodeElementList();
261 Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getNodeElementList();
262 bool gidsAreConsistentlyOrdered=
true;
265 if (rowGIDs[i] != colGIDs[i]) {
266 gidsAreConsistentlyOrdered=
false;
267 indexOfInconsistentGID=i;
271 TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered==
false, std::runtime_error,
272 "The ordering of the local GIDs in the row and column maps is not the same" 273 << std::endl <<
"at index " << indexOfInconsistentGID
274 <<
". Consistency is required, as all calculations are done with" 275 << std::endl <<
"local indexing.");
279 Teuchos::Array<local_ordinal_type> LI(MaxNumEntries);
280 Teuchos::Array<local_ordinal_type> UI(MaxNumEntries);
281 Teuchos::Array<scalar_type> LV(MaxNumEntries*blockMatSize);
282 Teuchos::Array<scalar_type> UV(MaxNumEntries*blockMatSize);
284 Teuchos::Array<scalar_type> diagValues(blockMatSize);
286 L_block_->setAllToScalar (STM::zero ());
287 U_block_->setAllToScalar (STM::zero ());
288 D_block_->setAllToScalar (STM::zero ());
290 RCP<const map_type> rowMap = L_block_->getRowMap ();
300 for (
size_t myRow=0; myRow<A.getNodeNumRows(); ++myRow) {
307 A.getLocalRowView(local_row, InI, InV, NumIn);
319 if (k == local_row) {
323 diagValues[jj] = this->Rthresh_ * InV[blockOffset+jj] + IFPACK2_SGN(InV[blockOffset+jj]) * this->Athresh_;
324 D_block_->replaceLocalValues(local_row, &InI[j], diagValues.getRawPtr(), 1);
327 TEUCHOS_TEST_FOR_EXCEPTION(
328 true, std::runtime_error,
"Ifpack2::RILUK::initAllValues: current " 329 "GID k = " << k <<
" < 0. I'm not sure why this is an error; it is " 330 "probably an artifact of the undocumented assumptions of the " 331 "original implementation (likely copied and pasted from Ifpack). " 332 "Nevertheless, the code I found here insisted on this being an error " 333 "state, so I will throw an exception here.");
335 else if (k < local_row) {
339 LV[LBlockOffset+jj] = InV[blockOffset+jj];
342 else if (Teuchos::as<size_t>(k) <= rowMap->getNodeNumElements()) {
346 UV[UBlockOffset+jj] = InV[blockOffset+jj];
358 diagValues[jj*(blockSize_+1)] = this->Athresh_;
359 D_block_->replaceLocalValues(local_row, &local_row, diagValues.getRawPtr(), 1);
363 L_block_->replaceLocalValues(local_row, &LI[0], &LV[0], NumL);
367 U_block_->replaceLocalValues(local_row, &UI[0], &UV[0], NumU);
371 this->isInitialized_ =
true;
374 template<
class MatrixType>
380 TEUCHOS_TEST_FOR_EXCEPTION(
381 A_block_.is_null (), std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 382 "The matrix is null. Please call setMatrix() with a nonnull input " 383 "before calling this method.");
389 typedef typename Tpetra::Details::GetLapackType<impl_scalar_type>::lapack_scalar_type LST;
390 typedef typename Tpetra::Details::GetLapackType<impl_scalar_type>::lapack_type lapack_type;
394 Teuchos::Time timer (
"RBILUK::compute");
396 this->isComputed_ =
false;
403 initAllValues (*A_block_);
409 const size_t MaxNumEntries =
410 L_block_->getNodeMaxNumRowEntries () + U_block_->getNodeMaxNumRowEntries () + 1;
417 Teuchos::Array<int> ipiv(blockSize_);
418 Teuchos::Array<LST> work(1);
420 size_t num_cols = U_block_->getColMap()->getNodeNumElements();
421 Teuchos::Array<int> colflag(num_cols);
423 Teuchos::Array<impl_scalar_type> diagMod(blockMatSize,STM::zero());
424 little_block_type diagModBlock(&diagMod[0], blockSize_, rowStride, colStride);
425 Teuchos::Array<impl_scalar_type> matTmpArray(blockMatSize,STM::zero());
426 little_block_type matTmp(&matTmpArray[0], blockSize_, rowStride, colStride);
427 Teuchos::Array<impl_scalar_type> multiplierArray(blockMatSize, STM::zero());
428 little_block_type multiplier(&multiplierArray[0], blockSize_, rowStride, colStride);
436 for (
size_t j = 0; j < num_cols; ++j) {
439 Teuchos::Array<local_ordinal_type> InI(MaxNumEntries, 0);
440 Teuchos::Array<scalar_type> InV(MaxNumEntries*blockMatSize,STM::zero());
447 NumIn = MaxNumEntries;
451 L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
455 little_block_type lmat(&valsL[matOffset],blockSize_,rowStride, colStride);
456 little_block_type lmatV(&InV[matOffset],blockSize_,rowStride, colStride);
458 InI[j] = colValsL[j];
461 little_block_type dmat = D_block_->getLocalBlock(local_row, local_row);
462 little_block_type dmatV(&InV[NumL*blockMatSize], blockSize_, rowStride, colStride);
464 InI[NumL] = local_row;
468 U_block_->getLocalRowView(local_row, colValsU, valsU, NumURead);
472 if (!(colValsU[j] < numLocalRows))
continue;
473 InI[NumL+1+j] = colValsU[j];
475 little_block_type umat(&valsU[blockMatSize*j], blockSize_, rowStride, colStride);
476 little_block_type umatV(&InV[matOffset], blockSize_, rowStride, colStride);
483 for (
size_t j = 0; j < NumIn; ++j) {
489 diagModBlock.fill(diagmod);
493 little_block_type currentVal(&InV[jj*blockMatSize], blockSize_, rowStride, colStride);
494 multiplier.assign(currentVal);
496 const little_block_type dmatInverse = D_block_->getLocalBlock(j,j);
497 blockMatOpts.square_matrix_matrix_multiply(reinterpret_cast<impl_scalar_type*>(currentVal.getRawPtr()), reinterpret_cast<impl_scalar_type*>(dmatInverse.getRawPtr()), reinterpret_cast<impl_scalar_type*>(matTmp.getRawPtr()), blockSize_);
498 currentVal.assign(matTmp);
502 U_block_->getLocalRowView(j, UUI, UUV, NumUU);
504 if (this->RelaxValue_ == STM::zero ()) {
506 if (!(UUI[k] < numLocalRows))
continue;
507 const int kk = colflag[UUI[k]];
509 little_block_type kkval(&InV[kk*blockMatSize], blockSize_, rowStride, colStride);
510 little_block_type uumat(&UUV[k*blockMatSize], blockSize_, rowStride, colStride);
511 blockMatOpts.square_matrix_matrix_multiply(reinterpret_cast<impl_scalar_type*>(multiplier.getRawPtr()), reinterpret_cast<impl_scalar_type*>(uumat.getRawPtr()), reinterpret_cast<impl_scalar_type*>(kkval.getRawPtr()), blockSize_, -STM::one(), STM::one());
517 if (!(UUI[k] < numLocalRows))
continue;
518 const int kk = colflag[UUI[k]];
519 little_block_type uumat(&UUV[k*blockMatSize], blockSize_, rowStride, colStride);
521 little_block_type kkval(&InV[kk*blockMatSize], blockSize_, rowStride, colStride);
522 blockMatOpts.square_matrix_matrix_multiply(reinterpret_cast<impl_scalar_type*>(multiplier.getRawPtr()), reinterpret_cast<impl_scalar_type*>(uumat.getRawPtr()), reinterpret_cast<impl_scalar_type*>(kkval.getRawPtr()), blockSize_, -STM::one(), STM::one());
525 blockMatOpts.square_matrix_matrix_multiply(reinterpret_cast<impl_scalar_type*>(multiplier.getRawPtr()), reinterpret_cast<impl_scalar_type*>(uumat.getRawPtr()), reinterpret_cast<impl_scalar_type*>(diagModBlock.getRawPtr()), blockSize_, -STM::one(), STM::one());
532 L_block_->replaceLocalValues(local_row, InI.getRawPtr(), InV.getRawPtr(), NumL);
537 if (this->RelaxValue_ != STM::zero ()) {
538 dmat.update(this->RelaxValue_, diagModBlock);
552 LST*
const d_raw =
reinterpret_cast<LST*
> (dmat.getRawPtr());
554 for (
int k = 0; k < blockSize_; ++k) {
558 lapack.GETRF(blockSize_, blockSize_, d_raw, blockSize_, ipiv.getRawPtr(), &lapackInfo);
559 TEUCHOS_TEST_FOR_EXCEPTION(
560 lapackInfo != 0, std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 561 "lapackInfo = " << lapackInfo <<
" which indicates an error in the factorization GETRF.");
564 lapack.GETRI(blockSize_, d_raw, blockSize_, ipiv.getRawPtr(), work.getRawPtr(), lwork, &lapackInfo);
565 TEUCHOS_TEST_FOR_EXCEPTION(
566 lapackInfo != 0, std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 567 "lapackInfo = " << lapackInfo <<
" which indicates an error in the matrix inverse GETRI.");
569 typedef typename Kokkos::Details::ArithTraits<impl_scalar_type>::mag_type ImplMagnitudeType;
570 ImplMagnitudeType worksize = Kokkos::Details::ArithTraits<impl_scalar_type>::magnitude(work[0]);
571 lwork =
static_cast<int>(worksize);
573 lapack.GETRI(blockSize_, d_raw, blockSize_, ipiv.getRawPtr(), work.getRawPtr(), lwork, &lapackInfo);
574 TEUCHOS_TEST_FOR_EXCEPTION(
575 lapackInfo != 0, std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 576 "lapackInfo = " << lapackInfo <<
" which indicates an error in the matrix inverse GETRI.");
580 little_block_type currentVal(&InV[(NumL+1+j)*blockMatSize], blockSize_, rowStride, colStride);
582 blockMatOpts.square_matrix_matrix_multiply(reinterpret_cast<impl_scalar_type*>(dmat.getRawPtr()), reinterpret_cast<impl_scalar_type*>(currentVal.getRawPtr()), reinterpret_cast<impl_scalar_type*>(matTmp.getRawPtr()), blockSize_);
583 currentVal.assign(matTmp);
588 U_block_->replaceLocalValues (local_row, &InI[NumL+1], &InV[blockMatSize*(NumL+1)], NumU);
592 for (
size_t j = 0; j < num_cols; ++j) {
599 this->isComputed_ =
true;
600 this->numCompute_ += 1;
601 this->computeTime_ += timer.totalElapsedTime ();
605 template<
class MatrixType>
608 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
609 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
610 Teuchos::ETransp mode,
615 typedef Tpetra::Experimental::BlockMultiVector<
scalar_type,
620 TEUCHOS_TEST_FOR_EXCEPTION(
621 A_block_.is_null (), std::runtime_error,
"Ifpack2::Experimental::RBILUK::apply: The matrix is " 622 "null. Please call setMatrix() with a nonnull input, then initialize() " 623 "and compute(), before calling this method.");
624 TEUCHOS_TEST_FOR_EXCEPTION(
626 "Ifpack2::Experimental::RBILUK::apply: If you have not yet called compute(), " 627 "you must call compute() before calling this method.");
628 TEUCHOS_TEST_FOR_EXCEPTION(
629 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
630 "Ifpack2::Experimental::RBILUK::apply: X and Y do not have the same number of columns. " 631 "X.getNumVectors() = " << X.getNumVectors ()
632 <<
" != Y.getNumVectors() = " << Y.getNumVectors () <<
".");
633 TEUCHOS_TEST_FOR_EXCEPTION(
634 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
635 "Ifpack2::Experimental::RBILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for " 636 "complex Scalar type. Please talk to the Ifpack2 developers to get this " 637 "fixed. There is a FIXME in this file about this very issue.");
639 const local_ordinal_type blockMatSize = blockSize_*blockSize_;
641 const local_ordinal_type rowStride = blockSize_;
642 const local_ordinal_type colStride = 1;
644 BMV yBlock (Y, * (A_block_->getGraph ()->getDomainMap ()), blockSize_);
645 const BMV xBlock (X, * (A_block_->getColMap ()), blockSize_);
647 Teuchos::Array<scalar_type> lclarray(blockSize_);
648 little_vec_type lclvec(&lclarray[0], blockSize_, colStride);
649 const scalar_type one = STM::one ();
650 const scalar_type zero = STM::zero ();
652 Teuchos::Time timer (
"RBILUK::apply");
654 Teuchos::TimeMonitor timeMon (timer);
655 if (alpha == one && beta == zero) {
656 if (mode == Teuchos::NO_TRANS) {
663 const local_ordinal_type numVectors = xBlock.getNumVectors();
664 BMV cBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
665 BMV rBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
666 for (local_ordinal_type imv = 0; imv < numVectors; ++imv)
668 for (
size_t i = 0; i < D_block_->getNodeNumRows(); ++i)
670 local_ordinal_type local_row = i;
671 little_vec_type xval = xBlock.getLocalBlock(local_row,imv);
672 little_vec_type cval = cBlock.getLocalBlock(local_row,imv);
675 local_ordinal_type NumL;
676 const local_ordinal_type * colValsL;
679 L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
681 for (local_ordinal_type j = 0; j < NumL; ++j)
683 local_ordinal_type col = colValsL[j];
684 little_vec_type prevVal = cBlock.getLocalBlock(col, imv);
686 const local_ordinal_type matOffset = blockMatSize*j;
687 little_block_type lij(&valsL[matOffset],blockSize_,rowStride, colStride);
689 cval.matvecUpdate(-one, lij, prevVal);
695 D_block_->applyBlock(cBlock, rBlock);
698 for (local_ordinal_type imv = 0; imv < numVectors; ++imv)
700 const local_ordinal_type numRows = D_block_->getNodeNumRows();
701 for (local_ordinal_type i = 0; i < numRows; ++i)
703 local_ordinal_type local_row = (numRows-1)-i;
704 little_vec_type rval = rBlock.getLocalBlock(local_row,imv);
705 little_vec_type yval = yBlock.getLocalBlock(local_row,imv);
708 local_ordinal_type NumU;
709 const local_ordinal_type * colValsU;
712 U_block_->getLocalRowView(local_row, colValsU, valsU, NumU);
714 for (local_ordinal_type j = 0; j < NumU; ++j)
716 local_ordinal_type col = colValsU[NumU-1-j];
717 little_vec_type prevVal = yBlock.getLocalBlock(col, imv);
719 const local_ordinal_type matOffset = blockMatSize*(NumU-1-j);
720 little_block_type uij(&valsU[matOffset], blockSize_, rowStride, colStride);
722 yval.matvecUpdate(-one, uij, prevVal);
728 TEUCHOS_TEST_FOR_EXCEPTION(
729 true, std::runtime_error,
730 "Ifpack2::Experimental::RBILUK::apply: transpose apply is not implemented for the block algorithm. ");
741 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
742 apply (X, Y_tmp, mode);
743 Y.update (alpha, Y_tmp, beta);
748 this->numApply_ += 1;
749 this->applyTime_ = timer.totalElapsedTime ();
753 template<
class MatrixType>
756 std::ostringstream os;
761 os <<
"\"Ifpack2::Experimental::RBILUK\": {";
762 os <<
"Initialized: " << (this->
isInitialized () ?
"true" :
"false") <<
", " 763 <<
"Computed: " << (this->
isComputed () ?
"true" :
"false") <<
", ";
767 if (A_block_.is_null ()) {
768 os <<
"Matrix: null";
771 os <<
"Global matrix dimensions: [" 772 << A_block_->getGlobalNumRows () <<
", " << A_block_->getGlobalNumCols () <<
"]" 773 <<
", Global nnz: " << A_block_->getGlobalNumEntries();
788 #define IFPACK2_EXPERIMENTAL_RBILUK_INSTANT(S,LO,GO,N) \ 789 template class Ifpack2::Experimental::RBILUK< Tpetra::Experimental::BlockCrsMatrix<S, LO, GO, N> >; \ 790 template class Ifpack2::Experimental::RBILUK< Tpetra::RowMatrix<S, LO, GO, N> >; Ifpack2 features that are experimental. Use at your own risk.
void initialize()
Initialize by computing the symbolic incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:179
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:179
virtual ~RBILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_Experimental_RBILUK_def.hpp:75
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const
Get the input matrix.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:174
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:131
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:170
Teuchos::RCP< crs_matrix_type > L_
The L (lower triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:571
void apply(const Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &X, Tpetra::MultiVector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, scalar_type alpha=Teuchos::ScalarTraits< scalar_type >::one(), scalar_type beta=Teuchos::ScalarTraits< scalar_type >::zero()) const
Apply the (inverse of the) incomplete factorization to X, resulting in Y.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:608
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:243
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:182
bool isComputed() const
Whether compute() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:344
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:176
const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:117
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:375
int getLevelOfFill() const
Get level of fill (the "k" in ILU(k)).
Definition: Ifpack2_RILUK_decl.hpp:501
Teuchos::RCP< Ifpack2::IlukGraph< Tpetra::CrsGraph< local_ordinal_type, global_ordinal_type, node_type > > > Graph_
The ILU(k) graph.
Definition: Ifpack2_RILUK_decl.hpp:559
File for utility functions.
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition: Ifpack2_IlukGraph.hpp:97
bool isInitialized() const
Whether initialize() has been called on this object.
Definition: Ifpack2_RILUK_decl.hpp:340
void setMatrix(const Teuchos::RCP< const block_crs_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:80
Definition: Ifpack2_Details_Amesos2Wrapper_decl.hpp:63
Teuchos::RCP< vec_type > D_
The diagonal entries of the ILU(k) factorization.
Definition: Ifpack2_RILUK_decl.hpp:575
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:103
Sparse matrix (Tpetra::RowMatrix subclass) with ghost rows.
Definition: Ifpack2_OverlappingRowMatrix_decl.hpp:59
Access only local rows and columns of a sparse matrix.
Definition: Ifpack2_LocalFilter_decl.hpp:160
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:191
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:754
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
ILU(k) factorization of a given Tpetra::Experimental::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:162
Teuchos::RCP< crs_matrix_type > U_
The U (upper triangular) factor of ILU(k).
Definition: Ifpack2_RILUK_decl.hpp:573