43 #ifndef IFPACK2_RELAXATION_DECL_HPP
44 #define IFPACK2_RELAXATION_DECL_HPP
48 #include "Ifpack2_Parameters.hpp"
49 #include "Tpetra_Vector.hpp"
50 #include "Teuchos_ScalarTraits.hpp"
51 #include "Tpetra_CrsMatrix_decl.hpp"
52 #include "Tpetra_Experimental_BlockCrsMatrix_decl.hpp"
53 #include <type_traits>
54 #include <KokkosKernels_Handle.hpp>
222 template<
class MatrixType>
225 typename MatrixType::local_ordinal_type,
226 typename MatrixType::global_ordinal_type,
227 typename MatrixType::node_type>,
229 typename MatrixType::local_ordinal_type,
230 typename MatrixType::global_ordinal_type,
231 typename MatrixType::node_type> >
250 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType
magnitude_type;
256 static_assert(std::is_same<MatrixType, row_matrix_type>::value,
"Ifpack2::Relaxation: Please use MatrixType = Tpetra::RowMatrix. This saves build times, library sizes, and executable sizes. Don't worry, this class still works with CrsMatrix and BlockCrsMatrix; those are both subclasses of RowMatrix.");
295 explicit Relaxation (
const Teuchos::RCP<const row_matrix_type>& A);
384 Teuchos::RCP<const Teuchos::ParameterList>
399 return isInitialized_;
444 setMatrix (
const Teuchos::RCP<const row_matrix_type>& A);
466 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
467 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
468 Teuchos::ETransp mode = Teuchos::NO_TRANS,
469 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
470 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero())
const;
473 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
477 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
496 applyMat (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
497 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
498 Teuchos::ETransp mode = Teuchos::NO_TRANS)
const;
505 Teuchos::RCP<const Teuchos::Comm<int> >
getComm()
const;
508 Teuchos::RCP<const row_matrix_type>
getMatrix ()
const;
572 describe (Teuchos::FancyOStream &out,
573 const Teuchos::EVerbosityLevel verbLevel =
574 Teuchos::Describable::verbLevel_default)
const;
581 typedef Teuchos::ScalarTraits<scalar_type> STS;
582 typedef Teuchos::ScalarTraits<magnitude_type> STM;
600 typedef typename crs_matrix_type::local_matrix_type local_matrix_type;
601 typedef typename local_matrix_type::StaticCrsGraphType::row_map_type lno_row_view_t;
602 typedef typename local_matrix_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
603 typedef typename local_matrix_type::values_type scalar_nonzero_view_t;
604 typedef typename local_matrix_type::StaticCrsGraphType::device_type TemporaryWorkSpace;
605 typedef typename local_matrix_type::StaticCrsGraphType::device_type PersistentWorkSpace;
606 typedef typename local_matrix_type::StaticCrsGraphType::execution_space MyExecSpace;
607 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
608 <
typename lno_row_view_t::const_value_type,
local_ordinal_type,
typename scalar_nonzero_view_t::value_type,
609 MyExecSpace, TemporaryWorkSpace,PersistentWorkSpace > mt_kernel_handle_type;
610 Teuchos::RCP<mt_kernel_handle_type> mtKernelHandle_;
630 void setParametersImpl (Teuchos::ParameterList& params);
633 void ApplyInverseJacobi(
634 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
635 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
638 void ApplyInverseJacobi_BlockCrsMatrix(
639 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
640 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
644 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
645 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
648 void ApplyInverseMTGS_CrsMatrix(
649 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
650 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
654 void ApplyInverseGS_RowMatrix(
655 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
656 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
660 ApplyInverseGS_CrsMatrix (
const crs_matrix_type& A,
661 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
662 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
666 ApplyInverseGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
667 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
668 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
671 void ApplyInverseSGS(
672 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
673 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
676 void ApplyInverseMTSGS_CrsMatrix(
677 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
678 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
681 const crs_matrix_type* crsMat,
682 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
683 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
684 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& D,
686 const Tpetra::ESweepDirection direction,
688 const bool zeroInitialGuess)
const;
691 void ApplyInverseSGS_RowMatrix(
692 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
693 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
697 ApplyInverseSGS_CrsMatrix (
const crs_matrix_type& A,
698 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
699 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
703 ApplyInverseSGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
704 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
705 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
707 void computeBlockCrs ();
710 void updateCachedMultiVector(
const Teuchos::RCP<
const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >& map,
size_t numVecs)
const;
723 mutable Teuchos::RCP<const Teuchos::ParameterList> validParams_;
726 Teuchos::RCP<const row_matrix_type> A_;
729 Teuchos::RCP<const Tpetra::Import<local_ordinal_type,global_ordinal_type,node_type> > Importer_;
731 Teuchos::RCP<Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> > Diagonal_;
733 mutable Teuchos::RCP<Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> > cachedMV_;
735 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
736 typename block_crs_matrix_type::device_type> block_diag_type;
737 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
738 typename block_crs_matrix_type::device_type,
739 Kokkos::MemoryUnmanaged> unmanaged_block_diag_type;
755 block_diag_type blockDiag_;
757 Teuchos::RCP<block_multivector_type> yBlockColumnPointMap_;
762 Details::RelaxationType PrecType_;
768 bool ZeroStartingSolution_;
778 bool fixTinyDiagEntries_;
780 bool checkDiagEntries_;
783 bool is_matrix_structurally_symmetric_;
786 bool ifpack2_dump_matrix_;
798 mutable int NumApply_;
800 double InitializeTime_;
804 mutable double ApplyTime_;
806 double ComputeFlops_;
808 mutable double ApplyFlops_;
815 size_t globalNumSmallDiagEntries_;
817 size_t globalNumZeroDiagEntries_;
819 size_t globalNumNegDiagEntries_;
831 Kokkos::View<size_t*, typename node_type::device_type> diagOffsets_;
838 bool savedDiagOffsets_;
840 bool hasBlockCrsMatrix_;
843 Teuchos::ArrayRCP<local_ordinal_type> localSmoothingIndices_;
850 #endif // IFPACK2_RELAXATION_DECL_HPP