14#ifndef IFPACK2_DETAILS_LINEARSOLVER_DECL_HPP
15#define IFPACK2_DETAILS_LINEARSOLVER_DECL_HPP
17#include "Ifpack2_ConfigDefs.hpp"
18#include "Trilinos_Details_LinearSolver.hpp"
20#include "Teuchos_Describable.hpp"
71template <
class SC,
class LO,
class GO,
class NT>
72class LinearSolver :
public Trilinos::Details::LinearSolver<Tpetra::MultiVector<SC, LO, GO, NT>,
73 Tpetra::Operator<SC, LO, GO, NT>,
74 typename Tpetra::MultiVector<SC, LO, GO, NT>::mag_type>,
75 virtual public Teuchos::Describable {
78 typedef Tpetra::Operator<SC, LO, GO, NT> OP;
79 typedef Tpetra::MultiVector<SC, LO, GO, NT> MV;
80 typedef typename Teuchos::ScalarTraits<SC>::magnitudeType magnitude_type;
81 typedef Tpetra::MultiVector<magnitude_type, LO, GO, NT> coord_type;
103 void setMatrix(
const Teuchos::RCP<const OP>& A);
108 void setCoord(
const Teuchos::RCP<const coord_type>& C);
111 Teuchos::RCP<const OP>
getMatrix()
const;
114 Teuchos::RCP<const coord_type>
getCoord()
const;
117 void solve(MV& X,
const MV& B);
134 const Teuchos::EVerbosityLevel
verbLevel =
135 Teuchos::Describable::verbLevel_default)
const;
139 Teuchos::RCP<prec_type> solver_;
141 std::string solverName_;
143 Teuchos::RCP<const OP> A_;
146 Teuchos::RCP<const coord_type> C_;
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
virtual ~LinearSolver()
Destructor (virtual for memory safety).
Definition Ifpack2_Details_LinearSolver_decl.hpp:97
Interface for all Ifpack2 preconditioners.
Definition Ifpack2_Preconditioner.hpp:74
Ifpack2 implementation details.
Preconditioners and smoothers for Tpetra sparse matrices.
Definition Ifpack2_AdditiveSchwarz_decl.hpp:40