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>
53 #include "Ifpack2_Details_GaussSeidel.hpp"
54 #include "Ifpack2_Details_InverseDiagonalKernel.hpp"
56 #ifndef DOXYGEN_SHOULD_SKIP_THIS
60 template<
class TpetraOperatorType>
61 class ScaledDampedResidual;
71 #endif // DOXYGEN_SHOULD_SKIP_THIS
236 template<
class MatrixType>
239 typename MatrixType::scalar_type,
240 typename MatrixType::local_ordinal_type,
241 typename MatrixType::global_ordinal_type,
242 typename MatrixType::node_type>,
244 Tpetra::RowMatrix<typename MatrixType::scalar_type,
245 typename MatrixType::local_ordinal_type,
246 typename MatrixType::global_ordinal_type,
247 typename MatrixType::node_type> >
281 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.");
406 bool supportsZeroStartingSolution() {
return true; }
419 return isInitialized_;
479 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
480 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
504 applyMat (
const Tpetra::MultiVector<
600 typedef typename Kokkos::ArithTraits<scalar_type>::val_type impl_scalar_type;
617 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
618 typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
620 typedef typename crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
621 typedef typename crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
629 typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
630 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
631 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
632 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
633 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type TemporaryWorkSpace;
634 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type PersistentWorkSpace;
635 typedef typename local_matrix_device_type::StaticCrsGraphType::execution_space MyExecSpace;
636 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
637 <
typename lno_row_view_t::const_value_type,
local_ordinal_type,
typename scalar_nonzero_view_t::value_type,
638 MyExecSpace, TemporaryWorkSpace,PersistentWorkSpace > mt_kernel_handle_type;
645 using SerialGaussSeidel = Ifpack2::Details::GaussSeidel<scalar_type, local_ordinal_type, global_ordinal_type, node_type>;
669 void ApplyInverseRichardson(
670 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
671 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
674 void ApplyInverseJacobi(
675 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
676 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
679 void ApplyInverseJacobi_BlockCrsMatrix(
680 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
681 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
684 void ApplyInverseMTGS_CrsMatrix(
685 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
686 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
687 Tpetra::ESweepDirection direction)
const;
690 void ApplyInverseSerialGS(
691 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
692 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
693 Tpetra::ESweepDirection direction)
const;
696 void ApplyInverseSerialGS_CrsMatrix (
const crs_matrix_type& A,
697 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
698 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
699 Tpetra::ESweepDirection direction)
const;
702 void ApplyInverseSerialGS_RowMatrix(
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,
705 Tpetra::ESweepDirection direction)
const;
708 void ApplyInverseSerialGS_BlockCrsMatrix(
709 const block_crs_matrix_type& A,
710 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
711 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
712 Tpetra::ESweepDirection direction);
715 void ApplyInverseMTSGS_CrsMatrix(
716 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
717 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
720 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
721 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
722 const Tpetra::ESweepDirection direction)
const;
724 void computeBlockCrs ();
727 void updateCachedMultiVector(
const Teuchos::RCP<
const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >& map,
size_t numVecs)
const;
754 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
755 typename block_crs_matrix_type::device_type> block_diag_type;
756 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
757 typename block_crs_matrix_type::device_type,
758 Kokkos::MemoryUnmanaged> unmanaged_block_diag_type;
774 block_diag_type blockDiag_;
781 Details::RelaxationType PrecType_ = Ifpack2::Details::JACOBI;
787 bool ZeroStartingSolution_ =
true;
789 bool DoBackwardGS_ =
false;
791 bool DoL1Method_ =
false;
797 bool fixTinyDiagEntries_ =
false;
799 bool checkDiagEntries_ =
false;
801 int clusterSize_ = 1;
803 int longRowThreshold_ = 0;
805 KokkosGraph::ColoringAlgorithm mtColoringAlgorithm_ = KokkosGraph::COLORING_DEFAULT;
808 int NumOuterSweeps_ = 1;
810 int NumInnerSweeps_ = 1;
812 bool InnerSpTrsv_ =
false;
816 bool CompactForm_ =
false;
819 bool is_matrix_structurally_symmetric_ =
false;
822 bool ifpack2_dump_matrix_ =
false;
826 bool isInitialized_ =
false;
828 bool IsComputed_ =
false;
830 int NumInitialize_ = 0;
834 bool TimerForApply_ =
true;
836 mutable int NumApply_ = 0;
838 double InitializeTime_ = 0.0;
840 double ComputeTime_ = 0.0;
842 mutable double ApplyTime_ = 0.0;
844 double ComputeFlops_ = 0.0;
846 mutable double ApplyFlops_ = 0.0;
853 size_t globalNumSmallDiagEntries_ = 0;
855 size_t globalNumZeroDiagEntries_ = 0;
857 size_t globalNumNegDiagEntries_ = 0;
869 Kokkos::View<size_t*, typename node_type::device_type> diagOffsets_;
876 bool savedDiagOffsets_ =
false;
878 bool hasBlockCrsMatrix_ =
false;
888 #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:525
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:93
Tpetra::Operator< scalar_type, local_ordinal_type, global_ordinal_type, node_type > op_type
Tpetra::Operator specialization used by this class.
Definition: Ifpack2_Relaxation_decl.hpp:279
double getComputeFlops() const
Total number of floating-point operations over all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:567
double getApplyFlops() const
Total number of floating-point operations over all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:573
void setZeroStartingSolution(bool zeroStartingSolution)
Set this preconditioner's parameters.
Definition: Ifpack2_Relaxation_decl.hpp:408
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:2295
Relaxation(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_Relaxation_def.hpp:214
void compute()
Compute the preconditioner ("numeric setup");.
Definition: Ifpack2_Relaxation_def.hpp:1025
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_Relaxation_def.hpp:580
MatrixType::node_type::device_type device_type
The Kokkos device type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:266
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix and vectors are distributed.
Definition: Ifpack2_Relaxation_def.hpp:482
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:193
int getNumApply() const
Total number of calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:543
double getComputeTime() const
Total time in seconds spent in all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:555
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:263
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:502
void setParameters(const Teuchos::ParameterList ¶ms)
Set the relaxation / preconditioner parameters.
Definition: Ifpack2_Relaxation_def.hpp:472
std::string description() const
A simple one-line description of this object.
Definition: Ifpack2_Relaxation_def.hpp:2197
int getNumInitialize() const
Total number of calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:531
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:426
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:515
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:260
double getInitializeTime() const
Total time in seconds spent in all calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:549
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:706
Teuchos::RCP< const row_matrix_type > getMatrix() const
The matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:493
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:257
virtual ~Relaxation()=default
Destructor.
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:254
static scalar_type zero()
int getNumCompute() const
Total number of calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:537
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_Relaxation_decl.hpp:418
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition: Ifpack2_Relaxation_decl.hpp:237
double getApplyTime() const
Total time in seconds spent in all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:561
void initialize()
Initialize the preconditioner ("symbolic setup").
Definition: Ifpack2_Relaxation_def.hpp:727
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:273
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:593
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_Relaxation_def.hpp:226
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Relaxation_decl.hpp:269