Ifpack2 Templated Preconditioning Package  Version 1.0
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Ifpack2_Details_ChebyshevKernel_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Ifpack2: Templated Object-Oriented Algebraic Preconditioner Package
4 //
5 // Copyright 2009 NTESS and the Ifpack2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
15 
16 #ifndef IFPACK2_DETAILS_CHEBYSHEVKERNEL_DECL_HPP
17 #define IFPACK2_DETAILS_CHEBYSHEVKERNEL_DECL_HPP
18 
19 #include "Ifpack2_config.h"
20 #include "Tpetra_CrsMatrix_fwd.hpp"
21 #include "Tpetra_MultiVector_fwd.hpp"
22 #include "Tpetra_Operator_fwd.hpp"
23 #include "Tpetra_Vector_fwd.hpp"
24 #include "Tpetra_Export_fwd.hpp"
25 #include "Tpetra_Import_fwd.hpp"
26 #include "Teuchos_RCP.hpp"
27 #include <memory>
28 
29 namespace Ifpack2 {
30 namespace Details {
31 
44 template<class TpetraOperatorType>
46 private:
47  using SC = typename TpetraOperatorType::scalar_type;
48  using LO = typename TpetraOperatorType::local_ordinal_type;
49  using GO = typename TpetraOperatorType::global_ordinal_type;
50  using NT = typename TpetraOperatorType::node_type;
51 
52  using crs_matrix_type = Tpetra::CrsMatrix<SC, LO, GO, NT>;
53  using multivector_type = Tpetra::MultiVector<SC, LO, GO, NT>;
54  using operator_type = Tpetra::Operator<SC, LO, GO, NT>;
55  using vector_type = Tpetra::Vector<SC, LO, GO, NT>;
56 
57 public:
59  const bool useNativeSpMV=false);
60 
61  void
62  setMatrix (const Teuchos::RCP<const operator_type>& A);
63 
64  void
65  compute (multivector_type& W,
66  const SC& alpha,
67  vector_type& D_inv,
68  multivector_type& B,
69  multivector_type& X,
70  const SC& beta);
71 
72 private:
73  using import_type = Tpetra::Import<LO, GO, NT>;
74  using export_type = Tpetra::Export<LO, GO, NT>;
75 
80  std::unique_ptr<vector_type> X_colMap_;
81  std::unique_ptr<multivector_type> V1_;
82 
83  Teuchos::RCP<vector_type> W_vec_, B_vec_, X_vec_;
84 
85  // External override to not fuse operations into a single kernel
86  // And use native blas/SpMV operations
87  bool useNativeSpMV_;
88 
89  // Do the Import, if needed, and return the column Map version of X.
90  vector_type&
91  importVector (vector_type& X_domMap);
92 
93  bool canFuse (const multivector_type& B) const;
94 
95  void
96  unfusedCase (multivector_type& W,
97  const SC& alpha,
98  vector_type& D_inv,
99  multivector_type& B,
100  const operator_type& A,
101  multivector_type& X,
102  const SC& beta);
103 
104  void
105  fusedCase (vector_type& W,
106  const SC& alpha,
107  vector_type& D_inv,
108  vector_type& B,
109  const crs_matrix_type& A,
110  vector_type& X,
111  const SC& beta);
112 };
113 
114 } // namespace Details
115 } // namespace Ifpack2
116 
117 #endif // IFPACK2_DETAILS_CHEBYSHEVKERNEL_DECL_HPP
Compute scaled damped residual for Chebyshev.
Definition: Ifpack2_Details_ChebyshevKernel_decl.hpp:45