48 #ifndef MUELU_DETAILS_LINEARSOLVERFACTORY_DEF_HPP 
   49 #define MUELU_DETAILS_LINEARSOLVERFACTORY_DEF_HPP 
   51 #include "MueLu_config.hpp" 
   54 #include <type_traits> 
   56 #ifdef HAVE_MUELU_EPETRA 
   57 #  include "Epetra_CrsMatrix.h" 
   59 #endif // HAVE_MUELU_EPETRA 
   62 #ifdef HAVE_MUELU_TPETRA 
   63 #  include "Tpetra_Operator.hpp" 
   65 #endif // HAVE_MUELU_TPETRA 
   70 template<
class MV, 
class OP, 
class NormType>
 
   96   void solve (MV& X, 
const MV& B);
 
  124 #if defined(HAVE_MUELU_SERIAL) and defined(HAVE_MUELU_EPETRA) 
  136     changedParams_(
false)
 
  148     const char prefix[] = 
"MueLu::Details::LinearSolver::setMatrix: ";
 
  152       if(solver_ != Teuchos::null)
 
  157         (
A_.
is_null(), std::runtime_error, prefix << 
"MueLu requires " 
  158          "an Epetra_CrsMatrix, but the matrix you provided is of a " 
  159          "different type.  Please provide an Epetra_CrsMatrix instead.");
 
  172     const char prefix[] = 
"MueLu::Details::LinearSolver::solve: ";
 
  174       (solver_.is_null (), std::runtime_error, prefix << 
"The solver does not " 
  175        "exist yet.  You must call numeric() before you may call this method.");
 
  177       (changedA_, std::runtime_error, prefix << 
"The matrix A has been reset " 
  178        "since the last call to numeric().  Please call numeric() again.");
 
  180       (changedParams_, std::runtime_error, prefix << 
"The parameters have been reset " 
  181        "since the last call to numeric().  Please call numeric() again.");
 
  183     int err = solver_->ApplyInverse(B, X);
 
  186       (err != 0, std::runtime_error, prefix << 
"EpetraOperator::ApplyInverse returned " 
  187        "nonzero error code " << err);
 
  193     if(solver_ != Teuchos::null && params != 
params_)
 
  194       changedParams_ = 
true;
 
  207     const char prefix[] = 
"MueLu::Details::LinearSolver::numeric: ";
 
  210     if(solver_ == Teuchos::null || changedA_ || changedParams_)
 
  213       changedParams_ = 
