Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_ScaledDampedResidual_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_SCALEDDAMPEDRESIDUAL_DECL_HPP
17#define IFPACK2_DETAILS_SCALEDDAMPEDRESIDUAL_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 ScaledDampedResidual(const Teuchos::RCP<const operator_type>& A);
59
60 void
61 setMatrix(const Teuchos::RCP<const operator_type>& A);
62
63 void
64 compute(multivector_type& W,
65 const SC& alpha,
66 vector_type& D_inv,
67 multivector_type& B,
68 multivector_type& X,
69 const SC& beta);
70
71 private:
72 using import_type = Tpetra::Import<LO, GO, NT>;
73 using export_type = Tpetra::Export<LO, GO, NT>;
74
75 Teuchos::RCP<const operator_type> A_op_;
76 Teuchos::RCP<const crs_matrix_type> A_crs_;
77 Teuchos::RCP<const import_type> imp_;
78 Teuchos::RCP<const export_type> exp_;
79 std::unique_ptr<vector_type> X_colMap_;
80 std::unique_ptr<multivector_type> V1_;
81
82 Teuchos::RCP<vector_type> W_vec_, B_vec_, X_vec_;
83
84 // Do the Import, if needed, and return the column Map version of X.
85 vector_type&
86 importVector(vector_type& X_domMap);
87
88 bool canFuse(const multivector_type& B) const;
89
90 void
91 unfusedCase(multivector_type& W,
92 const SC& alpha,
93 vector_type& D_inv,
94 multivector_type& B,
95 const operator_type& A,
96 multivector_type& X,
97 const SC& beta);
98
99 void
100 fusedCase(vector_type& W,
101 const SC& alpha,
102 vector_type& D_inv,
103 vector_type& B,
104 const crs_matrix_type& A,
105 vector_type& X,
106 const SC& beta);
107};
108
109} // namespace Details
110} // namespace Ifpack2
111
112#endif // IFPACK2_DETAILS_SCALEDDAMPEDRESIDUAL_DECL_HPP
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
Compute scaled damped residual for Chebyshev.
Definition Ifpack2_Details_ScaledDampedResidual_decl.hpp:45
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40