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> 55 #ifdef HAVE_IFPACK2_EXPERIMENTAL_KOKKOSKERNELS_FEATURES 56 #include <KokkosKernels_Handle.hpp> 225 template<
class MatrixType>
228 typename MatrixType::local_ordinal_type,
229 typename MatrixType::global_ordinal_type,
230 typename MatrixType::node_type>,
232 typename MatrixType::local_ordinal_type,
233 typename MatrixType::global_ordinal_type,
234 typename MatrixType::node_type> >
253 typedef typename Teuchos::ScalarTraits<scalar_type>::magnitudeType
magnitude_type;
256 typedef Tpetra::RowMatrix<scalar_type, local_ordinal_type,
259 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.");
298 explicit Relaxation (
const Teuchos::RCP<const row_matrix_type>& A);
384 void setParameters (
const Teuchos::ParameterList& params);
387 Teuchos::RCP<const Teuchos::ParameterList>
402 return isInitialized_;
447 setMatrix (
const Teuchos::RCP<const row_matrix_type>& A);
469 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
470 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
471 Teuchos::ETransp mode = Teuchos::NO_TRANS,
472 scalar_type alpha = Teuchos::ScalarTraits<scalar_type>::one(),
473 scalar_type beta = Teuchos::ScalarTraits<scalar_type>::zero())
const;
476 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
477 getDomainMap ()
const;
480 Teuchos::RCP<const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >
481 getRangeMap ()
const;
484 bool hasTransposeApply ()
const;
499 applyMat (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
500 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
501 Teuchos::ETransp mode = Teuchos::NO_TRANS)
const;
508 Teuchos::RCP<const Teuchos::Comm<int> > getComm()
const;
511 Teuchos::RCP<const row_matrix_type> getMatrix ()
const;
514 double getComputeFlops()
const;
517 double getApplyFlops()
const;
520 int getNumInitialize()
const;
523 int getNumCompute()
const;
526 int getNumApply()
const;
529 double getInitializeTime()
const;
532 double getComputeTime()
const;
535 double getApplyTime()
const;
538 size_t getNodeSmootherComplexity()
const;
550 std::string description ()
const;
575 describe (Teuchos::FancyOStream &out,
576 const Teuchos::EVerbosityLevel verbLevel =
577 Teuchos::Describable::verbLevel_default)
const;
584 typedef Teuchos::ScalarTraits<scalar_type> STS;
585 typedef Teuchos::ScalarTraits<magnitude_type> STM;
591 typedef Tpetra::CrsMatrix<scalar_type, local_ordinal_type,
592 global_ordinal_type, node_type> crs_matrix_type;
593 typedef Tpetra::Experimental::BlockCrsMatrix<scalar_type, local_ordinal_type,
594 global_ordinal_type, node_type> block_crs_matrix_type;
595 typedef Tpetra::Experimental::BlockMultiVector<scalar_type, local_ordinal_type,
596 global_ordinal_type, node_type> block_multivector_type;
598 #ifdef HAVE_IFPACK2_EXPERIMENTAL_KOKKOSKERNELS_FEATURES 604 typedef typename crs_matrix_type::local_matrix_type local_matrix_type;
605 typedef typename local_matrix_type::StaticCrsGraphType::row_map_type lno_row_view_t;
606 typedef typename local_matrix_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
607 typedef typename local_matrix_type::values_type scalar_nonzero_view_t;
608 typedef typename local_matrix_type::StaticCrsGraphType::device_type TemporaryWorkSpace;
609 typedef typename local_matrix_type::StaticCrsGraphType::device_type PersistentWorkSpace;
610 typedef typename local_matrix_type::StaticCrsGraphType::execution_space MyExecSpace;
611 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
612 <lno_row_view_t,lno_nonzero_view_t, scalar_nonzero_view_t,
613 MyExecSpace, TemporaryWorkSpace,PersistentWorkSpace > mt_kernel_handle_type;
614 Teuchos::RCP<mt_kernel_handle_type> mtKernelHandle_;
615 #endif // HAVE_IFPACK2_EXPERIMENTAL_KOKKOSKERNELS_FEATURES 635 void setParametersImpl (Teuchos::ParameterList& params);
638 void ApplyInverseJacobi(
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;
643 void ApplyInverseJacobi_BlockCrsMatrix(
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;
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;
653 void ApplyInverseMTGS_CrsMatrix(
654 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
655 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
659 void ApplyInverseGS_RowMatrix(
660 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
661 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
665 ApplyInverseGS_CrsMatrix (
const crs_matrix_type& A,
666 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
667 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
671 ApplyInverseGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
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);
676 void ApplyInverseSGS(
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 void ApplyInverseMTSGS_CrsMatrix(
682 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
683 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
686 const crs_matrix_type* crsMat,
687 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
688 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
689 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& D,
690 const scalar_type& dampingFactor,
691 const Tpetra::ESweepDirection direction,
693 const bool zeroInitialGuess)
const;
696 void ApplyInverseSGS_RowMatrix(
697 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
698 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
702 ApplyInverseSGS_CrsMatrix (
const crs_matrix_type& A,
703 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
704 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
708 ApplyInverseSGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
709 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
710 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
712 void computeBlockCrs ();
725 mutable Teuchos::RCP<const Teuchos::ParameterList> validParams_;
728 Teuchos::RCP<const row_matrix_type> A_;
730 Teuchos::RCP<Teuchos::Time> Time_;
732 Teuchos::RCP<const Tpetra::Import<local_ordinal_type,global_ordinal_type,node_type> > Importer_;
734 Teuchos::RCP<Tpetra::Vector<scalar_type,local_ordinal_type,global_ordinal_type,node_type> > Diagonal_;
737 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
738 typename block_crs_matrix_type::device_type> block_diag_type;
739 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
740 typename block_crs_matrix_type::device_type,
741 Kokkos::MemoryUnmanaged> unmanaged_block_diag_type;
757 block_diag_type blockDiag_;
759 Teuchos::RCP<block_multivector_type> yBlockColumnPointMap_;
764 Details::RelaxationType PrecType_;
766 scalar_type DampingFactor_;
770 bool ZeroStartingSolution_;
776 magnitude_type L1Eta_;
778 scalar_type MinDiagonalValue_;
780 bool fixTinyDiagEntries_;
782 bool checkDiagEntries_;
793 mutable int NumApply_;
795 double InitializeTime_;
799 mutable double ApplyTime_;
801 double ComputeFlops_;
803 mutable double ApplyFlops_;
806 magnitude_type globalMinMagDiagEntryMag_;
808 magnitude_type globalMaxMagDiagEntryMag_;
810 size_t globalNumSmallDiagEntries_;
812 size_t globalNumZeroDiagEntries_;
814 size_t globalNumNegDiagEntries_;
819 magnitude_type globalDiagNormDiff_;
826 Kokkos::View<size_t*, typename node_type::device_type> diagOffsets_;
833 bool savedDiagOffsets_;
835 bool hasBlockCrsMatrix_;
838 Teuchos::ArrayRCP<local_ordinal_type> localSmoothingIndices_;
845 #endif // IFPACK2_RELAXATION_DECL_HPP Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_Relaxation_decl.hpp:401
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:250
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:107
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:247
Declaration of interface for preconditioners that can change their matrix after construction.
bool isComputed() const
Return true if compute() has been called.
Definition: Ifpack2_Relaxation_decl.hpp:416
Definition: Ifpack2_Container.hpp:761
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:244
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:241
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition: Ifpack2_Relaxation_decl.hpp:226
Preconditioners and smoothers for Tpetra sparse matrices.
Definition: Ifpack2_AdditiveSchwarz_decl.hpp:72
void getValidParameters(Teuchos::ParameterList ¶ms)
Fills a list which contains all the parameters possibly used by Ifpack2.
Definition: Ifpack2_Parameters.cpp:50
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class.
Definition: Ifpack2_Relaxation_decl.hpp:257
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Relaxation_decl.hpp:253