46 #ifndef MUELU_PROJECTORSMOOTHER_DEF_HPP
47 #define MUELU_PROJECTORSMOOTHER_DEF_HPP
49 #include <Xpetra_Matrix.hpp>
50 #include <Xpetra_MultiVectorFactory.hpp>
56 #include "MueLu_Utilities.hpp"
61 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
63 : coarseSolver_(coarseSolver)
66 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
69 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
71 Factory::Input(currentLevel,
"A");
72 Factory::Input(currentLevel,
"Nullspace");
74 coarseSolver_->DeclareInput(currentLevel);
77 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
81 coarseSolver_->Setup(currentLevel);
84 this->GetOStream(
Warnings0) <<
"MueLu::ProjectorSmoother::Setup(): Setup() has already been called" << std::endl;
86 RCP<Matrix> A = Factory::Get< RCP<Matrix> > (currentLevel,
"A");
87 RCP<MultiVector> B = Factory::Get< RCP<MultiVector> >(currentLevel,
"Nullspace");
89 int m = B->getNumVectors();
99 for (
int i = 0; i < m; i++) {
100 Scalar rayleigh = br[i] / bb[i];
105 this->GetOStream(
Runtime0) <<
"Coarse level orth indices: " << selectedIndices << std::endl;
107 #if defined(HAVE_MUELU_TPETRA) && defined(HAVE_XPETRA_TPETRA)
108 #ifdef HAVE_MUELU_TPETRA_INST_INT_INT
114 for (
int i = 0; i < selectedIndices.
size(); i++) {
119 for (
int k = 0; k < i; k++) {
130 Borth_ =
rcp(static_cast<MultiVector*>(
new TpetraMultiVector(Borth)));
139 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
141 coarseSolver_->Apply(X, B, InitialGuessIsZero);
143 int m = Borth_->getNumVectors();
144 int n = X.getNumVectors();
147 for (
int i = 0; i < n; i++) {
151 for (
int k = 0; k < m; k++) {
160 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
165 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
167 std::ostringstream out;
172 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
180 #endif // MUELU_PROJECTORSMOOTHER_DEF_HPP
Important warning messages (one line)
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
std::string description() const
Return a simple one-line description of this object.
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > MV2TpetraMV(RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > const vec)
Helper utility to pull out the underlying Tpetra objects from an Xpetra object.
Timer to be used in factories. Similar to Monitor but with additional timers.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
One-liner description of what is happening.
virtual ~ProjectorSmoother()
Destructor.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
MueLu::DefaultScalar Scalar
Class that holds all level-specific information.
RCP< SmootherPrototype > Copy() const
bool IsSetup() const
Get the state of a smoother prototype.
void Setup(Level ¤tLevel)
Set up the direct solver.
#define MUELU_DESCRIBE
Helper macro for implementing Describable::describe() for BaseClass objects.
This class enables the elimination of the nullspace component of the solution through the use of proj...
void push_back(const value_type &x)
void Apply(MultiVector &X, const MultiVector &B, bool InitialGuessIsZero=false) const
Apply the direct solver. Solves the linear system AX=B using the constructed solver.
Exception throws to report errors in the internal logical of the program.
ProjectorSmoother(RCP< SmootherPrototype > coarseSolver)
Constructor.
void DeclareInput(Level ¤tLevel) const
Input.
virtual std::string description() const
Return a simple one-line description of this object.