43 #ifndef IFPACK2_RELAXATION_DECL_HPP
44 #define IFPACK2_RELAXATION_DECL_HPP
48 #include "Ifpack2_Parameters.hpp"
49 #include "Tpetra_Vector.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> >
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.");
399 return isInitialized_;
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,
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,
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;
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 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
682 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
683 const Tpetra::ESweepDirection direction)
const;
686 void ApplyInverseSGS_RowMatrix(
687 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
688 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
692 ApplyInverseSGS_CrsMatrix (
const crs_matrix_type& A,
693 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
694 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
698 ApplyInverseSGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
699 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
700 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
702 void computeBlockCrs ();
705 void updateCachedMultiVector(
const Teuchos::RCP<
const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >& map,
size_t numVecs)
const;
730 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
731 typename block_crs_matrix_type::device_type> block_diag_type;
732 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
733 typename block_crs_matrix_type::device_type,
734 Kokkos::MemoryUnmanaged> unmanaged_block_diag_type;
750 block_diag_type blockDiag_;
757 Details::RelaxationType PrecType_;
763 bool ZeroStartingSolution_;
773 bool fixTinyDiagEntries_;
775 bool checkDiagEntries_;
778 bool is_matrix_structurally_symmetric_;
781 bool ifpack2_dump_matrix_;
793 mutable int NumApply_;
795 double InitializeTime_;
799 mutable double ApplyTime_;
801 double ComputeFlops_;
803 mutable double ApplyFlops_;
810 size_t globalNumSmallDiagEntries_;
812 size_t globalNumZeroDiagEntries_;
814 size_t globalNumNegDiagEntries_;
826 Kokkos::View<size_t*, typename node_type::device_type> diagOffsets_;
833 bool savedDiagOffsets_;
835 bool hasBlockCrsMatrix_;
845 #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:409
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:451
double getApplyFlops() const
Total number of floating-point operations over all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:457
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:2497
Relaxation(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_Relaxation_def.hpp:177
void compute()
Compute the preconditioner.
Definition: Ifpack2_Relaxation_def.hpp:883
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_Relaxation_def.hpp:464
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix and vectors are distributed.
Definition: Ifpack2_Relaxation_def.hpp:366
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:157
int getNumApply() const
Total number of calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:427
double getComputeTime() const
Total time in seconds spent in all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:439
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:247
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:386
void setParameters(const Teuchos::ParameterList ¶ms)
Set the relaxation / preconditioner parameters.
Definition: Ifpack2_Relaxation_def.hpp:356
std::string description() const
A simple one-line description of this object.
Definition: Ifpack2_Relaxation_def.hpp:2442
int getNumInitialize() const
Total number of calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:415
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:413
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:399
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:244
double getInitializeTime() const
Total time in seconds spent in all calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:433
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 preconditioner to X, returning the result in Y.
Definition: Ifpack2_Relaxation_def.hpp:584
Teuchos::RCP< const row_matrix_type > getMatrix() const
The matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:377
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:241
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:238
virtual ~Relaxation()
Destructor.
Definition: Ifpack2_Relaxation_def.hpp:218
int getNumCompute() const
Total number of calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:421
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_Relaxation_decl.hpp:398
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition: Ifpack2_Relaxation_decl.hpp:223
double getApplyTime() const
Total time in seconds spent in all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:445
void initialize()
Initialize the preconditioner.
Definition: Ifpack2_Relaxation_def.hpp:605
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:254
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:477
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_Relaxation_def.hpp:223
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Relaxation_decl.hpp:250