Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
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
29namespace Ifpack2 {
30namespace Details {
31
44template <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:
58 ChebyshevKernel(const Teuchos::RCP<const operator_type>& A,
59 const bool useNativeSpMV = false);
60
61 void
62 setMatrix(const Teuchos::RCP<const operator_type>& A);
63
64 void setAuxiliaryVectors(size_t numVectors);
65
66 void
67 compute(multivector_type& W,
68 const SC& alpha,
69 vector_type& D_inv,
70 multivector_type& B,
71 multivector_type& X,
72 const SC& beta);
73
74 private:
75 using import_type = Tpetra::Import<LO, GO, NT>;
76 using export_type = Tpetra::Export<LO, GO, NT>;
77
78 Teuchos::RCP<const operator_type> A_op_;
79 Teuchos::RCP<const crs_matrix_type> A_crs_;
80 Teuchos::RCP<const import_type> imp_;
81 Teuchos::RCP<const export_type> exp_;
82 std::unique_ptr<multivector_type> X_colMap_;
83 std::unique_ptr<multivector_type> V1_;
84
85 Teuchos::RCP<vector_type> W_vec_, B_vec_, X_vec_;
86
87 // External override to not fuse operations into a single kernel
88 // And use native blas/SpMV operations
89 bool useNativeSpMV_;
90
91 // Do the Import, if needed, and return the column Map version of X.
92 multivector_type&
93 importVector(multivector_type& X_domMap);
94
95 bool canFuse(const multivector_type& B) const;
96
97 void
98 unfusedCase(multivector_type& W,
99 const SC& alpha,
100 vector_type& D_inv,
101 multivector_type& B,
102 const operator_type& A,
103 multivector_type& X,
104 const SC& beta);
105
106 void
107 fusedCase(multivector_type& W,
108 const SC& alpha,
109 multivector_type& D_inv,
110 multivector_type& B,
111 const crs_matrix_type& A,
112 multivector_type& X,
113 const SC& beta);
114};
115
116} // namespace Details
117} // namespace Ifpack2
118
119#endif // IFPACK2_DETAILS_CHEBYSHEVKERNEL_DECL_HPP
Compute scaled damped residual for Chebyshev.
Definition Ifpack2_Details_ChebyshevKernel_decl.hpp:45
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40