MueLu Version of the Day
Loading...
Searching...
No Matches
MueLu_Amesos2Smoother_decl.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// MueLu: A package for multigrid based preconditioning
4//
5// Copyright 2012 NTESS and the MueLu contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef MUELU_AMESOS2SMOOTHER_DECL_HPP
11#define MUELU_AMESOS2SMOOTHER_DECL_HPP
12
13#include "MueLu_ConfigDefs.hpp"
14
15#include <Teuchos_ParameterList.hpp>
16#include <Teuchos_SerialDenseMatrix.hpp>
17#include <Teuchos_LAPACK.hpp>
18
20
21#include "MueLu_SmootherPrototype.hpp"
24
25namespace Amesos2 {
26template <class OP, class MV>
27class Solver;
28}
29
30namespace MueLu {
31
32template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
34 public:
35 Projection(RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>>& Nullspace);
36
37 void projectOut(Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>& X);
38
39 Teuchos::RCP<Xpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node>> Nullspace_;
40
41 private:
42 Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> localMap_;
43};
44
54template <class Scalar = SmootherPrototype<>::scalar_type,
58class Amesos2Smoother : public SmootherPrototype<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
59#undef MUELU_AMESOS2SMOOTHER_SHORT
61
62 public:
64
65
70 Amesos2Smoother(const std::string& type = "", const Teuchos::ParameterList& paramList = Teuchos::ParameterList());
71
73 virtual ~Amesos2Smoother();
74
75 RCP<const ParameterList> GetValidParameterList() const;
76
78
80
81
82 void DeclareInput(Level& currentLevel) const;
83
85
87
88
93 void Setup(Level& currentLevel);
94
101 void Apply(MultiVector& X, const MultiVector& B, bool InitialGuessIsZero = false) const;
103
104 RCP<SmootherPrototype> Copy() const;
105
107
108
110 std::string description() const;
111
113 // using MueLu::Describable::describe; // overloading, not hiding
114 void print(Teuchos::FancyOStream& out, const VerbLevel verbLevel = Default) const;
115
117 size_t getNodeSmootherComplexity() const;
118
120
121 private:
122 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> Tpetra_CrsMatrix;
123 typedef Tpetra::MultiVector<Scalar, LocalOrdinal, GlobalOrdinal, Node> Tpetra_MultiVector;
124
126 std::string type_;
127
129 RCP<Amesos2::Solver<Tpetra_CrsMatrix, Tpetra_MultiVector>> prec_;
130
132 RCP<MultiVector> X_, B_;
133
134 RCP<Projection<Scalar, LocalOrdinal, GlobalOrdinal, Node>> projection_;
135
136}; // class Amesos2Smoother
137
138} // namespace MueLu
139
140#define MUELU_AMESOS2SMOOTHER_SHORT
141#endif // MUELU_AMESOS2SMOOTHER_DECL_HPP
MueLu::DefaultLocalOrdinal LocalOrdinal
MueLu::DefaultGlobalOrdinal GlobalOrdinal
MueLu::DefaultNode Node
Class that encapsulates Amesos2 direct solvers.
void DeclareInput(Level &currentLevel) const
Input.
size_t getNodeSmootherComplexity() const
Get a rough estimate of cost per iteration.
std::string type_
amesos2-specific key phrase that denote smoother type
RCP< SmootherPrototype > Copy() const
std::string description() const
Return a simple one-line description of this object.
void print(Teuchos::FancyOStream &out, const VerbLevel verbLevel=Default) const
Print the object with some verbosity level to an FancyOStream object.
Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Tpetra_CrsMatrix
RCP< Projection< Scalar, LocalOrdinal, GlobalOrdinal, Node > > projection_
RCP< Amesos2::Solver< Tpetra_CrsMatrix, Tpetra_MultiVector > > prec_
pointer to Amesos2 solver object
Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > Tpetra_MultiVector
RCP< const ParameterList > GetValidParameterList() const
Return a const parameter list of valid parameters that setParameterList() will accept.
void Setup(Level &currentLevel)
Set up the direct solver. This creates the underlying Amesos2 solver object according to the paramete...
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.
Class that holds all level-specific information.
void projectOut(Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X)
Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > localMap_
Teuchos::RCP< Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Nullspace_
Base class for smoother prototypes.
Namespace for MueLu classes and methods.