false;
 
  216         (
A_ == Teuchos::null, std::runtime_error, prefix << 
"The matrix has not been " 
  217          "set yet.  You must call setMatrix() with a nonnull matrix before you may " 
  218          "call this method.");
 
  232     if (solver_.is_null()) {
 
  233       return "\"MueLu::Details::LinearSolver\": {MV: Epetra_MultiVector, OP: Epetra_Operator, NormType: double}";
 
  236       return solver_->GetHierarchy()->description ();
 
  247     if (solver_.is_null()) {
 
  250         out << 
"\"MueLu::Details::LinearSolver\":" << endl;
 
  252         out << 
"MV: Epetra_MultiVector" << endl
 
  253             << 
"OP: Epetra_Operator" << endl
 
  254             << 
"NormType: double" << endl;
 
  258       solver_->GetHierarchy()->describe (out, verbLevel);
 
  269 #endif // HAVE_MUELU_EPETRA 
  271 #ifdef HAVE_MUELU_TPETRA 
  272 template<
class Scalar, 
class LO, 
class GO, 
class Node>
 
  274                    Tpetra::Operator<Scalar,LO,GO,Node>,
 
  277                                            Tpetra::Operator<Scalar,LO,GO,Node>, 
 
  278                                            typename Teuchos::ScalarTraits<Scalar>::magnitudeType>,
 
  287     changedParams_(false)
 
  301       if(solver_ != Teuchos::null)
 
  314   void solve (Tpetra::MultiVector<Scalar,LO,GO,Node>& X, 
const Tpetra::MultiVector<Scalar,LO,GO,Node>& B)
 
  317     const char prefix[] = 
"MueLu::Details::LinearSolver::solve: ";
 
  319       (solver_.is_null (), std::runtime_error, prefix << 
"The solver does not " 
  320        "exist yet.  You must call numeric() before you may call this method.");
 
  322       (changedA_, std::runtime_error, prefix << 
"The matrix A has been reset " 
  323        "since the last call to numeric().  Please call numeric() again.");
 
  324     TEUCHOS_TEST_FOR_EXCEPTION
 
  325       (changedParams_, std::runtime_error, prefix << 
"The parameters have been reset " 
  326        "since the last call to numeric().  Please call numeric() again.");
 
  328     solver_->apply(B, X);
 
  334     if(solver_ != Teuchos::null && params != 
params_)
 
  335       changedParams_ = 
true;
 
  348     const char prefix[] = 
"MueLu::Details::LinearSolver::numeric: ";
 
  351     if(solver_ == Teuchos::null || changedParams_)
 
  354         (
A_ == Teuchos::null, std::runtime_error, prefix << 
"The matrix has not been " 
  355          "set yet.  You must call setMatrix() with a nonnull matrix before you may " 
  356          "call this method.");
 
  368         (
A_ == Teuchos::null, std::runtime_error, prefix << 
"The matrix has not been " 
  369          "set yet.  You must call setMatrix() with a nonnull matrix before you may " 
  370          "call this method.");
 
  375       helperMat = rcp_dynamic_cast<
const Tpetra::CrsMatrix<Scalar,LO,GO,Node> >(
A_);
 
  377         (helperMat.
is_null(), std::runtime_error, prefix << 
"MueLu requires " 
  378          "a Tpetra::CrsMatrix, but the matrix you provided is of a " 
  379          "different type.  Please provide a Tpetra::CrsMatrix instead.");
 
  384     changedParams_ = 
false;
 
  391     if (solver_.is_null()) {
 
  392       std::ostringstream os;
 
  393       os << 
"\"MueLu::Details::LinearSolver\": {" 
  394          << 
"MV: " << TypeNameTraits<Tpetra::MultiVector<Scalar,LO,GO,Node> >::name()
 
  395          << 
"OP: " << TypeNameTraits<Tpetra::Operator<Scalar,LO,GO,Node> >::name()
 
  396          << 
"NormType: " << TypeNameTraits<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>::name()
 
  401       return solver_->GetHierarchy()->description ();
 
  413     if (solver_.is_null()) {
 
  416         out << 
"\"MueLu::Details::LinearSolver\":" << endl;
 
  418         out << 
"MV: " << TypeNameTraits<Tpetra::MultiVector<Scalar,LO,GO,Node> >::name() << endl
 
  419             << 
"OP: " << TypeNameTraits<Tpetra::Operator<Scalar,LO,GO,Node> >::name() << endl
 
  420             << 
"NormType: " << TypeNameTraits<typename Teuchos::ScalarTraits<Scalar>::magnitudeType>::name() << endl;
 
  424       solver_->GetHierarchy()->describe (out, verbLevel);
 
  435 #endif // HAVE_MUELU_TPETRA 
  437 template<
class MV, 
class OP, 
class NormType>
 
  446 template<
class MV, 
class OP, 
class NormType>
 
  451 #ifdef HAVE_TEUCHOSCORE_CXX11 
  452   typedef std::shared_ptr<MueLu::Details::LinearSolverFactory<MV, OP, NormType> > ptr_type;
 
  457 #endif // HAVE_TEUCHOSCORE_CXX11 
  460   Trilinos::Details::registerLinearSolverFactory<MV, OP, NormType> (
"MueLu", factory);
 
  473 #define MUELU_DETAILS_LINEARSOLVERFACTORY_INSTANT(SC, LO, GO, NT) \ 
  474   template class MueLu::Details::LinearSolverFactory<Tpetra::MultiVector<SC, LO, GO, NT>, \ 
  475                                                      Tpetra::Operator<SC, LO, GO, NT>, \ 
  476                                                      typename Tpetra::MultiVector<SC, LO, GO, NT>::mag_type>; 
  478 #endif // MUELU_DETAILS_LINEARSOLVERFACTORY_DEF_HPP 
Interface for a "factory" that creates MueLu solvers. 
 
Teuchos::RCP< const OP > getMatrix() const 
Get a pointer to this Solver's matrix. 
 
std::string description() const 
Implementation of Teuchos::Describable::description. 
 
void symbolic()
Set up any part of the solve that depends on the structure of the input matrix, but not its numerical...
 
Teuchos::RCP< const Tpetra::Operator< Scalar, LO, GO, Node > > getMatrix() const 
Get a pointer to this Solver's matrix. 
 
Teuchos::RCP< Teuchos::ParameterList > params_
 
Teuchos::RCP< const Tpetra::Operator< Scalar, LO, GO, Node > > A_
 
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set this solver's parameters. 
 
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
 
std::string description() const 
Implementation of Teuchos::Describable::description. 
 
virtual Teuchos::RCP< Trilinos::Details::LinearSolver< MV, OP, NormType > > getLinearSolver(const std::string &solverName)
Get an instance of a MueLu solver. 
 
Teuchos::RCP< TpetraOperator< Scalar, LO, GO, Node > > solver_
 
void numeric()
Set up any part of the solve that depends on both the structure and the numerical values of the input...
 
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
 
LinearSolver()
Constructor. 
 
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
Set this solver's parameters. 
 
virtual ~LinearSolver()
Destructor (virtual for memory safety). 
 
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Implementation of Teuchos::Describable::describe. 
 
Various adapters that will create a MueLu preconditioner that is a Tpetra::Operator. 
 
static const EVerbosityLevel verbLevel_default
 
LinearSolver()
Constructor. 
 
Teuchos::RCP< const OP > A_
 
Teuchos::RCP< MueLu::EpetraOperator > CreateEpetraPreconditioner(const Teuchos::RCP< Epetra_CrsMatrix > &inA, Teuchos::ParameterList ¶mListIn)
Helper function to create a MueLu preconditioner that can be used by Epetra.Given a EpetraCrs_Matrix...
 
void setMatrix(const Teuchos::RCP< const Tpetra::Operator< Scalar, LO, GO, Node > > &A)
Set the Solver's matrix. 
 
static void registerLinearSolverFactory()
Register this LinearSolverFactory with the central registry. 
 
Various adapters that will create a MueLu preconditioner that is an Epetra_Operator. 
 
virtual ~LinearSolver()
Destructor (virtual for memory safety). 
 
void symbolic()
Set up any part of the solve that depends on the structure of the input matrix, but not its numerical...
 
Teuchos::RCP< Teuchos::ParameterList > params_
 
void solve(MV &X, const MV &B)
Solve the linear system(s) AX=B. 
 
void numeric()
Set up any part of the solve that depends on both the structure and the numerical values of the input...
 
void ReuseTpetraPreconditioner(const Teuchos::RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &inA, MueLu::TpetraOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Op)
Helper function to reuse an existing MueLu preconditioner. 
 
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const 
Implementation of Teuchos::Describable::describe. 
 
Teuchos::RCP< MueLu::TpetraOperator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > CreateTpetraPreconditioner(const Teuchos::RCP< Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &inA, Teuchos::ParameterList &inParamList)
Helper function to create a MueLu or AMGX preconditioner that can be used by Tpetra.Given a Tpetra::Operator, this function returns a constructed MueLu preconditioner. 
 
void setMatrix(const Teuchos::RCP< const OP > &A)
Set the Solver's matrix. 
 
void solve(Tpetra::MultiVector< Scalar, LO, GO, Node > &X, const Tpetra::MultiVector< Scalar, LO, GO, Node > &B)
Solve the linear system(s) AX=B.