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.hpp"
52 #include "Tpetra_BlockCrsMatrix.hpp"
53 #include <type_traits>
54 #include <KokkosKernels_Handle.hpp>
57 #ifndef DOXYGEN_SHOULD_SKIP_THIS
61 template<
class TpetraOperatorType>
62 class ScaledDampedResidual;
72 #endif // DOXYGEN_SHOULD_SKIP_THIS
231 template<
class MatrixType>
234 typename MatrixType::scalar_type,
235 typename MatrixType::local_ordinal_type,
236 typename MatrixType::global_ordinal_type,
237 typename MatrixType::node_type>,
239 Tpetra::RowMatrix<typename MatrixType::scalar_type,
240 typename MatrixType::local_ordinal_type,
241 typename MatrixType::global_ordinal_type,
242 typename MatrixType::node_type> >
267 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.");
400 return isInitialized_;
460 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
461 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
485 applyMat (
const Tpetra::MultiVector<
597 typedef typename crs_matrix_type::local_matrix_type local_matrix_type;
598 typedef typename local_matrix_type::StaticCrsGraphType::row_map_type lno_row_view_t;
599 typedef typename local_matrix_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
600 typedef typename local_matrix_type::values_type scalar_nonzero_view_t;
601 typedef typename local_matrix_type::StaticCrsGraphType::device_type TemporaryWorkSpace;
602 typedef typename local_matrix_type::StaticCrsGraphType::device_type PersistentWorkSpace;
603 typedef typename local_matrix_type::StaticCrsGraphType::execution_space MyExecSpace;
604 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
605 <
typename lno_row_view_t::const_value_type,
local_ordinal_type,
typename scalar_nonzero_view_t::value_type,
606 MyExecSpace, TemporaryWorkSpace,PersistentWorkSpace > mt_kernel_handle_type;
630 void ApplyInverseJacobi(
631 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
632 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
635 void ApplyInverseJacobi_BlockCrsMatrix(
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;
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 ApplyInverseMTGS_CrsMatrix(
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 void ApplyInverseGS_RowMatrix(
652 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
653 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
657 ApplyInverseGS_CrsMatrix (
const crs_matrix_type& A,
658 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
659 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
663 ApplyInverseGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
664 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
665 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
668 void ApplyInverseSGS(
669 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
670 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
673 void ApplyInverseMTSGS_CrsMatrix(
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)
const;
678 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
679 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
680 const Tpetra::ESweepDirection direction)
const;
683 void ApplyInverseSGS_RowMatrix(
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;
689 ApplyInverseSGS_CrsMatrix (
const crs_matrix_type& A,
690 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
691 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
695 ApplyInverseSGS_BlockCrsMatrix (
const block_crs_matrix_type& A,
696 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
697 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y);
699 void computeBlockCrs ();
702 void updateCachedMultiVector(
const Teuchos::RCP<
const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >& map,
size_t numVecs)
const;
727 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
728 typename block_crs_matrix_type::device_type> block_diag_type;
729 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
730 typename block_crs_matrix_type::device_type,
731 Kokkos::MemoryUnmanaged> unmanaged_block_diag_type;
747 block_diag_type blockDiag_;
754 Details::RelaxationType PrecType_;
760 bool ZeroStartingSolution_;
770 bool fixTinyDiagEntries_;
772 bool checkDiagEntries_;
775 bool is_matrix_structurally_symmetric_;
778 bool ifpack2_dump_matrix_;
790 mutable int NumApply_;
792 double InitializeTime_;
796 mutable double ApplyTime_;
798 double ComputeFlops_;
800 mutable double ApplyFlops_;
807 size_t globalNumSmallDiagEntries_;
809 size_t globalNumZeroDiagEntries_;
811 size_t globalNumNegDiagEntries_;
823 Kokkos::View<size_t*, typename node_type::device_type> diagOffsets_;
830 bool savedDiagOffsets_;
832 bool hasBlockCrsMatrix_;
842 #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:444
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:486
double getApplyFlops() const
Total number of floating-point operations over all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:492
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:2577
Relaxation(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_Relaxation_def.hpp:216
void compute()
Compute the preconditioner ("numeric setup");.
Definition: Ifpack2_Relaxation_def.hpp:928
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_Relaxation_def.hpp:499
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix and vectors are distributed.
Definition: Ifpack2_Relaxation_def.hpp:401
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:196
int getNumApply() const
Total number of calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:462
double getComputeTime() const
Total time in seconds spent in all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:474
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:258
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:421
void setParameters(const Teuchos::ParameterList ¶ms)
Set the relaxation / preconditioner parameters.
Definition: Ifpack2_Relaxation_def.hpp:391
std::string description() const
A simple one-line description of this object.
Definition: Ifpack2_Relaxation_def.hpp:2522
int getNumInitialize() const
Total number of calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:450
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:407
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:434
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:255
double getInitializeTime() const
Total time in seconds spent in all calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:468
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:615
Teuchos::RCP< const row_matrix_type > getMatrix() const
The matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:412
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:252
virtual ~Relaxation()=default
Destructor.
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:249
int getNumCompute() const
Total number of calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:456
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_Relaxation_decl.hpp:399
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition: Ifpack2_Relaxation_decl.hpp:232
double getApplyTime() const
Total time in seconds spent in all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:480
void initialize()
Initialize the preconditioner ("symbolic setup").
Definition: Ifpack2_Relaxation_def.hpp:636
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:265
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:512
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_Relaxation_def.hpp:258
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Relaxation_decl.hpp:261