41 #ifndef IFPACK2_RELAXATION_DECL_HPP
42 #define IFPACK2_RELAXATION_DECL_HPP
46 #include "Ifpack2_Parameters.hpp"
47 #include "Tpetra_Vector.hpp"
49 #include "Tpetra_CrsMatrix.hpp"
50 #include "Tpetra_BlockCrsMatrix.hpp"
51 #include <type_traits>
52 #include <KokkosKernels_Handle.hpp>
55 #ifndef DOXYGEN_SHOULD_SKIP_THIS
59 template<
class TpetraOperatorType>
60 class ScaledDampedResidual;
70 #endif // DOXYGEN_SHOULD_SKIP_THIS
235 template<
class MatrixType>
238 typename MatrixType::scalar_type,
239 typename MatrixType::local_ordinal_type,
240 typename MatrixType::global_ordinal_type,
241 typename MatrixType::node_type>,
243 Tpetra::RowMatrix<typename MatrixType::scalar_type,
244 typename MatrixType::local_ordinal_type,
245 typename MatrixType::global_ordinal_type,
246 typename MatrixType::node_type> >
271 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.");
405 return isInitialized_;
465 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
466 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
490 applyMat (
const Tpetra::MultiVector<
602 typedef typename crs_matrix_type::local_matrix_type local_matrix_type;
603 typedef typename local_matrix_type::StaticCrsGraphType::row_map_type lno_row_view_t;
604 typedef typename local_matrix_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
605 typedef typename local_matrix_type::values_type scalar_nonzero_view_t;
606 typedef typename local_matrix_type::StaticCrsGraphType::device_type TemporaryWorkSpace;
607 typedef typename local_matrix_type::StaticCrsGraphType::device_type PersistentWorkSpace;
608 typedef typename local_matrix_type::StaticCrsGraphType::execution_space MyExecSpace;
609 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
610 <
typename lno_row_view_t::const_value_type,
local_ordinal_type,
typename scalar_nonzero_view_t::value_type,
611 MyExecSpace, TemporaryWorkSpace,PersistentWorkSpace > mt_kernel_handle_type;
635 void ApplyInverseRichardson(
636 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
637 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
640 void ApplyInverseJacobi(
641 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
642 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
645 void ApplyInverseJacobi_BlockCrsMatrix(
646 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
647 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
651 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
652 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
655 void ApplyInverseMTGS_CrsMatrix(
656 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
657 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
661 void ApplyInverseGS_RowMatrix(
662 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
663 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
667 ApplyInverseGS_CrsMatrix (
const crs_matrix_type& A,
668 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
669 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
673 ApplyInverseGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
674 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
675 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
678 void ApplyInverseSGS(
679 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
680 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
683 void ApplyInverseMTSGS_CrsMatrix(
684 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
685 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
688 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
689 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
690 const Tpetra::ESweepDirection direction)
const;
693 void ApplyInverseSGS_RowMatrix(
694 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
695 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
699 ApplyInverseSGS_CrsMatrix (
const crs_matrix_type& A,
700 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
701 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
705 ApplyInverseSGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
706 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
707 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
709 void computeBlockCrs ();
712 void updateCachedMultiVector(
const Teuchos::RCP<
const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >& map,
size_t numVecs)
const;
739 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
740 typename block_crs_matrix_type::device_type> block_diag_type;
741 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
742 typename block_crs_matrix_type::device_type,
743 Kokkos::MemoryUnmanaged> unmanaged_block_diag_type;
759 block_diag_type blockDiag_;
766 int NumInnerSweeps_ = 1;
768 Details::RelaxationType PrecType_ = Ifpack2::Details::JACOBI;
774 bool ZeroStartingSolution_ =
true;
776 bool DoBackwardGS_ =
false;
778 bool DoL1Method_ =
false;
784 bool fixTinyDiagEntries_ =
false;
786 bool checkDiagEntries_ =
false;
788 bool InnerSpTrsv_ =
false;
790 int clusterSize_ = 1;
793 bool is_matrix_structurally_symmetric_ =
false;
796 bool ifpack2_dump_matrix_ =
false;
800 bool isInitialized_ =
false;
802 bool IsComputed_ =
false;
804 int NumInitialize_ = 0;
808 mutable int NumApply_ = 0;
810 double InitializeTime_ = 0.0;
812 double ComputeTime_ = 0.0;
814 mutable double ApplyTime_ = 0.0;
816 double ComputeFlops_ = 0.0;
818 mutable double ApplyFlops_ = 0.0;
825 size_t globalNumSmallDiagEntries_ = 0;
827 size_t globalNumZeroDiagEntries_ = 0;
829 size_t globalNumNegDiagEntries_ = 0;
841 Kokkos::View<size_t*, typename node_type::device_type> diagOffsets_;
848 bool savedDiagOffsets_ =
false;
850 bool hasBlockCrsMatrix_ =
false;
860 #endif // IFPACK2_RELAXATION_DECL_HPP
bool hasTransposeApply() const
Whether apply() and applyMat() let you apply the transpose or conjugate transpose.
Definition: Ifpack2_Relaxation_def.hpp:440
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
double getComputeFlops() const
Total number of floating-point operations over all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:482
double getApplyFlops() const
Total number of floating-point operations over all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:488
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object's attributes to the given output stream.
Definition: Ifpack2_Relaxation_def.hpp:2681
Relaxation(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_Relaxation_def.hpp:215
void compute()
Compute the preconditioner ("numeric setup");.
Definition: Ifpack2_Relaxation_def.hpp:937
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_Relaxation_def.hpp:495
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix and vectors are distributed.
Definition: Ifpack2_Relaxation_def.hpp:397
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:194
int getNumApply() const
Total number of calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:458
double getComputeTime() const
Total time in seconds spent in all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:470
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:262
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getDomainMap() const
Returns the Tpetra::Map object associated with the domain of this operator.
Definition: Ifpack2_Relaxation_def.hpp:417
void setParameters(const Teuchos::ParameterList ¶ms)
Set the relaxation / preconditioner parameters.
Definition: Ifpack2_Relaxation_def.hpp:387
std::string description() const
A simple one-line description of this object.
Definition: Ifpack2_Relaxation_def.hpp:2622
int getNumInitialize() const
Total number of calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:446
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:107
bool isComputed() const
Return true if compute() has been called.
Definition: Ifpack2_Relaxation_decl.hpp:412
Teuchos::RCP< const Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > > getRangeMap() const
Returns the Tpetra::Map object associated with the range of this operator.
Definition: Ifpack2_Relaxation_def.hpp:430
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:259
double getInitializeTime() const
Total time in seconds spent in all calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:464
Declaration of interface for preconditioners that can change their matrix after construction.
static const EVerbosityLevel verbLevel_default
void applyMat(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) const
Apply the input matrix to X, returning the result in Y.
Definition: Ifpack2_Relaxation_def.hpp:616
Teuchos::RCP< const row_matrix_type > getMatrix() const
The matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:408
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:256
virtual ~Relaxation()=default
Destructor.
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:253
static scalar_type zero()
int getNumCompute() const
Total number of calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:452
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_Relaxation_decl.hpp:404
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition: Ifpack2_Relaxation_decl.hpp:236
double getApplyTime() const
Total time in seconds spent in all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:476
void initialize()
Initialize the preconditioner ("symbolic setup").
Definition: Ifpack2_Relaxation_def.hpp:637
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:269
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 preconditioner to X, returning the result in Y.
Definition: Ifpack2_Relaxation_def.hpp:508
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_Relaxation_def.hpp:227
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Relaxation_decl.hpp:265