14 #if defined(HAVE_MUELU_EPETRA) && defined(HAVE_MUELU_AMESOS)
16 #include <Epetra_LinearProblem.h>
18 #include <Amesos_config.h>
25 #include "MueLu_Utilities.hpp"
38 std::transform(
type_.begin(), ++
type_.begin(),
type_.begin(), ::toupper);
42 if (
type_ ==
"Amesos_umfpack")
type_ =
"Umfpack";
43 if (
type_ ==
"Superlu_dist")
type_ =
"Superludist";
49 std::string oldtype =
type_;
51 #if defined(HAVE_AMESOS_SUPERLU)
53 #elif defined(HAVE_AMESOS_KLU)
55 #elif defined(HAVE_AMESOS_SUPERLUDIST)
56 type_ =
"Superludist";
57 #elif defined(HAVE_AMESOS_UMFPACK)
60 this->
declareConstructionOutcome(
true,
"Amesos has been compiled without SuperLU_DIST, SuperLU, Umfpack or Klu. By default, MueLu tries" +
61 "to use one of these libraries. Amesos must be compiled with one of these solvers, " +
62 "or a valid Amesos solver has to be specified explicitly.");
66 this->
GetOStream(
Warnings0) <<
"MueLu::AmesosSmoother: \"" << oldtype <<
"\" is not available. Using \"" <<
type_ <<
"\" instead" << std::endl;
75 this->Input(currentLevel,
"A");
83 this->GetOStream(
Warnings0) <<
"MueLu::AmesosSmoother::Setup(): Setup() has already been called" << std::endl;
85 A_ = Factory::Get<RCP<Matrix> >(currentLevel,
"A");
89 linearProblem_->SetOperator(epA.get());
92 prec_ =
rcp(factory.
Create(type_, *linearProblem_));
97 if (A_->getRowMap()->isDistributed() ==
true && A_->getRowMap()->isContiguous() ==
false)
98 const_cast<ParameterList&>(this->GetParameterList()).set(
"Reindex",
true);
103 prec_->SetParameters(*precList);
105 const_cast<ParameterList&
>(paramList).setParameters(*precList);
107 int r = prec_->NumericFactorization();
113 template <
class Node>
123 linearProblem_->SetLHS(&epX);
124 linearProblem_->SetRHS(&nonconstB);
129 linearProblem_->SetLHS(0);
130 linearProblem_->SetRHS(0);
133 template <
class Node>
138 template <
class Node>
140 std::ostringstream out;
142 out <<
"{type = " << type_ <<
"}";
147 template <
class Node>
152 out0 <<
"Prec. type: " << type_ << std::endl;
155 out0 <<
"Parameter list: " << std::endl;
157 out << this->GetParameterList();
161 if (prec_ != Teuchos::null) {
162 prec_->PrintStatus();
163 prec_->PrintTiming();
166 if (verbLevel &
Debug) {
169 <<
"RCP<A_>: " << A_ << std::endl
170 <<
"RCP<linearProblem__>: " << linearProblem_ << std::endl
171 <<
"RCP<prec_>: " << prec_ << std::endl;
175 template <
class Node>
186 #if defined(HAVE_MUELU_EPETRA)
190 #endif // HAVE_MUELU_EPETRA && HAVE_MUELU_AMESOS
Important warning messages (one line)
Teuchos::FancyOStream & GetOStream(MsgType type, int thisProcRankOnly=0) const
Get an output stream for outputting the input message type.
void Apply(MultiVector &X, const MultiVector &B, bool=false) const
Apply the direct solver.
Print external lib objects.
static RCP< Epetra_CrsMatrix > Op2NonConstEpetraCrs(RCP< Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >> Op)
Timer to be used in factories. Similar to Monitor but with additional timers.
AmesosSmoother(const std::string &type="", const Teuchos::ParameterList ¶mList=Teuchos::ParameterList())
Constructor.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::string description() const
Return a simple one-line description of this object.
Print additional debugging information.
std::string tolower(const std::string &str)
void Setup(Level ¤tLevel)
Set up the direct solver. This creates the underlying Amesos solver object according to the parameter...
virtual void SetParameterList(const Teuchos::ParameterList ¶mList)
Set parameters from a parameter list and return with default values.
std::string type_
amesos-specific key phrase that denote smoother type
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
static RCP< Epetra_MultiVector > MV2NonConstEpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> vec)
void declareConstructionOutcome(bool fail, std::string msg)
Class that holds all level-specific information.
RCP< SmootherPrototype > Copy() const
bool IsSetup() const
Get the state of a smoother prototype.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
Amesos_BaseSolver * Create(const char *ClassType, const Epetra_LinearProblem &LinearProblem)
static RCP< const Epetra_MultiVector > MV2EpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node >> const vec)
Helper utility to pull out the underlying Epetra objects from an Xpetra object.
Print class parameters (more parameters, more verbose)
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
Exception throws to report errors in the internal logical of the program.
Description of what is happening (more verbose)
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
virtual std::string description() const
Return a simple one-line description of this object.
std::string toString(const T &t)
bool Query(const char *ClassType)
void DeclareInput(Level ¤tLevel) const
Input.
Class that encapsulates Amesos direct solvers.