43 #ifndef IFPACK2_CONTAINER_HPP
44 #define IFPACK2_CONTAINER_HPP
49 #include "Ifpack2_ConfigDefs.hpp"
50 #include "Tpetra_RowMatrix.hpp"
51 #include "Teuchos_Describable.hpp"
52 #include <Ifpack2_Partitioner.hpp>
53 #include <Tpetra_Map.hpp>
54 #include <Tpetra_Experimental_BlockCrsMatrix_decl.hpp>
55 #include <Teuchos_ParameterList.hpp>
56 #include <Teuchos_Time.hpp>
58 #include <type_traits>
60 #ifndef DOXYGEN_SHOULD_SKIP_THIS
65 #endif // DOXYGEN_SHOULD_SKIP_THIS
113 template<
class MatrixType>
118 typedef typename MatrixType::scalar_type scalar_type;
119 typedef typename MatrixType::local_ordinal_type local_ordinal_type;
120 typedef typename MatrixType::global_ordinal_type global_ordinal_type;
121 typedef typename MatrixType::node_type node_type;
122 typedef Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> mv_type;
123 typedef Tpetra::Vector<scalar_type, local_ordinal_type, global_ordinal_type, node_type> vector_type;
124 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
125 typedef Teuchos::ScalarTraits<scalar_type> STS;
126 typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
128 typedef Tpetra::Experimental::BlockCrsMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> block_crs_matrix_type;
129 typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type, global_ordinal_type, node_type> row_matrix_type;
131 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
132 "Ifpack2::Container: Please use MatrixType = Tpetra::RowMatrix.");
135 typedef typename Kokkos::Details::ArithTraits<scalar_type>::val_type
impl_scalar_type;
139 typedef typename mv_type::dual_view_type::t_host HostView;
149 Container (
const Teuchos::RCP<const row_matrix_type>& matrix,
150 const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions,
151 const Teuchos::RCP<const import_type>& importer,
153 scalar_type DampingFactor) :
161 using Teuchos::Array;
162 using Teuchos::ArrayView;
170 Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(
inputMatrix_);
190 Container (
const Teuchos::RCP<const row_matrix_type>& matrix,
191 const Teuchos::Array<local_ordinal_type>& localRows) :
207 const block_crs_matrix_type* global_bcrs =
208 Teuchos::rcp_dynamic_cast<const block_crs_matrix_type>(
inputMatrix_).get();
214 for(local_ordinal_type i = 0; i < localRows.size(); i++)
236 Teuchos::ArrayView<const local_ordinal_type>
getLocalRows(
int blockIndex)
const
238 return Teuchos::ArrayView<const local_ordinal_type>
253 void setBlockSizes(
const Teuchos::Array<Teuchos::Array<local_ordinal_type> >& partitions)
257 local_ordinal_type totalBlockRows = 0;
263 local_ordinal_type rowsInBlock = partitions[i].size();
266 totalBlockRows += rowsInBlock;
270 local_ordinal_type iter = 0;
280 void getMatDiag()
const
296 void DoJacobi(HostView& X, HostView& Y,
int stride)
const;
297 void DoOverlappingJacobi(HostView& X, HostView& Y, HostView& W,
int stride)
const;
298 void DoGaussSeidel(HostView& X, HostView& Y, HostView& Y2,
int stride)
const;
299 void DoSGS(HostView& X, HostView& Y, HostView& Y2,
int stride)
const;
302 virtual void setParameters (
const Teuchos::ParameterList& List) = 0;
327 Teuchos::ETransp mode = Teuchos::NO_TRANS,
328 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
329 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero())
const = 0;
342 { TEUCHOS_TEST_FOR_EXCEPT_MSG(
true,
"Not implemented."); }
347 HostView XView = X.template getLocalView<Kokkos::HostSpace>();
348 HostView YView = Y.template getLocalView<Kokkos::HostSpace>();
349 this->
apply (XView, YView, 0, X.getStride());
378 Teuchos::ETransp mode = Teuchos::NO_TRANS,
379 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
380 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero())
const = 0;
387 HostView XView = X.template getLocalView<Kokkos::HostSpace>();
388 HostView YView = Y.template getLocalView<Kokkos::HostSpace>();
389 HostView WView = W.template getLocalView<Kokkos::HostSpace>();
393 virtual void clearBlocks();
396 virtual std::ostream&
print (std::ostream& os)
const = 0;
418 mutable Teuchos::RCP<vector_type>
Diag_;
426 Teuchos::RCP<const Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type>>
Importer_;
440 template <
class MatrixType>
444 return obj.
print (os);
447 template <
class MatrixType>
448 void Container<MatrixType>::DoJacobi(HostView& X, HostView& Y,
int stride)
const
450 const scalar_type one = STS::one();
453 size_t numVecs = X.extent(1);
455 for (local_ordinal_type i = 0; i < numBlocks_; i++)
458 if(blockRows_[i] != 1 || hasBlockCrs_)
460 if(blockRows_[i] == 0 )
462 apply(X, Y, i, stride, Teuchos::NO_TRANS, DampingFactor_, one);
466 local_ordinal_type LRID = partitions_[partitionIndices_[i]];
468 HostView diagView = Diag_->template getLocalView<Kokkos::HostSpace>();
469 impl_scalar_type d = (impl_scalar_type) one / diagView(LRID, 0);
470 for(
size_t nv = 0; nv < numVecs; nv++)
472 impl_scalar_type x = X(LRID, nv);
479 template <
class MatrixType>
480 void Container<MatrixType>::DoOverlappingJacobi(HostView& X, HostView& Y, HostView& W,
int stride)
const
483 for(local_ordinal_type i = 0; i < numBlocks_; i++)
486 if(blockRows_[i] == 0)
488 if(blockRows_[i] != 1)
489 weightedApply(X, Y, W, i, stride, Teuchos::NO_TRANS, DampingFactor_, STS::one());
493 template<
class MatrixType>
494 void Container<MatrixType>::DoGaussSeidel(HostView& X, HostView& Y, HostView& Y2,
int stride)
const
496 using Teuchos::Array;
497 using Teuchos::ArrayRCP;
498 using Teuchos::ArrayView;
502 using Teuchos::rcpFromRef;
504 const scalar_type one = STS::one();
505 const size_t Length = inputMatrix_->getNodeMaxNumRowEntries();
506 auto numVecs = X.extent(1);
507 Array<scalar_type> Values;
508 Array<local_ordinal_type> Indices;
509 Indices.resize(Length);
511 Values.resize(bcrsBlockSize_ * bcrsBlockSize_ * Length);
516 HostView Resid(
"", X.extent(0), X.extent(1));
517 for(local_ordinal_type i = 0; i < numBlocks_; i++)
519 if(blockRows_[i] > 1 || hasBlockCrs_)
521 if (blockRows_[i] == 0)
524 ArrayView<const local_ordinal_type> localRows = getLocalRows (i);
525 const size_t localNumRows = blockRows_[i];
526 for(
size_t j = 0; j < localNumRows; j++)
528 const local_ordinal_type LID = localRows[j];
530 inputMatrix_->getLocalRowCopy(LID, Indices(), Values(), NumEntries);
531 for(
size_t m = 0; m < numVecs; m++)
535 for (
int localR = 0; localR < bcrsBlockSize_; localR++)
536 Resid(LID * bcrsBlockSize_ + localR, m) = X(LID * bcrsBlockSize_ + localR, m);
537 for (
size_t k = 0; k < NumEntries; ++k)
539 const local_ordinal_type col = Indices[k];
540 for (
int localR = 0; localR < bcrsBlockSize_; localR++)
542 for(
int localC = 0; localC < bcrsBlockSize_; localC++)
544 Resid(LID * bcrsBlockSize_ + localR, m) -=
545 Values[k * bcrsBlockSize_ * bcrsBlockSize_ + localR + localC * bcrsBlockSize_]
546 * Y2(col * bcrsBlockSize_ + localC, m);
553 Resid(LID, m) = X(LID, m);
554 for (
size_t k = 0; k < NumEntries; ++k)
556 const local_ordinal_type col = Indices[k];
557 Resid(LID, m) -= Values[k] * Y2(col, m);
568 apply(Resid, Y2, i, stride, Teuchos::NO_TRANS, DampingFactor_, one);
570 else if(blockRows_[i] == 1)
574 local_ordinal_type LRID = partitionIndices_[i];
576 HostView diagView = Diag_->template getLocalView<Kokkos::HostSpace>();
577 impl_scalar_type d = (impl_scalar_type) one / diagView(LRID, 0);
578 for(
size_t m = 0; m < numVecs; m++)
580 impl_scalar_type x = X(LRID, m);
581 impl_scalar_type newy = x * d;
588 auto numMyRows = inputMatrix_->getNodeNumRows();
589 for (
size_t m = 0; m < numVecs; ++m)
591 for (
size_t i = 0; i < numMyRows * bcrsBlockSize_; ++i)
599 template<
class MatrixType>
600 void Container<MatrixType>::DoSGS(HostView& X, HostView& Y, HostView& Y2,
int stride)
const
602 using Teuchos::Array;
603 using Teuchos::ArrayRCP;
604 using Teuchos::ArrayView;
609 using Teuchos::rcpFromRef;
610 const scalar_type one = STS::one();
611 auto numVecs = X.extent(1);
612 const size_t Length = inputMatrix_->getNodeMaxNumRowEntries();
613 Array<scalar_type> Values;
614 Array<local_ordinal_type> Indices(Length);
615 Values.resize(bcrsBlockSize_ * bcrsBlockSize_ * Length);
620 HostView Resid(
"", X.extent(0), X.extent(1));
622 for(local_ordinal_type i = 0; i < numBlocks_; i++)
624 if(blockRows_[i] != 1 || hasBlockCrs_)
626 if(blockRows_[i] == 0)
629 ArrayView<const local_ordinal_type> localRows = getLocalRows(i);
630 for(local_ordinal_type j = 0; j < blockRows_[i]; j++)
632 const local_ordinal_type LID = localRows[j];
634 inputMatrix_->getLocalRowCopy(LID, Indices(), Values(), NumEntries);
636 for(
size_t m = 0; m < numVecs; m++)
640 for(
int localR = 0; localR < bcrsBlockSize_; localR++)
641 Resid(LID * bcrsBlockSize_ + localR, m) = X(LID * bcrsBlockSize_ + localR, m);
642 for(
size_t k = 0; k < NumEntries; ++k)
644 const local_ordinal_type col = Indices[k];
645 for (
int localR = 0; localR < bcrsBlockSize_; localR++)
647 for(
int localC = 0; localC < bcrsBlockSize_; localC++)
648 Resid(LID * bcrsBlockSize_ + localR, m) -=
649 Values[k * (bcrsBlockSize_ * bcrsBlockSize_) + (localR + localC * bcrsBlockSize_)]
650 * Y2(col * bcrsBlockSize_ + localC, m);
656 Resid(LID, m) = X(LID, m);
657 for(
size_t k = 0; k < NumEntries; k++)
659 local_ordinal_type col = Indices[k];
660 Resid(LID, m) -= Values[k] * Y2(col, m);
671 apply(Resid, Y2, i, stride, Teuchos::NO_TRANS, DampingFactor_, one);
676 local_ordinal_type LRID = partitions_[partitionIndices_[i]];
678 HostView diagView = Diag_->template getLocalView<Kokkos::HostSpace>();
679 impl_scalar_type d = (impl_scalar_type) one / diagView(LRID, 0);
680 for(
size_t m = 0; m < numVecs; m++)
682 impl_scalar_type x = X(LRID, m);
683 impl_scalar_type newy = x * d;
695 for(local_ordinal_type i = numBlocks_; i > 0; --i)
697 if(hasBlockCrs_ || blockRows_[i-1] != 1)
699 if(blockRows_[i - 1] == 0)
702 ArrayView<const local_ordinal_type> localRows = getLocalRows(i-1);
703 for(local_ordinal_type j = 0; j < blockRows_[i-1]; j++)
705 const local_ordinal_type LID = localRows[j];
707 inputMatrix_->getLocalRowCopy(LID, Indices(), Values(), NumEntries);
709 for (
size_t m = 0; m < numVecs; m++)
713 for(
int localR = 0; localR < bcrsBlockSize_; localR++)
714 Resid(LID * bcrsBlockSize_ + localR, m) = X(LID * bcrsBlockSize_ + localR, m);
715 for(
size_t k = 0; k < NumEntries; ++k)
717 const local_ordinal_type col = Indices[k];
718 for(
int localR = 0; localR < bcrsBlockSize_; localR++)
720 for(
int localC = 0; localC < bcrsBlockSize_; localC++)
721 Resid(LID*bcrsBlockSize_+localR, m) -=
722 Values[k * (bcrsBlockSize_ * bcrsBlockSize_) + (localR + localC * bcrsBlockSize_)]
723 * Y2(col * bcrsBlockSize_ + localC, m);
729 Resid(LID, m) = X(LID, m);
730 for(
size_t k = 0; k < NumEntries; ++k)
732 local_ordinal_type col = Indices[k];
733 Resid(LID, m) -= Values[k] * Y2(col, m);
744 apply(Resid, Y2, i - 1, stride, Teuchos::NO_TRANS, DampingFactor_, one);
751 auto numMyRows = inputMatrix_->getNodeNumRows();
752 for (
size_t m = 0; m < numVecs; ++m)
754 for (
size_t i = 0; i < numMyRows * bcrsBlockSize_; ++i)
762 template<
class MatrixType>
763 void Container<MatrixType>::clearBlocks()
768 partitionIndices_.clear();
769 Diag_ = Teuchos::null;
780 template<
class MatrixType>
784 static std::string name () {
785 return std::string (
"Ifpack2::Container<") +
786 TypeNameTraits<MatrixType>::name () +
" >";
789 static std::string concreteName (const ::Ifpack2::Container<MatrixType>&) {
796 #endif // IFPACK2_CONTAINER_HPP