10 #ifndef IFPACK2_RELAXATION_DECL_HPP
11 #define IFPACK2_RELAXATION_DECL_HPP
15 #include "Ifpack2_Parameters.hpp"
16 #include "Tpetra_Vector.hpp"
18 #include "Tpetra_CrsMatrix.hpp"
19 #include "Tpetra_BlockCrsMatrix.hpp"
20 #include <type_traits>
21 #include <KokkosKernels_Handle.hpp>
22 #include "Ifpack2_Details_GaussSeidel.hpp"
23 #include "Ifpack2_Details_InverseDiagonalKernel.hpp"
25 #ifndef DOXYGEN_SHOULD_SKIP_THIS
29 template<
class TpetraOperatorType>
30 class ScaledDampedResidual;
40 #endif // DOXYGEN_SHOULD_SKIP_THIS
205 template<
class MatrixType>
208 typename MatrixType::scalar_type,
209 typename MatrixType::local_ordinal_type,
210 typename MatrixType::global_ordinal_type,
211 typename MatrixType::node_type>,
213 Tpetra::RowMatrix<typename MatrixType::scalar_type,
214 typename MatrixType::local_ordinal_type,
215 typename MatrixType::global_ordinal_type,
216 typename MatrixType::node_type> >
250 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.");
375 bool supportsZeroStartingSolution() {
return true; }
388 return isInitialized_;
448 apply (
const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
449 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y,
473 applyMat (
const Tpetra::MultiVector<
569 typedef typename Kokkos::ArithTraits<scalar_type>::val_type impl_scalar_type;
586 typedef Tpetra::Map<local_ordinal_type, global_ordinal_type, node_type> map_type;
587 typedef Tpetra::Import<local_ordinal_type, global_ordinal_type, node_type> import_type;
589 typedef typename crs_matrix_type::nonconst_local_inds_host_view_type nonconst_local_inds_host_view_type;
590 typedef typename crs_matrix_type::nonconst_values_host_view_type nonconst_values_host_view_type;
598 typedef typename crs_matrix_type::local_matrix_device_type local_matrix_device_type;
599 typedef typename local_matrix_device_type::StaticCrsGraphType::row_map_type lno_row_view_t;
600 typedef typename local_matrix_device_type::StaticCrsGraphType::entries_type lno_nonzero_view_t;
601 typedef typename local_matrix_device_type::values_type scalar_nonzero_view_t;
602 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type TemporaryWorkSpace;
603 typedef typename local_matrix_device_type::StaticCrsGraphType::device_type PersistentWorkSpace;
604 typedef typename local_matrix_device_type::StaticCrsGraphType::execution_space MyExecSpace;
605 typedef typename KokkosKernels::Experimental::KokkosKernelsHandle
606 <
typename lno_row_view_t::const_value_type,
local_ordinal_type,
typename scalar_nonzero_view_t::value_type,
607 MyExecSpace, TemporaryWorkSpace,PersistentWorkSpace > mt_kernel_handle_type;
614 using SerialGaussSeidel = Ifpack2::Details::GaussSeidel<scalar_type, local_ordinal_type, global_ordinal_type, node_type>;
638 void ApplyInverseRichardson(
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(
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 ApplyInverseJacobi_BlockCrsMatrix(
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>& B,
655 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
656 Tpetra::ESweepDirection direction)
const;
659 void ApplyInverseSerialGS(
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,
662 Tpetra::ESweepDirection direction)
const;
665 void ApplyInverseSerialGS_CrsMatrix (
const crs_matrix_type& A,
666 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
667 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
668 Tpetra::ESweepDirection direction)
const;
671 void ApplyInverseSerialGS_RowMatrix(
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,
674 Tpetra::ESweepDirection direction)
const;
677 void ApplyInverseSerialGS_BlockCrsMatrix(
678 const block_crs_matrix_type& A,
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,
681 Tpetra::ESweepDirection direction);
684 void ApplyInverseMTSGS_CrsMatrix(
685 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
686 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& Y)
const;
689 const Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& B,
690 Tpetra::MultiVector<scalar_type,local_ordinal_type,global_ordinal_type,node_type>& X,
691 const Tpetra::ESweepDirection direction)
const;
693 void computeBlockCrs ();
696 void updateCachedMultiVector(
const Teuchos::RCP<
const Tpetra::Map<local_ordinal_type,global_ordinal_type,node_type> >& map,
size_t numVecs)
const;
723 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
724 typename block_crs_matrix_type::device_type> block_diag_type;
725 typedef Kokkos::View<
typename block_crs_matrix_type::impl_scalar_type***,
726 typename block_crs_matrix_type::device_type,
727 Kokkos::MemoryUnmanaged> unmanaged_block_diag_type;
743 block_diag_type blockDiag_;
750 Details::RelaxationType PrecType_ = Ifpack2::Details::JACOBI;
756 bool ZeroStartingSolution_ =
true;
758 bool DoBackwardGS_ =
false;
760 bool DoL1Method_ =
false;
766 bool fixTinyDiagEntries_ =
false;
768 bool checkDiagEntries_ =
false;
770 int clusterSize_ = 1;
772 int longRowThreshold_ = 0;
774 KokkosGraph::ColoringAlgorithm mtColoringAlgorithm_ = KokkosGraph::COLORING_DEFAULT;
777 int NumOuterSweeps_ = 1;
779 int NumInnerSweeps_ = 1;
781 bool InnerSpTrsv_ =
false;
785 bool CompactForm_ =
false;
788 bool is_matrix_structurally_symmetric_ =
false;
791 bool ifpack2_dump_matrix_ =
false;
795 bool isInitialized_ =
false;
797 bool IsComputed_ =
false;
799 int NumInitialize_ = 0;
803 bool TimerForApply_ =
true;
805 mutable int NumApply_ = 0;
807 double InitializeTime_ = 0.0;
809 double ComputeTime_ = 0.0;
811 mutable double ApplyTime_ = 0.0;
813 double ComputeFlops_ = 0.0;
815 mutable double ApplyFlops_ = 0.0;
822 size_t globalNumSmallDiagEntries_ = 0;
824 size_t globalNumZeroDiagEntries_ = 0;
826 size_t globalNumNegDiagEntries_ = 0;
838 Kokkos::View<size_t*, typename node_type::device_type> diagOffsets_;
845 bool savedDiagOffsets_ =
false;
847 bool hasBlockCrsMatrix_ =
false;
857 #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:494
Mix-in interface for preconditioners that can change their matrix after construction.
Definition: Ifpack2_Details_CanChangeMatrix.hpp:60
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:248
double getComputeFlops() const
Total number of floating-point operations over all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:536
double getApplyFlops() const
Total number of floating-point operations over all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:542
void setZeroStartingSolution(bool zeroStartingSolution)
Set this preconditioner's parameters.
Definition: Ifpack2_Relaxation_decl.hpp:377
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:2264
Relaxation(const Teuchos::RCP< const row_matrix_type > &A)
Constructor.
Definition: Ifpack2_Relaxation_def.hpp:183
void compute()
Compute the preconditioner ("numeric setup");.
Definition: Ifpack2_Relaxation_def.hpp:994
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Definition: Ifpack2_Relaxation_def.hpp:549
MatrixType::node_type::device_type device_type
The Kokkos device type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:235
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
The communicator over which the matrix and vectors are distributed.
Definition: Ifpack2_Relaxation_def.hpp:451
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:162
int getNumApply() const
Total number of calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:512
double getComputeTime() const
Total time in seconds spent in all calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:524
MatrixType::node_type node_type
The Node type used by the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:232
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:471
void setParameters(const Teuchos::ParameterList ¶ms)
Set the relaxation / preconditioner parameters.
Definition: Ifpack2_Relaxation_def.hpp:441
std::string description() const
A simple one-line description of this object.
Definition: Ifpack2_Relaxation_def.hpp:2166
int getNumInitialize() const
Total number of calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:500
Interface for all Ifpack2 preconditioners.
Definition: Ifpack2_Preconditioner.hpp:74
bool isComputed() const
Return true if compute() has been called.
Definition: Ifpack2_Relaxation_decl.hpp:395
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:484
MatrixType::global_ordinal_type global_ordinal_type
The type of global indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:229
double getInitializeTime() const
Total time in seconds spent in all calls to initialize().
Definition: Ifpack2_Relaxation_def.hpp:518
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:675
Teuchos::RCP< const row_matrix_type > getMatrix() const
The matrix to be preconditioned.
Definition: Ifpack2_Relaxation_def.hpp:462
MatrixType::local_ordinal_type local_ordinal_type
The type of local indices in the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:226
virtual ~Relaxation()=default
Destructor.
MatrixType::scalar_type scalar_type
The type of the entries of the input MatrixType.
Definition: Ifpack2_Relaxation_decl.hpp:223
static scalar_type zero()
int getNumCompute() const
Total number of calls to compute().
Definition: Ifpack2_Relaxation_def.hpp:506
bool isInitialized() const
Returns true if the preconditioner has been successfully initialized.
Definition: Ifpack2_Relaxation_decl.hpp:387
Relaxation preconditioners for Tpetra::RowMatrix and Tpetra::CrsMatrix sparse matrices.
Definition: Ifpack2_Relaxation_decl.hpp:206
double getApplyTime() const
Total time in seconds spent in all calls to apply().
Definition: Ifpack2_Relaxation_def.hpp:530
void initialize()
Initialize the preconditioner ("symbolic setup").
Definition: Ifpack2_Relaxation_def.hpp:696
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:242
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:562
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Return a list of all the parameters that this class accepts.
Definition: Ifpack2_Relaxation_def.hpp:195
Teuchos::ScalarTraits< scalar_type >::magnitudeType magnitude_type
The type of the magnitude (absolute value) of a matrix entry.
Definition: Ifpack2_Relaxation_decl.hpp:238