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> 57 template<
class MatrixType>
61 A_block_(
Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(Matrix_in))
64 template<
class MatrixType>
71 template<
class MatrixType>
75 template<
class MatrixType>
85 if (A.
getRawPtr () != A_block_.getRawPtr ())
87 this->isAllocated_ =
false;
88 this->isInitialized_ =
false;
89 this->isComputed_ =
false;
90 this->Graph_ = Teuchos::null;
91 L_block_ = Teuchos::null;
92 U_block_ = Teuchos::null;
93 D_block_ = Teuchos::null;
100 template<
class MatrixType>
101 const typename RBILUK<MatrixType>::block_crs_matrix_type&
105 L_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getL: The L factor " 106 "is null. Please call initialize() (and preferably also compute()) " 107 "before calling this method. If the input matrix has not yet been set, " 108 "you must first call setMatrix() with a nonnull input matrix before you " 109 "may call initialize() or compute().");
114 template<
class MatrixType>
115 const typename RBILUK<MatrixType>::block_crs_matrix_type&
119 D_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getD: The D factor " 120 "(of diagonal entries) is null. Please call initialize() (and " 121 "preferably also compute()) before calling this method. If the input " 122 "matrix has not yet been set, you must first call setMatrix() with a " 123 "nonnull input matrix before you may call initialize() or compute().");
128 template<
class MatrixType>
129 const typename RBILUK<MatrixType>::block_crs_matrix_type&
133 U_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getU: The U factor " 134 "is null. Please call initialize() (and preferably also compute()) " 135 "before calling this method. If the input matrix has not yet been set, " 136 "you must first call setMatrix() with a nonnull input matrix before you " 137 "may call initialize() or compute().");
141 template<
class MatrixType>
147 if (! this->isAllocated_) {
159 L_block_ =
rcp(
new block_crs_matrix_type(*this->Graph_->getL_Graph (), blockSize_) );
160 U_block_ =
rcp(
new block_crs_matrix_type(*this->Graph_->getU_Graph (), blockSize_) );
161 D_block_ =
rcp(
new block_crs_matrix_type(*(Ifpack2::Details::computeDiagonalGraph(*(this->Graph_->getOverlapGraph()))),
163 L_block_->setAllToScalar (STM::zero ());
164 U_block_->setAllToScalar (STM::zero ());
165 D_block_->setAllToScalar (STM::zero ());
168 this->isAllocated_ =
true;
171 template<
class MatrixType>
177 template<
class MatrixType>
182 using Teuchos::rcp_dynamic_cast;
183 const char prefix[] =
"Ifpack2::Experimental::RBILUK::initialize: ";
194 if (A_block_.is_null ()) {
199 RCP<const LocalFilter<row_matrix_type> > filteredA =
202 (filteredA.is_null (), std::runtime_error, prefix <<
203 "Cannot cast to filtered matrix.");
204 RCP<const OverlappingRowMatrix<row_matrix_type> > overlappedA =
206 if (! overlappedA.is_null ()) {
207 A_block_ = rcp_dynamic_cast<
const block_crs_matrix_type>(overlappedA->getUnderlyingMatrix());
210 A_block_ = rcp_dynamic_cast<
const block_crs_matrix_type>(filteredA->getUnderlyingMatrix());
215 (A_block_.is_null (), std::runtime_error, prefix <<
"The matrix (A_block_, " 216 "the BlockCrsMatrix) is null. Please call setMatrix() with a nonnull " 217 "input before calling this method.");
219 (! A_block_->isFillComplete (), std::runtime_error, prefix <<
"The matrix " 220 "(A_block_, the BlockCrsMatrix) is not fill complete. You may not invoke " 221 "initialize() or compute() with this matrix until the matrix is fill " 222 "complete. Note: BlockCrsMatrix is fill complete if and only if its " 223 "underlying graph is fill complete.");
225 blockSize_ = A_block_->getBlockSize();
238 this->isInitialized_ =
false;
239 this->isAllocated_ =
false;
240 this->isComputed_ =
false;
241 this->Graph_ = Teuchos::null;
247 RCP<const crs_graph_type> matrixCrsGraph = Teuchos::rcpFromRef(A_block_->getCrsGraph() );
249 this->LevelOfFill_, 0));
251 this->Graph_->initialize ();
252 allocate_L_and_U_blocks ();
253 initAllValues (*A_block_);
256 this->isInitialized_ =
true;
257 this->numInitialize_ += 1;
262 template<
class MatrixType>
268 typedef Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> map_type;
270 local_ordinal_type NumIn = 0, NumL = 0, NumU = 0;
271 bool DiagFound =
false;
272 size_t NumNonzeroDiags = 0;
273 size_t MaxNumEntries = A.getNodeMaxNumRowEntries();
274 local_ordinal_type blockMatSize = blockSize_*blockSize_;
281 bool gidsAreConsistentlyOrdered=
true;
282 global_ordinal_type indexOfInconsistentGID=0;
283 for (global_ordinal_type i=0; i<rowGIDs.
size(); ++i) {
284 if (rowGIDs[i] != colGIDs[i]) {
285 gidsAreConsistentlyOrdered=
false;
286 indexOfInconsistentGID=i;
291 "The ordering of the local GIDs in the row and column maps is not the same" 292 << std::endl <<
"at index " << indexOfInconsistentGID
293 <<
". Consistency is required, as all calculations are done with" 294 << std::endl <<
"local indexing.");
305 L_block_->setAllToScalar (STM::zero ());
306 U_block_->setAllToScalar (STM::zero ());
307 D_block_->setAllToScalar (STM::zero ());
312 const_cast<block_crs_matrix_type&
> (A).
template sync<Kokkos::HostSpace> ();
313 L_block_->template sync<Kokkos::HostSpace> ();
314 U_block_->template sync<Kokkos::HostSpace> ();
315 D_block_->template sync<Kokkos::HostSpace> ();
317 L_block_->template modify<Kokkos::HostSpace> ();
318 U_block_->template modify<Kokkos::HostSpace> ();
319 D_block_->template modify<Kokkos::HostSpace> ();
321 RCP<const map_type> rowMap = L_block_->getRowMap ();
331 for (
size_t myRow=0; myRow<A.getNodeNumRows(); ++myRow) {
332 local_ordinal_type local_row = myRow;
336 const local_ordinal_type * InI = 0;
337 scalar_type * InV = 0;
338 A.getLocalRowView(local_row, InI, InV, NumIn);
346 for (local_ordinal_type j = 0; j < NumIn; ++j) {
347 const local_ordinal_type k = InI[j];
348 const local_ordinal_type blockOffset = blockMatSize*j;
350 if (k == local_row) {
353 for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
354 diagValues[jj] = this->Rthresh_ * InV[blockOffset+jj] + IFPACK2_SGN(InV[blockOffset+jj]) * this->Athresh_;
355 D_block_->replaceLocalValues(local_row, &InI[j], diagValues.getRawPtr(), 1);
359 true, std::runtime_error,
"Ifpack2::RILUK::initAllValues: current " 360 "GID k = " << k <<
" < 0. I'm not sure why this is an error; it is " 361 "probably an artifact of the undocumented assumptions of the " 362 "original implementation (likely copied and pasted from Ifpack). " 363 "Nevertheless, the code I found here insisted on this being an error " 364 "state, so I will throw an exception here.");
366 else if (k < local_row) {
368 const local_ordinal_type LBlockOffset = NumL*blockMatSize;
369 for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
370 LV[LBlockOffset+jj] = InV[blockOffset+jj];
373 else if (Teuchos::as<size_t>(k) <= rowMap->getNodeNumElements()) {
375 const local_ordinal_type UBlockOffset = NumU*blockMatSize;
376 for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
377 UV[UBlockOffset+jj] = InV[blockOffset+jj];
388 for (local_ordinal_type jj = 0; jj < blockSize_; ++jj)
389 diagValues[jj*(blockSize_+1)] = this->Athresh_;
390 D_block_->replaceLocalValues(local_row, &local_row, diagValues.getRawPtr(), 1);
394 L_block_->replaceLocalValues(local_row, &LI[0], &LV[0], NumL);
398 U_block_->replaceLocalValues(local_row, &UI[0], &UV[0], NumU);
405 typedef typename block_crs_matrix_type::device_type device_type;
406 const_cast<block_crs_matrix_type&
> (A).
template sync<device_type> ();
407 L_block_->template sync<device_type> ();
408 U_block_->template sync<device_type> ();
409 D_block_->template sync<device_type> ();
411 this->isInitialized_ =
true;
420 template<
class LittleBlockType>
421 struct GetManagedView {
422 static_assert (Kokkos::Impl::is_view<LittleBlockType>::value,
423 "LittleBlockType must be a Kokkos::View.");
424 typedef Kokkos::View<
typename LittleBlockType::non_const_data_type,
425 typename LittleBlockType::array_layout,
426 typename LittleBlockType::device_type> managed_non_const_type;
427 static_assert (static_cast<int> (managed_non_const_type::rank) ==
428 static_cast<int> (LittleBlockType::rank),
429 "managed_non_const_type::rank != LittleBlock::rank. " 430 "Please report this bug to the Ifpack2 developers.");
435 template<
class MatrixType>
438 typedef impl_scalar_type IST;
439 const char prefix[] =
"Ifpack2::Experimental::RBILUK::compute: ";
445 (A_block_.is_null (), std::runtime_error, prefix <<
"The matrix (A_block_, " 446 "the BlockCrsMatrix) is null. Please call setMatrix() with a nonnull " 447 "input before calling this method.");
449 (! A_block_->isFillComplete (), std::runtime_error, prefix <<
"The matrix " 450 "(A_block_, the BlockCrsMatrix) is not fill complete. You may not invoke " 451 "initialize() or compute() with this matrix until the matrix is fill " 452 "complete. Note: BlockCrsMatrix is fill complete if and only if its " 453 "underlying graph is fill complete.");
455 if (! this->isInitialized ()) {
462 if (! A_block_.is_null ()) {
464 Teuchos::rcp_const_cast<block_crs_matrix_type> (A_block_);
465 A_nc->template sync<Kokkos::HostSpace> ();
467 L_block_->template sync<Kokkos::HostSpace> ();
468 U_block_->template sync<Kokkos::HostSpace> ();
469 D_block_->template sync<Kokkos::HostSpace> ();
471 L_block_->template modify<Kokkos::HostSpace> ();
472 U_block_->template modify<Kokkos::HostSpace> ();
473 D_block_->template modify<Kokkos::HostSpace> ();
478 this->isComputed_ =
false;
485 initAllValues (*A_block_);
491 const size_t MaxNumEntries =
492 L_block_->getNodeMaxNumRowEntries () + U_block_->getNodeMaxNumRowEntries () + 1;
502 Kokkos::View<
int*, Kokkos::HostSpace,
503 Kokkos::MemoryUnmanaged> ipiv (ipiv_teuchos.
getRawPtr (), blockSize_);
505 Kokkos::View<IST*, Kokkos::HostSpace,
506 Kokkos::MemoryUnmanaged> work (work_teuchos.getRawPtr (), blockSize_);
508 size_t num_cols = U_block_->getColMap()->getNodeNumElements();
511 typename GetManagedView<little_block_type>::managed_non_const_type diagModBlock (
"diagModBlock", blockSize_, blockSize_);
512 typename GetManagedView<little_block_type>::managed_non_const_type matTmp (
"matTmp", blockSize_, blockSize_);
513 typename GetManagedView<little_block_type>::managed_non_const_type multiplier (
"multiplier", blockSize_, blockSize_);
521 for (
size_t j = 0; j < num_cols; ++j) {
532 NumIn = MaxNumEntries;
536 L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
540 little_block_type lmat((
typename little_block_type::value_type*) &valsL[matOffset],blockSize_,rowStride);
541 little_block_type lmatV((
typename little_block_type::value_type*) &InV[matOffset],blockSize_,rowStride);
543 Tpetra::Experimental::COPY (lmat, lmatV);
544 InI[j] = colValsL[j];
547 little_block_type dmat = D_block_->getLocalBlock(local_row, local_row);
548 little_block_type dmatV((
typename little_block_type::value_type*) &InV[NumL*blockMatSize], blockSize_, rowStride);
550 Tpetra::Experimental::COPY (dmat, dmatV);
551 InI[NumL] = local_row;
555 U_block_->getLocalRowView(local_row, colValsU, valsU, NumURead);
559 if (!(colValsU[j] < numLocalRows))
continue;
560 InI[NumL+1+j] = colValsU[j];
562 little_block_type umat((
typename little_block_type::value_type*) &valsU[blockMatSize*j], blockSize_, rowStride);
563 little_block_type umatV((
typename little_block_type::value_type*) &InV[matOffset], blockSize_, rowStride);
565 Tpetra::Experimental::COPY (umat, umatV);
571 for (
size_t j = 0; j < NumIn; ++j) {
576 Kokkos::deep_copy (diagModBlock, diagmod);
580 little_block_type currentVal((
typename little_block_type::value_type*) &InV[jj*blockMatSize], blockSize_, rowStride);
582 Tpetra::Experimental::COPY (currentVal, multiplier);
584 const little_block_type dmatInverse = D_block_->getLocalBlock(j,j);
586 Tpetra::Experimental::GEMM (
"N",
"N", STS::one (), currentVal, dmatInverse,
587 STS::zero (), matTmp);
590 Tpetra::Experimental::COPY (matTmp, currentVal);
594 U_block_->getLocalRowView(j, UUI, UUV, NumUU);
596 if (this->RelaxValue_ == STM::zero ()) {
598 if (!(UUI[k] < numLocalRows))
continue;
599 const int kk = colflag[UUI[k]];
601 little_block_type kkval((
typename little_block_type::value_type*) &InV[kk*blockMatSize], blockSize_, rowStride);
602 little_block_type uumat((
typename little_block_type::value_type*) &UUV[k*blockMatSize], blockSize_, rowStride);
603 Tpetra::Experimental::GEMM (
"N",
"N",
magnitude_type(-STM::one ()), multiplier, uumat,
611 if (!(UUI[k] < numLocalRows))
continue;
612 const int kk = colflag[UUI[k]];
613 little_block_type uumat((
typename little_block_type::value_type*) &UUV[k*blockMatSize], blockSize_, rowStride);
615 little_block_type kkval((
typename little_block_type::value_type*) &InV[kk*blockMatSize], blockSize_, rowStride);
616 Tpetra::Experimental::GEMM (
"N",
"N",
magnitude_type(-STM::one ()), multiplier, uumat,
621 Tpetra::Experimental::GEMM (
"N",
"N",
magnitude_type(-STM::one ()), multiplier, uumat,
622 STM::one (), diagModBlock);
634 Tpetra::Experimental::COPY (dmatV, dmat);
636 if (this->RelaxValue_ != STM::zero ()) {
638 Tpetra::Experimental::AXPY (this->RelaxValue_, diagModBlock, dmat);
652 for (
int k = 0; k < blockSize_; ++k) {
656 Tpetra::Experimental::GETF2 (dmat, ipiv, lapackInfo);
659 lapackInfo != 0, std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 660 "lapackInfo = " << lapackInfo <<
" which indicates an error in the factorization GETRF.");
662 Tpetra::Experimental::GETRI (dmat, ipiv, work, lapackInfo);
665 lapackInfo != 0, std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 666 "lapackInfo = " << lapackInfo <<
" which indicates an error in the matrix inverse GETRI.");
670 little_block_type currentVal((
typename little_block_type::value_type*) &InV[(NumL+1+j)*blockMatSize], blockSize_, rowStride);
672 Tpetra::Experimental::GEMM (
"N",
"N", STS::one (), dmat, currentVal,
673 STS::zero (), matTmp);
676 Tpetra::Experimental::COPY (matTmp, currentVal);
681 U_block_->replaceLocalValues (local_row, &InI[NumL+1], &InV[blockMatSize*(NumL+1)], NumU);
685 for (
size_t j = 0; j < num_cols; ++j) {
694 typedef typename block_crs_matrix_type::device_type device_type;
695 if (! A_block_.is_null ()) {
697 Teuchos::rcp_const_cast<block_crs_matrix_type> (A_block_);
698 A_nc->template sync<device_type> ();
700 L_block_->template sync<device_type> ();
701 U_block_->template sync<device_type> ();
702 D_block_->template sync<device_type> ();
705 this->isComputed_ =
true;
706 this->numCompute_ += 1;
707 this->computeTime_ += timer.totalElapsedTime ();
711 template<
class MatrixType>
714 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
715 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
721 typedef Tpetra::Experimental::BlockMultiVector<
scalar_type,
727 A_block_.is_null (), std::runtime_error,
"Ifpack2::Experimental::RBILUK::apply: The matrix is " 728 "null. Please call setMatrix() with a nonnull input, then initialize() " 729 "and compute(), before calling this method.");
731 ! this->isComputed (), std::runtime_error,
732 "Ifpack2::Experimental::RBILUK::apply: If you have not yet called compute(), " 733 "you must call compute() before calling this method.");
735 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
736 "Ifpack2::Experimental::RBILUK::apply: X and Y do not have the same number of columns. " 737 "X.getNumVectors() = " << X.getNumVectors ()
738 <<
" != Y.getNumVectors() = " << Y.getNumVectors () <<
".");
740 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
741 "Ifpack2::Experimental::RBILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for " 742 "complex Scalar type. Please talk to the Ifpack2 developers to get this " 743 "fixed. There is a FIXME in this file about this very issue.");
749 BMV yBlock (Y, * (A_block_->getGraph ()->getDomainMap ()), blockSize_);
750 const BMV xBlock (X, * (A_block_->getColMap ()), blockSize_);
753 little_vec_type lclvec((
typename little_vec_type::value_type*)&lclarray[0], blockSize_);
760 if (alpha == one && beta ==
zero) {
761 if (mode == Teuchos::NO_TRANS) {
769 BMV cBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
770 BMV rBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
773 for (
size_t i = 0; i < D_block_->getNodeNumRows(); ++i)
776 little_vec_type xval = xBlock.getLocalBlock(local_row,imv);
777 little_vec_type cval = cBlock.getLocalBlock(local_row,imv);
779 Tpetra::Experimental::COPY (xval, cval);
785 L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
790 little_vec_type prevVal = cBlock.getLocalBlock(col, imv);
793 little_block_type lij((
typename little_block_type::value_type*) &valsL[matOffset],blockSize_,rowStride);
796 Tpetra::Experimental::GEMV (-one, lij, prevVal, cval);
802 D_block_->applyBlock(cBlock, rBlock);
811 little_vec_type rval = rBlock.getLocalBlock(local_row,imv);
812 little_vec_type yval = yBlock.getLocalBlock(local_row,imv);
814 Tpetra::Experimental::COPY (rval, yval);
820 U_block_->getLocalRowView(local_row, colValsU, valsU, NumU);
825 little_vec_type prevVal = yBlock.getLocalBlock(col, imv);
828 little_block_type uij((
typename little_block_type::value_type*) &valsU[matOffset], blockSize_, rowStride);
831 Tpetra::Experimental::GEMV (-one, uij, prevVal, yval);
838 true, std::runtime_error,
839 "Ifpack2::Experimental::RBILUK::apply: transpose apply is not implemented for the block algorithm. ");
850 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
851 apply (X, Y_tmp, mode);
852 Y.update (alpha, Y_tmp, beta);
857 this->numApply_ += 1;
862 template<
class MatrixType>
865 std::ostringstream os;
870 os <<
"\"Ifpack2::Experimental::RBILUK\": {";
871 os <<
"Initialized: " << (this->isInitialized () ?
"true" :
"false") <<
", " 872 <<
"Computed: " << (this->isComputed () ?
"true" :
"false") <<
", ";
874 os <<
"Level-of-fill: " << this->getLevelOfFill() <<
", ";
876 if (A_block_.is_null ()) {
877 os <<
"Matrix: null";
880 os <<
"Global matrix dimensions: [" 881 << A_block_->getGlobalNumRows () <<
", " << A_block_->getGlobalNumCols () <<
"]" 882 <<
", Global nnz: " << A_block_->getGlobalNumEntries();
897 #define IFPACK2_EXPERIMENTAL_RBILUK_INSTANT(S,LO,GO,N) \ 898 template class Ifpack2::Experimental::RBILUK< Tpetra::Experimental::BlockCrsMatrix<S, LO, GO, N> >; \ 899 template class Ifpack2::Experimental::RBILUK< Tpetra::RowMatrix<S, LO, GO, N> >; const block_crs_matrix_type & getDBlock() const
Return the diagonal entries of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:116
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:178
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const
Get the input matrix.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:173
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:145
virtual ~RBILUK()
Destructor (declared virtual for memory safety).
Definition: Ifpack2_Experimental_RBILUK_def.hpp:72
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:136
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:151
ILU(k) factorization of a given Tpetra::RowMatrix.
Definition: Ifpack2_RILUK_decl.hpp:254
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:148
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:714
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:142
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:436
LinearOp zero(const VectorSpace &vs)
File for utility functions.
const block_crs_matrix_type & getUBlock() const
Return the U factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:130
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition: Ifpack2_IlukGraph.hpp:97
double totalElapsedTime(bool readCurrentTime=false) const
void setMatrix(const Teuchos::RCP< const block_crs_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:77
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:863
Sparse matrix (Tpetra::RowMatrix subclass) with ghost rows.
Definition: Ifpack2_OverlappingRowMatrix_decl.hpp:59
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:102
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:157
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:128