14#ifndef IFPACK2_DETAILS_LINEARSOLVER_DEF_HPP
15#define IFPACK2_DETAILS_LINEARSOLVER_DEF_HPP
18#include "Tpetra_MultiVector.hpp"
32template <
class SC,
class LO,
class GO,
class NT>
38 using Teuchos::rcp_dynamic_cast;
39 const char prefix[] =
"Ifpack2::Details::LinearSolver: ";
41 prefix <<
"Input solver is NULL.");
43 typedef Tpetra::RowMatrix<SC, LO, GO, NT> row_matrix_type;
44 typedef ::Ifpack2::Details::CanChangeMatrix<row_matrix_type>
mixin_type;
47 "solver does not implement the setMatrix() feature. Only Ifpack2 solvers "
48 "that inherit from Ifpack2::Details::CanChangeMatrix implement this feature.");
51template <
class SC,
class LO,
class GO,
class NT>
53 setMatrix(
const Teuchos::RCP<const OP>& A) {
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: ";
69 "if not null, must be a Tpetra::RowMatrix.");
72 "This should never happen! Please report this bug to the Ifpack2 "
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.");
85template <
class SC,
class LO,
class GO,
class NT>
87 setCoord(
const Teuchos::RCP<const coord_type>& C) {
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: ";
94 "This should never happen! Please report this bug to the Ifpack2 "
99 "implement the setCoord() feature. We should never get here! "
100 "Please report this bug to the Ifpack2 developers.");
106template <
class SC,
class LO,
class GO,
class NT>
107Teuchos::RCP<const typename LinearSolver<SC, LO, GO, NT>::OP>
113template <
class SC,
class LO,
class GO,
class NT>
114Teuchos::RCP<const typename LinearSolver<SC, LO, GO, NT>::coord_type>
120template <
class SC,
class LO,
class GO,
class NT>
122 solve(MV& X,
const MV& B) {
123 const char prefix[] =
"Ifpack2::Details::LinearSolver::solve: ";
125 "This should never happen. Please report this bug to the Ifpack2 "
128 "set yet. You must call setMatrix() with a nonnull matrix before you "
129 "may call this method.");
130 solver_->apply(B, X);
133template <
class SC,
class LO,
class GO,
class NT>
136 solver_->setParameters(*
params);
139template <
class SC,
class LO,
class GO,
class NT>
142 const char prefix[] =
"Ifpack2::Details::LinearSolver::symbolic: ";
144 "This should never happen. Please report this bug to the Ifpack2 "
147 "set yet. You must call setMatrix() with a nonnull matrix before you "
148 "may call this method.");
149 solver_->initialize();
152template <
class SC,
class LO,
class GO,
class NT>
155 const char prefix[] =
"Ifpack2::Details::LinearSolver::numeric: ";
157 "This should never happen. Please report this bug to the Ifpack2 "
160 "set yet. You must call setMatrix() with a nonnull matrix before you "
161 "may call this method.");
165template <
class SC,
class LO,
class GO,
class NT>
169 const char prefix[] =
"Ifpack2::Details::LinearSolver::description: ";
171 "This should never happen. Please report this bug to the Ifpack2 "
173 return solver_->description();
176template <
class SC,
class LO,
class GO,
class NT>
179 const Teuchos::EVerbosityLevel
verbLevel)
const {
180 const char prefix[] =
"Ifpack2::Details::LinearSolver::describe: ";
182 "This should never happen. Please report this bug to the Ifpack2 "
193#define IFPACK2_DETAILS_LINEARSOLVER_INSTANT(SC, LO, GO, NT) \
194 template class Ifpack2::Details::LinearSolver<SC, LO, GO, NT>;
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 > ¶ms)
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