43 #ifndef IFPACK2_EXPERIMENTAL_CRSRBILUK_DEF_HPP 
   44 #define IFPACK2_EXPERIMENTAL_CRSRBILUK_DEF_HPP 
   46 #include "Tpetra_BlockMultiVector.hpp" 
   47 #include "Tpetra_BlockView.hpp" 
   48 #include "Ifpack2_OverlappingRowMatrix.hpp" 
   49 #include "Ifpack2_LocalFilter.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" 
   64 namespace Experimental {
 
   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&
 
  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&
 
  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&
 
  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>
 
  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 =
 
  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());
 
  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.");
 
  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();
 
  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;
 
  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_;
 
  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;
 
  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.");
 
  316   L_block_->setAllToScalar (STM::zero ()); 
 
  317   U_block_->setAllToScalar (STM::zero ());
 
  318   D_block_->setAllToScalar (STM::zero ()); 
 
  324   const_cast<block_crs_matrix_type&
> (A).sync_host ();
 
  325   L_block_->sync_host ();
 
  326   U_block_->sync_host ();
 
  327   D_block_->sync_host ();
 
  329   L_block_->modify_host ();
 
  330   U_block_->modify_host ();
 
  331   D_block_->modify_host ();
 
  333   RCP<const map_type> rowMap = L_block_->getRowMap ();
 
  343   for (
size_t myRow=0; myRow<A.getNodeNumRows(); ++myRow) {
 
  344     local_ordinal_type local_row = myRow;
 
  348     const local_ordinal_type * InI = 0;
 
  349     scalar_type * InV = 0;
 
  350     A.getLocalRowView(local_row, InI, InV, NumIn);
 
  358     for (local_ordinal_type j = 0; j < NumIn; ++j) {
 
  359       const local_ordinal_type k = InI[j];
 
  360       const local_ordinal_type blockOffset = blockMatSize*j;
 
  362       if (k == local_row) {
 
  365         for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
 
  366           diagValues[jj] = this->Rthresh_ * InV[blockOffset+jj] + IFPACK2_SGN(InV[blockOffset+jj]) * this->Athresh_;
 
  367         D_block_->replaceLocalValues(local_row, &InI[j], diagValues.getRawPtr(), 1);
 
  371           true, std::runtime_error, 
"Ifpack2::RILUK::initAllValues: current " 
  372           "GID k = " << k << 
" < 0.  I'm not sure why this is an error; it is " 
  373           "probably an artifact of the undocumented assumptions of the " 
  374           "original implementation (likely copied and pasted from Ifpack).  " 
  375           "Nevertheless, the code I found here insisted on this being an error " 
  376           "state, so I will throw an exception here.");
 
  378       else if (k < local_row) {
 
  380         const local_ordinal_type LBlockOffset = NumL*blockMatSize;
 
  381         for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
 
  382           LV[LBlockOffset+jj] = InV[blockOffset+jj];
 
  385       else if (Teuchos::as<size_t>(k) <= rowMap->getNodeNumElements()) {
 
  387         const local_ordinal_type UBlockOffset = NumU*blockMatSize;
 
  388         for (local_ordinal_type jj = 0; jj < blockMatSize; ++jj)
 
  389           UV[UBlockOffset+jj] = InV[blockOffset+jj];
 
  400       for (local_ordinal_type jj = 0; jj < blockSize_; ++jj)
 
  401         diagValues[jj*(blockSize_+1)] = this->Athresh_;
 
  402       D_block_->replaceLocalValues(local_row, &local_row, diagValues.getRawPtr(), 1);
 
  406       L_block_->replaceLocalValues(local_row, &LI[0], &LV[0], NumL);
 
  410       U_block_->replaceLocalValues(local_row, &UI[0], &UV[0], NumU);
 
  417     typedef typename block_crs_matrix_type::device_type device_type;
 
  418     const_cast<block_crs_matrix_type&
> (A).
template sync<device_type> ();
 
  419     L_block_->template sync<device_type> ();
 
  420     U_block_->template sync<device_type> ();
 
  421     D_block_->template sync<device_type> ();
 
  423   this->isInitialized_ = 
true;
 
  432 template<
class LittleBlockType>
 
  433 struct GetManagedView {
 
  434   static_assert (Kokkos::Impl::is_view<LittleBlockType>::value,
 
  435                  "LittleBlockType must be a Kokkos::View.");
 
  436   typedef Kokkos::View<
typename LittleBlockType::non_const_data_type,
 
  437                        typename LittleBlockType::array_layout,
 
  438                        typename LittleBlockType::device_type> managed_non_const_type;
 
  439   static_assert (static_cast<int> (managed_non_const_type::rank) ==
 
  440                  static_cast<int> (LittleBlockType::rank),
 
  441                  "managed_non_const_type::rank != LittleBlock::rank.  " 
  442                  "Please report this bug to the Ifpack2 developers.");
 
  447 template<
class MatrixType>
 
  450   typedef impl_scalar_type IST;
 
  451   const char prefix[] = 
"Ifpack2::Experimental::RBILUK::compute: ";
 
  457     (A_block_.is_null (), std::runtime_error, prefix << 
"The matrix (A_block_, " 
  458      "the BlockCrsMatrix) is null.  Please call setMatrix() with a nonnull " 
  459      "input before calling this method.");
 
  461     (! A_block_->isFillComplete (), std::runtime_error, prefix << 
"The matrix " 
  462      "(A_block_, the BlockCrsMatrix) is not fill complete.  You may not invoke " 
  463      "initialize() or compute() with this matrix until the matrix is fill " 
  464      "complete.  Note: BlockCrsMatrix is fill complete if and only if its " 
  465      "underlying graph is fill complete.");
 
  467   if (! this->isInitialized ()) {
 
  474   if (! A_block_.is_null ()) {
 
  476       Teuchos::rcp_const_cast<block_crs_matrix_type> (A_block_);
 
  479   L_block_->sync_host ();
 
  480   U_block_->sync_host ();
 
  481   D_block_->sync_host ();
 
  483   L_block_->modify_host ();
 
  484   U_block_->modify_host ();
 
  485   D_block_->modify_host ();
 
  490     this->isComputed_ = 
false;
 
  497     initAllValues (*A_block_);
 
  503     const size_t MaxNumEntries =
 
  504       L_block_->getNodeMaxNumRowEntries () + U_block_->getNodeMaxNumRowEntries () + 1;
 
  514     Kokkos::View<
int*, Kokkos::HostSpace,
 
  515       Kokkos::MemoryUnmanaged> ipiv (ipiv_teuchos.
getRawPtr (), blockSize_);
 
  517     Kokkos::View<IST*, Kokkos::HostSpace,
 
  518       Kokkos::MemoryUnmanaged> work (work_teuchos.getRawPtr (), blockSize_);
 
  520     size_t num_cols = U_block_->getColMap()->getNodeNumElements();
 
  523     typename GetManagedView<little_block_type>::managed_non_const_type diagModBlock (
"diagModBlock", blockSize_, blockSize_);
 
  524     typename GetManagedView<little_block_type>::managed_non_const_type matTmp (
"matTmp", blockSize_, blockSize_);
 
  525     typename GetManagedView<little_block_type>::managed_non_const_type multiplier (
"multiplier", blockSize_, blockSize_);
 
  533     for (
size_t j = 0; j < num_cols; ++j) {
 
  544       NumIn = MaxNumEntries;
 
  548       L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
 
  552         little_block_type lmat((
typename little_block_type::value_type*) &valsL[matOffset],blockSize_,rowStride);
 
  553         little_block_type lmatV((
typename little_block_type::value_type*) &InV[matOffset],blockSize_,rowStride);
 
  555         Tpetra::COPY (lmat, lmatV);
 
  556         InI[j] = colValsL[j];
 
  559       little_block_type dmat = D_block_->getLocalBlock(local_row, local_row);
 
  560       little_block_type dmatV((
typename little_block_type::value_type*) &InV[NumL*blockMatSize], blockSize_, rowStride);
 
  562       Tpetra::COPY (dmat, dmatV);
 
  563       InI[NumL] = local_row;
 
  567       U_block_->getLocalRowView(local_row, colValsU, valsU, NumURead);
 
  571         if (!(colValsU[j] < numLocalRows)) 
continue;
 
  572         InI[NumL+1+j] = colValsU[j];
 
  574         little_block_type umat((
typename little_block_type::value_type*) &valsU[blockMatSize*j], blockSize_, rowStride);
 
  575         little_block_type umatV((
typename little_block_type::value_type*) &InV[matOffset], blockSize_, rowStride);
 
  577         Tpetra::COPY (umat, umatV);
 
  583       for (
size_t j = 0; j < NumIn; ++j) {
 
  587 #ifndef IFPACK2_RBILUK_INITIAL 
  591             diagModBlock(i,j) = 0;
 
  596       Kokkos::deep_copy (diagModBlock, diagmod);
 
  601         little_block_type currentVal((
typename little_block_type::value_type*) &InV[jj*blockMatSize], blockSize_, rowStride); 
 
  603         Tpetra::COPY (currentVal, multiplier);
 
  605         const little_block_type dmatInverse = D_block_->getLocalBlock(j,j);
 
  607 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 
  608         KokkosBatched::Experimental::SerialGemm
 
  609           <KokkosBatched::Experimental::Trans::NoTranspose,
 
  610           KokkosBatched::Experimental::Trans::NoTranspose,
 
  611           KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
 
  612           (STS::one (), currentVal, dmatInverse, STS::zero (), matTmp);
 
  614         Tpetra::GEMM (
"N", 
"N", STS::one (), currentVal, dmatInverse,
 
  615                       STS::zero (), matTmp);
 
  619         Tpetra::COPY (matTmp, currentVal);
 
  623         U_block_->getLocalRowView(j, UUI, UUV, NumUU);
 
  625         if (this->RelaxValue_ == STM::zero ()) {
 
  627             if (!(UUI[k] < numLocalRows)) 
continue;
 
  628             const int kk = colflag[UUI[k]];
 
  630               little_block_type kkval((
typename little_block_type::value_type*) &InV[kk*blockMatSize], blockSize_, rowStride);
 
  631               little_block_type uumat((
typename little_block_type::value_type*) &UUV[k*blockMatSize], blockSize_, rowStride);
 
  632 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 
  633         KokkosBatched::Experimental::SerialGemm
 
  634           <KokkosBatched::Experimental::Trans::NoTranspose,
 
  635           KokkosBatched::Experimental::Trans::NoTranspose,
 
  636           KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
 
  637           ( 
magnitude_type(-STM::one ()), multiplier, uumat, STM::one (), kkval);
 
  639               Tpetra::GEMM (
"N", 
"N", 
magnitude_type(-STM::one ()), multiplier, uumat,
 
  648             if (!(UUI[k] < numLocalRows)) 
continue;
 
  649             const int kk = colflag[UUI[k]];
 
  650             little_block_type uumat((
typename little_block_type::value_type*) &UUV[k*blockMatSize], blockSize_, rowStride);
 
  652               little_block_type kkval((
typename little_block_type::value_type*) &InV[kk*blockMatSize], blockSize_, rowStride);
 
  653 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 
  654         KokkosBatched::Experimental::SerialGemm
 
  655           <KokkosBatched::Experimental::Trans::NoTranspose,
 
  656           KokkosBatched::Experimental::Trans::NoTranspose,
 
  657           KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
 
  658           (
magnitude_type(-STM::one ()), multiplier, uumat, STM::one (), kkval);
 
  660               Tpetra::GEMM (
"N", 
"N", 
magnitude_type(-STM::one ()), multiplier, uumat,
 
  666 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 
  667         KokkosBatched::Experimental::SerialGemm
 
  668           <KokkosBatched::Experimental::Trans::NoTranspose,
 
  669           KokkosBatched::Experimental::Trans::NoTranspose,
 
  670           KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
 
  671           (
magnitude_type(-STM::one ()), multiplier, uumat, STM::one (), diagModBlock);
 
  673               Tpetra::GEMM (
"N", 
"N", 
magnitude_type(-STM::one ()), multiplier, uumat,
 
  674                             STM::one (), diagModBlock);
 
  687       Tpetra::COPY (dmatV, dmat);
 
  689       if (this->RelaxValue_ != STM::zero ()) {
 
  691         Tpetra::AXPY (this->RelaxValue_, diagModBlock, dmat);
 
  705         for (
int k = 0; k < blockSize_; ++k) {
 
  709         Tpetra::GETF2 (dmat, ipiv, lapackInfo);
 
  712           lapackInfo != 0, std::runtime_error, 
"Ifpack2::Experimental::RBILUK::compute: " 
  713           "lapackInfo = " << lapackInfo << 
" which indicates an error in the factorization GETRF.");
 
  715         Tpetra::GETRI (dmat, ipiv, work, lapackInfo);
 
  718           lapackInfo != 0, std::runtime_error, 
"Ifpack2::Experimental::RBILUK::compute: " 
  719           "lapackInfo = " << lapackInfo << 
" which indicates an error in the matrix inverse GETRI.");
 
  723         little_block_type currentVal((
typename little_block_type::value_type*) &InV[(NumL+1+j)*blockMatSize], blockSize_, rowStride); 
 
  725 #ifndef IFPACK2_RBILUK_INITIAL_NOKK 
  726         KokkosBatched::Experimental::SerialGemm
 
  727           <KokkosBatched::Experimental::Trans::NoTranspose,
 
  728           KokkosBatched::Experimental::Trans::NoTranspose,
 
  729           KokkosBatched::Experimental::Algo::Gemm::Blocked>::invoke
 
  730           (STS::one (), dmat, currentVal, STS::zero (), matTmp);
 
  732         Tpetra::GEMM (
"N", 
"N", STS::one (), dmat, currentVal,
 
  733                                     STS::zero (), matTmp);
 
  737         Tpetra::COPY (matTmp, currentVal);
 
  742         U_block_->replaceLocalValues (local_row, &InI[NumL+1], &InV[blockMatSize*(NumL+1)], NumU);
 
  745 #ifndef IFPACK2_RBILUK_INITIAL 
  747       for (
size_t j = 0; j < NumIn; ++j) {
 
  748         colflag[InI[j]] = -1;
 
  752       for (
size_t j = 0; j < num_cols; ++j) {
 
  762     typedef typename block_crs_matrix_type::device_type device_type;
 
  763     if (! A_block_.is_null ()) {
 
  765         Teuchos::rcp_const_cast<block_crs_matrix_type> (A_block_);
 
  766       A_nc->template sync<device_type> ();
 
  768     L_block_->template sync<device_type> ();
 
  769     U_block_->template sync<device_type> ();
 
  770     D_block_->template sync<device_type> ();
 
  773   this->isComputed_ = 
true;
 
  774   this->numCompute_ += 1;
 
  779 template<
class MatrixType>
 
  782 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
 
  783        Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
 
  791   typedef Tpetra::MultiVector<scalar_type,
 
  792     local_ordinal_type, global_ordinal_type, 
node_type> MV;
 
  795     A_block_.is_null (), std::runtime_error, 
"Ifpack2::Experimental::RBILUK::apply: The matrix is " 
  796     "null.  Please call setMatrix() with a nonnull input, then initialize() " 
  797     "and compute(), before calling this method.");
 
  799     ! this->isComputed (), std::runtime_error,
 
  800     "Ifpack2::Experimental::RBILUK::apply: If you have not yet called compute(), " 
  801     "you must call compute() before calling this method.");
 
  802   TEUCHOS_TEST_FOR_EXCEPTION(
 
  803     X.getNumVectors () != Y.getNumVectors (), std::invalid_argument,
 
  804     "Ifpack2::Experimental::RBILUK::apply: X and Y do not have the same number of columns.  " 
  805     "X.getNumVectors() = " << X.getNumVectors ()
 
  806     << 
" != Y.getNumVectors() = " << Y.getNumVectors () << 
".");
 
  807   TEUCHOS_TEST_FOR_EXCEPTION(
 
  809     "Ifpack2::Experimental::RBILUK::apply: mode = Teuchos::CONJ_TRANS is not implemented for " 
  810     "complex Scalar type.  Please talk to the Ifpack2 developers to get this " 
  811     "fixed.  There is a FIXME in this file about this very issue.");
 
  813   const local_ordinal_type blockMatSize = blockSize_*blockSize_;
 
  815   const local_ordinal_type rowStride = blockSize_;
 
  817   BMV yBlock (Y, * (A_block_->getGraph ()->getDomainMap ()), blockSize_);
 
  818   const BMV xBlock (X, * (A_block_->getColMap ()), blockSize_);
 
  821   little_vec_type lclvec((
typename little_vec_type::value_type*)&lclarray[0], blockSize_);
 
  822   const scalar_type one = STM::one ();
 
  823   const scalar_type zero = STM::zero ();
 
  828     if (alpha == one && beta == zero) {
 
  836         const local_ordinal_type numVectors = xBlock.getNumVectors();
 
  837         BMV cBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
 
  838         BMV rBlock (* (A_block_->getGraph ()->getDomainMap ()), blockSize_, numVectors);
 
  839         for (local_ordinal_type imv = 0; imv < numVectors; ++imv)
 
  841           for (
size_t i = 0; i < D_block_->getNodeNumRows(); ++i)
 
  843             local_ordinal_type local_row = i;
 
  844             little_vec_type xval = xBlock.getLocalBlock(local_row,imv);
 
  845             little_vec_type cval = cBlock.getLocalBlock(local_row,imv);
 
  847             Tpetra::COPY (xval, cval);
 
  849             local_ordinal_type NumL;
 
  850             const local_ordinal_type * colValsL;
 
  853             L_block_->getLocalRowView(local_row, colValsL, valsL, NumL);
 
  855             for (local_ordinal_type j = 0; j < NumL; ++j)
 
  857               local_ordinal_type col = colValsL[j];
 
  858               little_vec_type prevVal = cBlock.getLocalBlock(col, imv);
 
  860               const local_ordinal_type matOffset = blockMatSize*j;
 
  861               little_block_type lij((
typename little_block_type::value_type*) &valsL[matOffset],blockSize_,rowStride);
 
  864               Tpetra::GEMV (-one, lij, prevVal, cval);
 
  870         D_block_->applyBlock(cBlock, rBlock);
 
  873         for (local_ordinal_type imv = 0; imv < numVectors; ++imv)
 
  875           const local_ordinal_type numRows = D_block_->getNodeNumRows();
 
  876           for (local_ordinal_type i = 0; i < numRows; ++i)
 
  878             local_ordinal_type local_row = (numRows-1)-i;
 
  879             little_vec_type rval = rBlock.getLocalBlock(local_row,imv);
 
  880             little_vec_type yval = yBlock.getLocalBlock(local_row,imv);
 
  882             Tpetra::COPY (rval, yval);
 
  884             local_ordinal_type NumU;
 
  885             const local_ordinal_type * colValsU;
 
  888             U_block_->getLocalRowView(local_row, colValsU, valsU, NumU);
 
  890             for (local_ordinal_type j = 0; j < NumU; ++j)
 
  892               local_ordinal_type col = colValsU[NumU-1-j];
 
  893               little_vec_type prevVal = yBlock.getLocalBlock(col, imv);
 
  895               const local_ordinal_type matOffset = blockMatSize*(NumU-1-j);
 
  896               little_block_type uij((
typename little_block_type::value_type*) &valsU[matOffset], blockSize_, rowStride);
 
  899               Tpetra::GEMV (-one, uij, prevVal, yval);
 
  905         TEUCHOS_TEST_FOR_EXCEPTION(
 
  906           true, std::runtime_error,
 
  907           "Ifpack2::Experimental::RBILUK::apply: transpose apply is not implemented for the block algorithm. ");
 
  918         MV Y_tmp (Y.getMap (), Y.getNumVectors ());
 
  919         apply (X, Y_tmp, mode);
 
  920         Y.update (alpha, Y_tmp, beta);
 
  925   this->numApply_ += 1;
 
  930 template<
class MatrixType>
 
  933   std::ostringstream os;
 
  938   os << 
"\"Ifpack2::Experimental::RBILUK\": {";
 
  939   os << 
"Initialized: " << (this->isInitialized () ? 
"true" : 
"false") << 
", " 
  940      << 
"Computed: " << (this->isComputed () ? 
"true" : 
"false") << 
", ";
 
  942   os << 
"Level-of-fill: " << this->getLevelOfFill() << 
", ";
 
  944   if (A_block_.is_null ()) {
 
  945     os << 
"Matrix: null";
 
  948     os << 
"Global matrix dimensions: [" 
  949        << A_block_->getGlobalNumRows () << 
", " << A_block_->getGlobalNumCols () << 
"]" 
  950        << 
", Global nnz: " << A_block_->getGlobalNumEntries();
 
  965 #define IFPACK2_EXPERIMENTAL_RBILUK_INSTANT(S,LO,GO,N)                            \ 
  966   template class Ifpack2::Experimental::RBILUK< Tpetra::BlockCrsMatrix<S, LO, GO, N> >; \ 
  967   template class Ifpack2::Experimental::RBILUK< Tpetra::RowMatrix<S, LO, GO, N> >; 
void initialize()
Initialize by computing the symbolic incomplete factorization. 
Definition: Ifpack2_Experimental_RBILUK_def.hpp:187
 
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
 
Teuchos::RCP< const block_crs_matrix_type > getBlockMatrix() const 
Get the input matrix. 
Definition: Ifpack2_Experimental_RBILUK_def.hpp:182
 
const block_crs_matrix_type & getUBlock() const 
Return the U factor of the ILU factorization. 
Definition: Ifpack2_Experimental_RBILUK_def.hpp:139
 
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType. 
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:136
 
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:782
 
#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:243
 
MatrixType::node_type node_type
The Node type used by the input MatrixType. 
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:148
 
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType. 
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:142
 
const block_crs_matrix_type & getDBlock() const 
Return the diagonal entries of the ILU factorization. 
Definition: Ifpack2_Experimental_RBILUK_def.hpp:125
 
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:448
 
File for utility functions. 
 
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
 
const block_crs_matrix_type & getLBlock() const 
Return the L factor of the ILU factorization. 
Definition: Ifpack2_Experimental_RBILUK_def.hpp:111
 
Sparse matrix (Tpetra::RowMatrix subclass) with ghost rows. 
Definition: Ifpack2_OverlappingRowMatrix_decl.hpp:58
 
double totalElapsedTime(bool readCurrentTime=false) const 
 
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
 
std::string description() const 
A one-line description of this object. 
Definition: Ifpack2_Experimental_RBILUK_def.hpp:931
 
ILU(k) factorization of a given Tpetra::BlockCrsMatrix. 
Definition: Ifpack2_Experimental_RBILUK_decl.hpp:128