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 #define IFPACK2_RBILUK_INITIAL_NOKK 56 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 57 #include <KokkosBatched_Gemm_Decl.hpp> 58 #include <KokkosBatched_Gemm_Serial_Impl.hpp> 59 #include <KokkosBatched_Util.hpp> 66 template<
class MatrixType>
70 A_block_(
Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(Matrix_in))
73 template<
class MatrixType>
80 template<
class MatrixType>
84 template<
class MatrixType>
94 if (A.getRawPtr () != A_block_.getRawPtr ())
96 this->isAllocated_ =
false;
97 this->isInitialized_ =
false;
98 this->isComputed_ =
false;
99 this->Graph_ = Teuchos::null;
100 L_block_ = Teuchos::null;
101 U_block_ = Teuchos::null;
102 D_block_ = Teuchos::null;
109 template<
class MatrixType>
110 const typename RBILUK<MatrixType>::block_crs_matrix_type&
113 TEUCHOS_TEST_FOR_EXCEPTION(
114 L_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getL: The L factor " 115 "is null. Please call initialize() (and preferably also compute()) " 116 "before calling this method. If the input matrix has not yet been set, " 117 "you must first call setMatrix() with a nonnull input matrix before you " 118 "may call initialize() or compute().");
123 template<
class MatrixType>
124 const typename RBILUK<MatrixType>::block_crs_matrix_type&
127 TEUCHOS_TEST_FOR_EXCEPTION(
128 D_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getD: The D factor " 129 "(of diagonal entries) is null. Please call initialize() (and " 130 "preferably also compute()) before calling this method. If the input " 131 "matrix has not yet been set, you must first call setMatrix() with a " 132 "nonnull input matrix before you may call initialize() or compute().");
137 template<
class MatrixType>
138 const typename RBILUK<MatrixType>::block_crs_matrix_type&
141 TEUCHOS_TEST_FOR_EXCEPTION(
142 U_block_.is_null (), std::runtime_error,
"Ifpack2::RILUK::getU: The U factor " 143 "is null. Please call initialize() (and preferably also compute()) " 144 "before calling this method. If the input matrix has not yet been set, " 145 "you must first call setMatrix() with a nonnull input matrix before you " 146 "may call initialize() or compute().");
150 template<
class MatrixType>
156 if (! this->isAllocated_) {
168 L_block_ = rcp(
new block_crs_matrix_type(*this->Graph_->getL_Graph (), blockSize_) );
169 U_block_ = rcp(
new block_crs_matrix_type(*this->Graph_->getU_Graph (), blockSize_) );
170 D_block_ = rcp(
new block_crs_matrix_type(*(Ifpack2::Details::computeDiagonalGraph(*(this->Graph_->getOverlapGraph()))),
172 L_block_->setAllToScalar (STM::zero ());
173 U_block_->setAllToScalar (STM::zero ());
174 D_block_->setAllToScalar (STM::zero ());
177 this->isAllocated_ =
true;
180 template<
class MatrixType>
181 Teuchos::RCP<const typename RBILUK<MatrixType>::block_crs_matrix_type>
186 template<
class MatrixType>
191 using Teuchos::rcp_dynamic_cast;
192 const char prefix[] =
"Ifpack2::Experimental::RBILUK::initialize: ";
203 if (A_block_.is_null ()) {
208 RCP<const LocalFilter<row_matrix_type> > filteredA =
210 TEUCHOS_TEST_FOR_EXCEPTION
211 (filteredA.is_null (), std::runtime_error, prefix <<
212 "Cannot cast to filtered matrix.");
213 RCP<const OverlappingRowMatrix<row_matrix_type> > overlappedA =
215 if (! overlappedA.is_null ()) {
216 A_block_ = rcp_dynamic_cast<
const block_crs_matrix_type>(overlappedA->getUnderlyingMatrix());
219 A_block_ = rcp_dynamic_cast<
const block_crs_matrix_type>(filteredA->getUnderlyingMatrix());
223 TEUCHOS_TEST_FOR_EXCEPTION
224 (A_block_.is_null (), std::runtime_error, prefix <<
"The matrix (A_block_, " 225 "the BlockCrsMatrix) is null. Please call setMatrix() with a nonnull " 226 "input before calling this method.");
227 TEUCHOS_TEST_FOR_EXCEPTION
228 (! A_block_->isFillComplete (), std::runtime_error, prefix <<
"The matrix " 229 "(A_block_, the BlockCrsMatrix) is not fill complete. You may not invoke " 230 "initialize() or compute() with this matrix until the matrix is fill " 231 "complete. Note: BlockCrsMatrix is fill complete if and only if its " 232 "underlying graph is fill complete.");
234 blockSize_ = A_block_->getBlockSize();
236 Teuchos::Time timer (
"RBILUK::initialize");
238 Teuchos::TimeMonitor timeMon (timer);
247 this->isInitialized_ =
false;
248 this->isAllocated_ =
false;
249 this->isComputed_ =
false;
250 this->Graph_ = Teuchos::null;
256 RCP<const crs_graph_type> matrixCrsGraph = Teuchos::rcpFromRef(A_block_->getCrsGraph() );
258 this->LevelOfFill_, 0));
260 this->Graph_->initialize ();
261 allocate_L_and_U_blocks ();
262 #ifdef IFPACK2_RBILUK_INITIAL 263 initAllValues (*A_block_);
267 this->isInitialized_ =
true;
268 this->numInitialize_ += 1;
269 this->initializeTime_ += timer.totalElapsedTime ();
273 template<
class MatrixType>
279 typedef Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> map_type;
281 local_ordinal_type NumIn = 0, NumL = 0, NumU = 0;
282 bool DiagFound =
false;
283 size_t NumNonzeroDiags = 0;
284 size_t MaxNumEntries = A.getNodeMaxNumRowEntries();
285 local_ordinal_type blockMatSize = blockSize_*blockSize_;
290 Teuchos::ArrayView<const global_ordinal_type> rowGIDs = A.getRowMap()->getNodeElementList();
291 Teuchos::ArrayView<const global_ordinal_type> colGIDs = A.getColMap()->getNodeElementList();
292 bool gidsAreConsistentlyOrdered=
true;
293 global_ordinal_type indexOfInconsistentGID=0;
294 for (global_ordinal_type i=0; i<rowGIDs.size(); ++i) {
295 if (rowGIDs[i] != colGIDs[i]) {
296 gidsAreConsistentlyOrdered=
false;
297 indexOfInconsistentGID=i;
301 TEUCHOS_TEST_FOR_EXCEPTION(gidsAreConsistentlyOrdered==
false, std::runtime_error,
302 "The ordering of the local GIDs in the row and column maps is not the same" 303 << std::endl <<
"at index " << indexOfInconsistentGID
304 <<
". Consistency is required, as all calculations are done with" 305 << std::endl <<
"local indexing.");
309 Teuchos::Array<local_ordinal_type> LI(MaxNumEntries);
310 Teuchos::Array<local_ordinal_type> UI(MaxNumEntries);
311 Teuchos::Array<scalar_type> LV(MaxNumEntries*blockMatSize);
312 Teuchos::Array<scalar_type> UV(MaxNumEntries*blockMatSize);
314 Teuchos::Array<scalar_type> diagValues(blockMatSize);
316 L_block_->setAllToScalar (STM::zero ());
317 U_block_->setAllToScalar (STM::zero ());
318 D_block_->setAllToScalar (STM::zero ());
323 const_cast<block_crs_matrix_type&
> (A).
template sync<Kokkos::HostSpace> ();
324 L_block_->template sync<Kokkos::HostSpace> ();
325 U_block_->template sync<Kokkos::HostSpace> ();
326 D_block_->template sync<Kokkos::HostSpace> ();
328 L_block_->template modify<Kokkos::HostSpace> ();
329 U_block_->template modify<Kokkos::HostSpace> ();
330 D_block_->template modify<Kokkos::HostSpace> ();
332 RCP<const map_type> rowMap = L_block_->getRowMap ();
342 for (
size_t myRow=0; myRow<A.getNodeNumRows(); ++myRow) {
343 local_ordinal_type local_row = myRow;
347 const local_ordinal_type * InI = 0;
348 scalar_type * InV = 0;
349 A.getLocalRowView(local_row, InI, InV, NumIn);
357 for (local_ordinal_type j = 0; j < NumIn; ++j) {
358 const local_ordinal_type k = InI[j];
359 const local_ordinal_type blockOffset = blockMatSize*j;
361 if (k == local_row) {
364 for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
365 diagValues[jj] = this->Rthresh_ * InV[blockOffset+jj] + IFPACK2_SGN(InV[blockOffset+jj]) * this->Athresh_;
366 D_block_->replaceLocalValues(local_row, &InI[j], diagValues.getRawPtr(), 1);
369 TEUCHOS_TEST_FOR_EXCEPTION(
370 true, std::runtime_error,
"Ifpack2::RILUK::initAllValues: current " 371 "GID k = " << k <<
" < 0. I'm not sure why this is an error; it is " 372 "probably an artifact of the undocumented assumptions of the " 373 "original implementation (likely copied and pasted from Ifpack). " 374 "Nevertheless, the code I found here insisted on this being an error " 375 "state, so I will throw an exception here.");
377 else if (k < local_row) {
379 const local_ordinal_type LBlockOffset = NumL*blockMatSize;
380 for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
381 LV[LBlockOffset+jj] = InV[blockOffset+jj];
384 else if (Teuchos::as<size_t>(k) <= rowMap->getNodeNumElements()) {
386 const local_ordinal_type UBlockOffset = NumU*blockMatSize;
387 for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
388 UV[UBlockOffset+jj] = InV[blockOffset+jj];
399 for (local_ordinal_type jj = 0; jj < blockSize_; ++jj)
400 diagValues[jj*(blockSize_+1)] = this->Athresh_;
401 D_block_->replaceLocalValues(local_row, &local_row, diagValues.getRawPtr(), 1);
405 L_block_->replaceLocalValues(local_row, &LI[0], &LV[0], NumL);
409 U_block_->replaceLocalValues(local_row, &UI[0], &UV[0], NumU);
416 typedef typename block_crs_matrix_type::device_type device_type;
417 const_cast<block_crs_matrix_type&
> (A).
template sync<device_type> ();
418 L_block_->template sync<device_type> ();
419 U_block_->template sync<device_type> ();
420 D_block_->template sync<device_type> ();
422 this->isInitialized_ =
true;
431 template<
class LittleBlockType>
432 struct GetManagedView {
433 static_assert (Kokkos::Impl::is_view<LittleBlockType>::value,
434 "LittleBlockType must be a Kokkos::View.");
435 typedef Kokkos::View<
typename LittleBlockType::non_const_data_type,
436 typename LittleBlockType::array_layout,
437 typename LittleBlockType::device_type> managed_non_const_type;
438 static_assert (static_cast<int> (managed_non_const_type::rank) ==
439 static_cast<int> (LittleBlockType::rank),
440 "managed_non_const_type::rank != LittleBlock::rank. " 441 "Please report this bug to the Ifpack2 developers.");
446 template<
class MatrixType>
449 typedef impl_scalar_type IST;
450 const char prefix[] =
"Ifpack2::Experimental::RBILUK::compute: ";
455 TEUCHOS_TEST_FOR_EXCEPTION
456 (A_block_.is_null (), std::runtime_error, prefix <<
"The matrix (A_block_, " 457 "the BlockCrsMatrix) is null. Please call setMatrix() with a nonnull " 458 "input before calling this method.");
459 TEUCHOS_TEST_FOR_EXCEPTION
460 (! A_block_->isFillComplete (), std::runtime_error, prefix <<
"The matrix " 461 "(A_block_, the BlockCrsMatrix) is not fill complete. You may not invoke " 462 "initialize() or compute() with this matrix until the matrix is fill " 463 "complete. Note: BlockCrsMatrix is fill complete if and only if its " 464 "underlying graph is fill complete.");
466 if (! this->isInitialized ()) {
473 if (! A_block_.is_null ()) {
474 Teuchos::RCP<block_crs_matrix_type> A_nc =
475 Teuchos::rcp_const_cast<block_crs_matrix_type> (A_block_);
476 A_nc->template sync<Kokkos::HostSpace> ();
478 L_block_->template sync<Kokkos::HostSpace> ();
479 U_block_->template sync<Kokkos::HostSpace> ();
480 D_block_->template sync<Kokkos::HostSpace> ();
482 L_block_->template modify<Kokkos::HostSpace> ();
483 U_block_->template modify<Kokkos::HostSpace> ();
484 D_block_->template modify<Kokkos::HostSpace> ();
486 Teuchos::Time timer (
"RBILUK::compute");
488 Teuchos::TimeMonitor timeMon (timer);
489 this->isComputed_ =
false;
496 initAllValues (*A_block_);
502 const size_t MaxNumEntries =
503 L_block_->getNodeMaxNumRowEntries () + U_block_->getNodeMaxNumRowEntries () + 1;
512 Teuchos::Array<int> ipiv_teuchos(blockSize_);
513 Kokkos::View<
int*, Kokkos::HostSpace,
514 Kokkos::MemoryUnmanaged> ipiv (ipiv_teuchos.getRawPtr (), blockSize_);
515 Teuchos::Array<IST> work_teuchos(blockSize_);
516 Kokkos::View<IST*, Kokkos::HostSpace,
517 Kokkos::MemoryUnmanaged> work (work_teuchos.getRawPtr (), blockSize_);
519 size_t num_cols = U_block_->getColMap()->getNodeNumElements();
520 Teuchos::Array<int> colflag(num_cols);
522 typename GetManagedView<little_block_type>::managed_non_const_type diagModBlock (
"diagModBlock", blockSize_, blockSize_);
523 typename GetManagedView<little_block_type>::managed_non_const_type matTmp (
"matTmp", blockSize_, blockSize_);
524 typename GetManagedView<little_block_type>::managed_non_const_type multiplier (
"multiplier", blockSize_, blockSize_);
532 for (
size_t j = 0; j < num_cols; ++j) {
535 Teuchos::Array<local_ordinal_type> InI(MaxNumEntries, 0);
536 Teuchos::Array<scalar_type> InV(MaxNumEntries*blockMatSize,STM::zero());
543 NumIn = MaxNumEntries;
547 L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
551 little_block_type lmat((
typename little_block_type::value_type*) &valsL[matOffset],blockSize_,rowStride);
552 little_block_type lmatV((
typename little_block_type::value_type*) &InV[matOffset],blockSize_,rowStride);
554 Tpetra::Experimental::COPY (lmat, lmatV);
555 InI[j] = colValsL[j];
558 little_block_type dmat = D_block_->getLocalBlock(local_row, local_row);
559 little_block_type dmatV((
typename little_block_type::value_type*) &InV[NumL*blockMatSize], blockSize_, rowStride);
561 Tpetra::Experimental::COPY (dmat, dmatV);
562 InI[NumL] = local_row;
566 U_block_->getLocalRowView(local_row, colValsU, valsU, NumURead);
570 if (!(colValsU[j] < numLocalRows))
continue;
571 InI[NumL+1+j] = colValsU[j];
573 little_block_type umat((
typename little_block_type::value_type*) &valsU[blockMatSize*j], blockSize_, rowStride);
574 little_block_type umatV((
typename little_block_type::value_type*) &InV[matOffset], blockSize_, rowStride);
576 Tpetra::Experimental::COPY (umat, umatV);
582 for (
size_t j = 0; j < NumIn; ++j) {
586 #ifndef IFPACK2_RBILUK_INITIAL 590 diagModBlock(i,j) = 0;
595 Kokkos::deep_copy (diagModBlock, diagmod);
600 little_block_type currentVal((
typename little_block_type::value_type*) &InV[jj*blockMatSize], blockSize_, rowStride);
602 Tpetra::Experimental::COPY (currentVal, multiplier);
604 const little_block_type dmatInverse = D_block_->getLocalBlock(j,j);
606 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 607 KokkosBatched::Experimental::SerialGemm
608 <KokkosBatched::Experimental::Trans::NoTranspose,
609 KokkosBatched::Experimental::Trans::NoTranspose,
610 KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
611 (STS::one (), currentVal, dmatInverse, STS::zero (), matTmp);
613 Tpetra::Experimental::GEMM (
"N",
"N", STS::one (), currentVal, dmatInverse,
614 STS::zero (), matTmp);
618 Tpetra::Experimental::COPY (matTmp, currentVal);
622 U_block_->getLocalRowView(j, UUI, UUV, NumUU);
624 if (this->RelaxValue_ == STM::zero ()) {
626 if (!(UUI[k] < numLocalRows))
continue;
627 const int kk = colflag[UUI[k]];
629 little_block_type kkval((
typename little_block_type::value_type*) &InV[kk*blockMatSize], blockSize_, rowStride);
630 little_block_type uumat((
typename little_block_type::value_type*) &UUV[k*blockMatSize], blockSize_, rowStride);
631 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 632 KokkosBatched::Experimental::SerialGemm
633 <KokkosBatched::Experimental::Trans::NoTranspose,
634 KokkosBatched::Experimental::Trans::NoTranspose,
635 KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
636 (
magnitude_type(-STM::one ()), multiplier, uumat, STM::one (), kkval);
638 Tpetra::Experimental::GEMM (
"N",
"N",
magnitude_type(-STM::one ()), multiplier, uumat,
647 if (!(UUI[k] < numLocalRows))
continue;
648 const int kk = colflag[UUI[k]];
649 little_block_type uumat((
typename little_block_type::value_type*) &UUV[k*blockMatSize], blockSize_, rowStride);
651 little_block_type kkval((
typename little_block_type::value_type*) &InV[kk*blockMatSize], blockSize_, rowStride);
652 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 653 KokkosBatched::Experimental::SerialGemm
654 <KokkosBatched::Experimental::Trans::NoTranspose,
655 KokkosBatched::Experimental::Trans::NoTranspose,
656 KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
657 (
magnitude_type(-STM::one ()), multiplier, uumat, STM::one (), kkval);
659 Tpetra::Experimental::GEMM (
"N",
"N",
magnitude_type(-STM::one ()), multiplier, uumat,
665 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 666 KokkosBatched::Experimental::SerialGemm
667 <KokkosBatched::Experimental::Trans::NoTranspose,
668 KokkosBatched::Experimental::Trans::NoTranspose,
669 KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
670 (
magnitude_type(-STM::one ()), multiplier, uumat, STM::one (), diagModBlock);
672 Tpetra::Experimental::GEMM (
"N",
"N",
magnitude_type(-STM::one ()), multiplier, uumat,
673 STM::one (), diagModBlock);
682 L_block_->replaceLocalValues (local_row, InI.getRawPtr (), InV.getRawPtr (), NumL);
686 Tpetra::Experimental::COPY (dmatV, dmat);
688 if (this->RelaxValue_ != STM::zero ()) {
690 Tpetra::Experimental::AXPY (this->RelaxValue_, diagModBlock, dmat);
704 for (
int k = 0; k < blockSize_; ++k) {
708 Tpetra::Experimental::GETF2 (dmat, ipiv, lapackInfo);
710 TEUCHOS_TEST_FOR_EXCEPTION(
711 lapackInfo != 0, std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 712 "lapackInfo = " << lapackInfo <<
" which indicates an error in the factorization GETRF.");
714 Tpetra::Experimental::GETRI (dmat, ipiv, work, lapackInfo);
716 TEUCHOS_TEST_FOR_EXCEPTION(
717 lapackInfo != 0, std::runtime_error,
"Ifpack2::Experimental::RBILUK::compute: " 718 "lapackInfo = " << lapackInfo <<
" which indicates an error in the matrix inverse GETRI.");
722 little_block_type currentVal((
typename little_block_type::value_type*) &InV[(NumL+1+j)*blockMatSize], blockSize_, rowStride);
724 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 725 KokkosBatched::Experimental::SerialGemm
726 <KokkosBatched::Experimental::Trans::NoTranspose,
727 KokkosBatched::Experimental::Trans::NoTranspose,
728 KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
729 (STS::one (), dmat, currentVal, STS::zero (), matTmp);
731 Tpetra::Experimental::GEMM (
"N",
"N", STS::one (), dmat, currentVal,
732 STS::zero (), matTmp);
736 Tpetra::Experimental::COPY (matTmp, currentVal);
741 U_block_->replaceLocalValues (local_row, &InI[NumL+1], &InV[blockMatSize*(NumL+1)], NumU);
744 #ifndef IFPACK2_RBILUK_INITIAL 746 for (
size_t j = 0; j < NumIn; ++j) {
747 colflag[InI[j]] = -1;
751 for (
size_t j = 0; j < num_cols; ++j) {
761 typedef typename block_crs_matrix_type::device_type device_type;
762 if (! A_block_.is_null ()) {
763 Teuchos::RCP<block_crs_matrix_type> A_nc =
764 Teuchos::rcp_const_cast<block_crs_matrix_type> (A_block_);
765 A_nc->template sync<device_type> ();
767 L_block_->template sync<device_type> ();
768 U_block_->template sync<device_type> ();
769 D_block_->template sync<device_type> ();
772 this->isComputed_ =
true;
773 this->numCompute_ += 1;
774 this->computeTime_ += timer.totalElapsedTime ();
778 template<
class MatrixType>
781 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
782 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
783 Teuchos::ETransp mode,
788 typedef Tpetra::Experimental::BlockMultiVector<
scalar_type,
793 TEUCHOS_TEST_FOR_EXCEPTION(
794 A_block_.is_null (), std::runtime_error,
"Ifpack2::Experimental::RBILUK::apply: The matrix is " 795 "null. Please call setMatrix() with a nonnull input, then initialize() " 796 "and compute(), before calling this method.");
797 TEUCHOS_TEST_FOR_EXCEPTION(
798 ! this->isComputed (), std::runtime_error,
799 "Ifpack2::Experimental::RBILUK::apply: If you have not yet called compute(), " 800 "you must call compute() before calling this method.");
801 TEUCHOS_TEST_FOR_EXCEPTION(
802 X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
803 "Ifpack2::Experimental::RBILUK::apply: X and Y do not have the same number of columns. " 804 "X.getNumVectors() = " << X.getNumVectors ()
805 <<
" != Y.getNumVectors() = " << Y.getNumVectors () <<
".");
806 TEUCHOS_TEST_FOR_EXCEPTION(
807 STS::isComplex && mode == Teuchos::CONJ_TRANS, std::logic_error,
808 "Ifpack2::Experimental::RBILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for " 809 "complex Scalar type. Please talk to the Ifpack2 developers to get this " 810 "fixed. There is a FIXME in this file about this very issue.");
816 BMV yBlock (Y, * (A_block_->getGraph ()->getDomainMap ()), blockSize_);
817 const BMV xBlock (X, * (A_block_->getColMap ()), blockSize_);
819 Teuchos::Array<scalar_type> lclarray(blockSize_);
820 little_vec_type lclvec((
typename little_vec_type::value_type*)&lclarray[0], blockSize_);
824 Teuchos::Time timer (
"RBILUK::apply");
826 Teuchos::TimeMonitor timeMon (timer);
827 if (alpha == one && beta == zero) {
828 if (mode == Teuchos::NO_TRANS) {
836 BMV cBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
837 BMV rBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
840 for (
size_t i = 0; i < D_block_->getNodeNumRows(); ++i)
843 little_vec_type xval = xBlock.getLocalBlock(local_row,imv);
844 little_vec_type cval = cBlock.getLocalBlock(local_row,imv);
846 Tpetra::Experimental::COPY (xval, cval);
852 L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
857 little_vec_type prevVal = cBlock.getLocalBlock(col, imv);
860 little_block_type lij((
typename little_block_type::value_type*) &valsL[matOffset],blockSize_,rowStride);
863 Tpetra::Experimental::GEMV (-one, lij, prevVal, cval);
869 D_block_->applyBlock(cBlock, rBlock);
878 little_vec_type rval = rBlock.getLocalBlock(local_row,imv);
879 little_vec_type yval = yBlock.getLocalBlock(local_row,imv);
881 Tpetra::Experimental::COPY (rval, yval);
887 U_block_->getLocalRowView(local_row, colValsU, valsU, NumU);
892 little_vec_type prevVal = yBlock.getLocalBlock(col, imv);
895 little_block_type uij((
typename little_block_type::value_type*) &valsU[matOffset], blockSize_, rowStride);
898 Tpetra::Experimental::GEMV (-one, uij, prevVal, yval);
904 TEUCHOS_TEST_FOR_EXCEPTION(
905 true, std::runtime_error,
906 "Ifpack2::Experimental::RBILUK::apply: transpose apply is not implemented for the block algorithm. ");
917 MV Y_tmp (Y.getMap (), Y.getNumVectors ());
918 apply (X, Y_tmp, mode);
919 Y.update (alpha, Y_tmp, beta);
924 this->numApply_ += 1;
925 this->applyTime_ = timer.totalElapsedTime ();
929 template<
class MatrixType>
932 std::ostringstream os;
937 os <<
"\"Ifpack2::Experimental::RBILUK\": {";
938 os <<
"Initialized: " << (this->isInitialized () ?
"true" :
"false") <<
", " 939 <<
"Computed: " << (this->isComputed () ?
"true" :
"false") <<
", ";
941 os <<
"Level-of-fill: " << this->getLevelOfFill() <<
", ";
943 if (A_block_.is_null ()) {
944 os <<
"Matrix: null";
947 os <<
"Global matrix dimensions: [" 948 << A_block_->getGlobalNumRows () <<
", " << A_block_->getGlobalNumCols () <<
"]" 949 <<
", Global nnz: " << A_block_->getGlobalNumEntries();
964 #define IFPACK2_EXPERIMENTAL_RBILUK_INSTANT(S,LO,GO,N) \ 965 template class Ifpack2::Experimental::RBILUK< Tpetra::BlockCrsMatrix<S, LO, GO, N> >; \ 966 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:125
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:187
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const
Get the input matrix.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:182
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:81
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:136
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:243
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:781
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:142
void compute()
Compute the (numeric) incomplete factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:447
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:139
Construct a level filled graph for use in computing an ILU(k) incomplete factorization.
Definition: Ifpack2_IlukGraph.hpp:97
void setMatrix(const Teuchos::RCP< const block_crs_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:86
Definition: Ifpack2_Container.hpp:774
std::string description() const
A one-line description of this object.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:930
Sparse matrix (Tpetra::RowMatrix subclass) with ghost rows.
Definition: Ifpack2_OverlappingRowMatrix_decl.hpp:58
const block_crs_matrix_type & getLBlock() const
Return the L factor of the ILU factorization.
Definition: Ifpack2_Experimental_RBILUK_def.hpp:111
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::BlockCrsMatrix.
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:128