Ifpack2 Templated Preconditioning Package Version 1.0
Loading...
Searching...
No Matches
Ifpack2_Details_LinearSolver_def.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
13
14#ifndef IFPACK2_DETAILS_LINEARSOLVER_DEF_HPP
15#define IFPACK2_DETAILS_LINEARSOLVER_DEF_HPP
16
18#include "Tpetra_MultiVector.hpp"
20
21// Ifpack2: key is for Ifpack2's factory to have subordinate
22// factories. That way, each package still has one factory, but we
23// don't have to worry about intrapackage circular dependencies (e.g.,
24// relating to AdditiveSchwarz). There are two approaches:
25//
26// 1. Reuse existing Ifpack2::Details::OneLevelFactory
27// 2. Have each Ifpack2 solver register itself with Ifpack2's factory
28
29namespace Ifpack2 {
30namespace Details {
31
32template <class SC, class LO, class GO, class NT>
34 LinearSolver(const Teuchos::RCP<prec_type>& solver, const std::string& solverName)
35 : solver_(solver)
36 , solverName_(solverName) {
37 using Teuchos::RCP;
38 using Teuchos::rcp_dynamic_cast;
39 const char prefix[] = "Ifpack2::Details::LinearSolver: ";
40 TEUCHOS_TEST_FOR_EXCEPTION(solver.is_null(), std::invalid_argument,
41 prefix << "Input solver is NULL.");
42
43 typedef Tpetra::RowMatrix<SC, LO, GO, NT> row_matrix_type;
44 typedef ::Ifpack2::Details::CanChangeMatrix<row_matrix_type> mixin_type;
46 TEUCHOS_TEST_FOR_EXCEPTION(innerSolver.is_null(), std::invalid_argument, prefix << "The input "
47 "solver does not implement the setMatrix() feature. Only Ifpack2 solvers "
48 "that inherit from Ifpack2::Details::CanChangeMatrix implement this feature.");
49}
50
51template <class SC, class LO, class GO, class NT>
53 setMatrix(const Teuchos::RCP<const OP>& A) {
54 using Teuchos::RCP;
55 using Teuchos::rcp_dynamic_cast;
56 typedef Tpetra::RowMatrix<SC, LO, GO, NT> row_matrix_type;
57 typedef ::Ifpack2::Details::CanChangeMatrix<row_matrix_type> mixin_type;
58 const char prefix[] = "Ifpack2::Details::LinearSolver::setMatrix: ";
59
60 // It's OK for the input matrix to be null. Ifpack2 solvers may
61 // interpret this as a hint to clear out their state. It's also a
62 // way to keep the preconditioner around, but disassociate it from a
63 // particular matrix. (The code that uses the preconditioner might
64 // not want to or be able to recreate it.)
66 if (!A.is_null()) {
68 TEUCHOS_TEST_FOR_EXCEPTION(A_row.is_null(), std::invalid_argument, prefix << "The input matrix A, "
69 "if not null, must be a Tpetra::RowMatrix.");
70 }
71 TEUCHOS_TEST_FOR_EXCEPTION(solver_.is_null(), std::logic_error, prefix << "Solver is NULL. "
72 "This should never happen! Please report this bug to the Ifpack2 "
73 "developers.");
74
76 TEUCHOS_TEST_FOR_EXCEPTION(innerSolver.is_null(), std::logic_error, prefix << "The solver does not "
77 "implement the setMatrix() feature. Only input preconditioners that "
78 "inherit from Ifpack2::Details::CanChangeMatrix implement this. We should"
79 " never get here! Please report this bug to the Ifpack2 developers.");
81
82 A_ = A; // keep a pointer to A, so that getMatrix() works
83}
84
85template <class SC, class LO, class GO, class NT>
87 setCoord(const Teuchos::RCP<const coord_type>& C) {
88 using Teuchos::RCP;
89 using Teuchos::rcp_dynamic_cast;
90 typedef Tpetra::RowMatrix<SC, LO, GO, NT> row_matrix_type;
91 const char prefix[] = "Ifpack2::Details::LinearSolver::setCoord: ";
92
93 TEUCHOS_TEST_FOR_EXCEPTION(solver_.is_null(), std::logic_error, prefix << "Solver is NULL. "
94 "This should never happen! Please report this bug to the Ifpack2 "
95 "developers.");
96
98 TEUCHOS_TEST_FOR_EXCEPTION(innerSolver.is_null(), std::logic_error, prefix << "The solver does not "
99 "implement the setCoord() feature. We should never get here! "
100 "Please report this bug to the Ifpack2 developers.");
102
103 C_ = C;
104}
105
106template <class SC, class LO, class GO, class NT>
107Teuchos::RCP<const typename LinearSolver<SC, LO, GO, NT>::OP>
109 getMatrix() const {
110 return A_; // may be null
111}
112
113template <class SC, class LO, class GO, class NT>
114Teuchos::RCP<const typename LinearSolver<SC, LO, GO, NT>::coord_type>
116 getCoord() const {
117 return C_;
118}
119
120template <class SC, class LO, class GO, class NT>
122 solve(MV& X, const MV& B) {
123 const char prefix[] = "Ifpack2::Details::LinearSolver::solve: ";
124 TEUCHOS_TEST_FOR_EXCEPTION(solver_.is_null(), std::logic_error, prefix << "The solver is NULL! "
125 "This should never happen. Please report this bug to the Ifpack2 "
126 "developers.");
127 TEUCHOS_TEST_FOR_EXCEPTION(A_.is_null(), std::runtime_error, prefix << "The matrix has not been "
128 "set yet. You must call setMatrix() with a nonnull matrix before you "
129 "may call this method.");
130 solver_->apply(B, X);
131}
132
133template <class SC, class LO, class GO, class NT>
135 setParameters(const Teuchos::RCP<Teuchos::ParameterList>& params) {
136 solver_->setParameters(*params);
137}
138
139template <class SC, class LO, class GO, class NT>
141 symbolic() {
142 const char prefix[] = "Ifpack2::Details::LinearSolver::symbolic: ";
143 TEUCHOS_TEST_FOR_EXCEPTION(solver_.is_null(), std::logic_error, prefix << "The solver is NULL! "
144 "This should never happen. Please report this bug to the Ifpack2 "
145 "developers.");
146 TEUCHOS_TEST_FOR_EXCEPTION(A_.is_null(), std::runtime_error, prefix << "The matrix has not been "
147 "set yet. You must call setMatrix() with a nonnull matrix before you "
148 "may call this method.");
149 solver_->initialize();
150}
151
152template <class SC, class LO, class GO, class NT>
154 numeric() {
155 const char prefix[] = "Ifpack2::Details::LinearSolver::numeric: ";
156 TEUCHOS_TEST_FOR_EXCEPTION(solver_.is_null(), std::logic_error, prefix << "The solver is NULL! "
157 "This should never happen. Please report this bug to the Ifpack2 "
158 "developers.");
159 TEUCHOS_TEST_FOR_EXCEPTION(A_.is_null(), std::runtime_error, prefix << "The matrix has not been "
160 "set yet. You must call setMatrix() with a nonnull matrix before you "
161 "may call this method.");
162 solver_->compute();
163}
164
165template <class SC, class LO, class GO, class NT>
166std::string
168 description() const {
169 const char prefix[] = "Ifpack2::Details::LinearSolver::description: ";
170 TEUCHOS_TEST_FOR_EXCEPTION(solver_.is_null(), std::logic_error, prefix << "The solver is NULL! "
171 "This should never happen. Please report this bug to the Ifpack2 "
172 "developers.");
173 return solver_->description();
174}
175
176template <class SC, class LO, class GO, class NT>
178 describe(Teuchos::FancyOStream& out,
179 const Teuchos::EVerbosityLevel verbLevel) const {
180 const char prefix[] = "Ifpack2::Details::LinearSolver::describe: ";
181 TEUCHOS_TEST_FOR_EXCEPTION(solver_.is_null(), std::logic_error, prefix << "The solver is NULL! "
182 "This should never happen. Please report this bug to the Ifpack2 "
183 "developers.");
184 solver_->describe(out, verbLevel);
185}
186
187} // namespace Details
188} // namespace Ifpack2
189
190// Explicit template instantiation macro for LinearSolver. This is
191// generally not for users! It is used by automatically generated
192// code, and perhaps by expert Trilinos developers.
193#define IFPACK2_DETAILS_LINEARSOLVER_INSTANT(SC, LO, GO, NT) \
194 template class Ifpack2::Details::LinearSolver<SC, LO, GO, NT>;
195
196#endif // IFPACK2_DETAILS_LINEARSOLVER_DEF_HPP
Declaration of interface for preconditioners that can change their matrix after construction.
Declaration of MDF interface.
Ifpack2's implementation of Trilinos::Details::LinearSolver interface.
Definition Ifpack2_Details_LinearSolver_decl.hpp:75
Teuchos::RCP< const coord_type > getCoord() const
Get the coordinates associated with the solver's matrix rows.
Definition Ifpack2_Details_LinearSolver_def.hpp:116
void setCoord(const Teuchos::RCP< const coord_type > &C)
Set the matrix rows' coordinates.
Definition Ifpack2_Details_LinearSolver_def.hpp:87
void numeric()
Precompute for matrix values' changes.
Definition Ifpack2_Details_LinearSolver_def.hpp:154
void symbolic()
Precompute for matrix structure changes.
Definition Ifpack2_Details_LinearSolver_def.hpp:141
void solve(MV &X, const MV &B)
Solve the linear system AX=B for X.
Definition Ifpack2_Details_LinearSolver_def.hpp:122
void setMatrix(const Teuchos::RCP< const OP > &A)
Set the solver's matrix.
Definition Ifpack2_Details_LinearSolver_def.hpp:53
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params)
Set the solver's parameters.
Definition Ifpack2_Details_LinearSolver_def.hpp:135
std::string description() const
Implementation of Teuchos::Describable::description.
Definition Ifpack2_Details_LinearSolver_def.hpp:168
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Implementation of Teuchos::Describable::describe.
Definition Ifpack2_Details_LinearSolver_def.hpp:178
Teuchos::RCP< const OP > getMatrix() const
Get the solver's matrix.
Definition Ifpack2_Details_LinearSolver_def.hpp:109
LinearSolver(const Teuchos::RCP< prec_type > &solver, const std::string &solverName)
Constructor.
Definition Ifpack2_Details_LinearSolver_def.hpp:34
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40