10 #ifndef IFPACK2_DIAGONAL_DECL_HPP 
   11 #define IFPACK2_DIAGONAL_DECL_HPP 
   15 #include "Tpetra_CrsMatrix_decl.hpp" 
   16 #include <type_traits> 
   37 template <
class MatrixType>
 
   39                                                         typename MatrixType::local_ordinal_type,
 
   40                                                         typename MatrixType::global_ordinal_type,
 
   41                                                         typename MatrixType::node_type>,
 
   43                                                                                     typename MatrixType::local_ordinal_type,
 
   44                                                                                     typename MatrixType::global_ordinal_type,
 
   45                                                                                     typename MatrixType::node_type> > {
 
   47   typedef typename MatrixType::scalar_type scalar_type;
 
   48   typedef typename MatrixType::local_ordinal_type local_ordinal_type;
 
   49   typedef typename MatrixType::global_ordinal_type global_ordinal_type;
 
   50   typedef typename MatrixType::node_type node_type;
 
   54   typedef Tpetra::RowMatrix<scalar_type,
 
   60   static_assert(std::is_same<MatrixType, row_matrix_type>::value, 
"Ifpack2::Diagonal: The template parameter MatrixType must be a Tpetra::RowMatrix specialization.  Please don't use Tpetra::CrsMatrix (a subclass of Tpetra::RowMatrix) here anymore.  The constructor can take either a RowMatrix or a CrsMatrix just fine.");
 
   63   typedef Tpetra::CrsMatrix<scalar_type,
 
   69   typedef Tpetra::Vector<scalar_type,
 
   75   typedef Tpetra::Map<local_ordinal_type,
 
  121     return isInitialized_;
 
  171   apply(
const Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& X,
 
  172         Tpetra::MultiVector<scalar_type, local_ordinal_type, global_ordinal_type, node_type>& Y,
 
  253   typedef Kokkos::View<size_t*, typename node_type::device_type> offsets_type;
 
  254   offsets_type offsets_;
 
  256   double initializeTime_;
 
  258   mutable double applyTime_;
 
  262   mutable int numApply_;
 
  284 template <
class MatrixType, 
class VectorType>
 
  286                                                  typename MatrixType::local_ordinal_type,
 
  287                                                  typename MatrixType::global_ordinal_type,
 
  288                                                  typename MatrixType::node_type> > >
 
  290   typedef Tpetra::RowMatrix<
typename MatrixType::scalar_type,
 
  291                             typename MatrixType::local_ordinal_type,
 
  292                             typename MatrixType::global_ordinal_type,
 
  293                             typename MatrixType::node_type>
 
std::string description() const 
Return a one-line description of this object. 
Definition: Ifpack2_Diagonal_def.hpp:253
Tpetra::RowMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > row_matrix_type
Tpetra::RowMatrix specialization used by this class. 
Definition: Ifpack2_Diagonal_decl.hpp:58
Mix-in interface for preconditioners that can change their matrix after construction. 
Definition: Ifpack2_Details_CanChangeMatrix.hpp:60
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, putting the result in Y. 
Definition: Ifpack2_Diagonal_def.hpp:201
Diagonal preconditioner. 
Definition: Ifpack2_Diagonal_decl.hpp:38
bool isInitialized() const 
Returns true if the preconditioner has been successfully initialized. 
Definition: Ifpack2_Diagonal_decl.hpp:120
double getApplyFlops() const 
Return the number of flops for the application of the preconditioner. 
Teuchos::RCP< const map_type > getDomainMap() const 
The Tpetra::Map representing this operator's domain. 
Definition: Ifpack2_Diagonal_def.hpp:60
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
double getComputeFlops() const 
Return the number of flops in the computation phase. 
Diagonal(const Teuchos::RCP< const row_matrix_type > &A)
Constructor that takes a Tpetra::RowMatrix. 
Definition: Ifpack2_Diagonal_def.hpp:19
void compute()
Compute the preconditioner. 
Definition: Ifpack2_Diagonal_def.hpp:158
Interface for all Ifpack2 preconditioners. 
Definition: Ifpack2_Preconditioner.hpp:74
virtual ~Diagonal()
Destructor. 
Definition: Ifpack2_Diagonal_def.hpp:56
Teuchos::RCP< Ifpack2::Diagonal< Tpetra::RowMatrix< typename MatrixType::scalar_type, typename MatrixType::local_ordinal_type, typename MatrixType::global_ordinal_type, typename MatrixType::node_type > > > createDiagonalPreconditioner(const Teuchos::RCP< const VectorType > &invdiag)
Definition: Ifpack2_Diagonal_decl.hpp:289
double getComputeTime() const 
Return the time spent in compute(). 
Definition: Ifpack2_Diagonal_def.hpp:243
int getNumCompute() const 
Return the number of calls to compute(). 
Definition: Ifpack2_Diagonal_def.hpp:228
Declaration of interface for preconditioners that can change their matrix after construction. 
double getInitializeTime() const 
Return the time spent in initialize(). 
Definition: Ifpack2_Diagonal_def.hpp:238
static const EVerbosityLevel verbLevel_default
void setParameters(const Teuchos::ParameterList ¶ms)
Sets parameters on this object. 
Definition: Ifpack2_Diagonal_def.hpp:98
int getNumInitialize() const 
Return the number of calls to initialize(). 
Definition: Ifpack2_Diagonal_def.hpp:223
Teuchos::RCP< const row_matrix_type > getMatrix() const 
Return the communicator associated with this matrix operator. 
Definition: Ifpack2_Diagonal_decl.hpp:195
double getApplyTime() const 
Return the time spent in apply(). 
Definition: Ifpack2_Diagonal_def.hpp:248
virtual void setMatrix(const Teuchos::RCP< const row_matrix_type > &A)
Change the matrix to be preconditioned. 
Definition: Ifpack2_Diagonal_def.hpp:110
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Print the object with some verbosity level to an FancyOStream object. 
Definition: Ifpack2_Diagonal_def.hpp:279
Teuchos::RCP< const map_type > getRangeMap() const 
The Tpetra::Map representing this operator's range. 
Definition: Ifpack2_Diagonal_def.hpp:79
void initialize()
Initialize. 
Definition: Ifpack2_Diagonal_def.hpp:118
Tpetra::Map< local_ordinal_type, global_ordinal_type, node_type > map_type
Tpetra::Map specialization used by this class. 
Definition: Ifpack2_Diagonal_decl.hpp:78
Tpetra::CrsMatrix< scalar_type, local_ordinal_type, global_ordinal_type, node_type > crs_matrix_type
Tpetra::CrsMatrix specialization used by this class. 
Definition: Ifpack2_Diagonal_decl.hpp:60
int getNumApply() const 
Return the number of calls to apply(). 
Definition: Ifpack2_Diagonal_def.hpp:233
Tpetra::Vector< scalar_type, local_ordinal_type, global_ordinal_type, node_type > vector_type
Tpetra::Vector specialization used by this class. 
Definition: Ifpack2_Diagonal_decl.hpp:73
bool isComputed() const 
Return true if compute() has been called. 
Definition: Ifpack2_Diagonal_decl.hpp:128