OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
MueLuBottomSolver.h
Go to the documentation of this file.
1 #ifndef MUELU_BOTTOM_SOLVER_H
2 #define MUELU_BOTTOM_SOLVER_H
3 
4 #include "BottomSolver.h"
5 
6 #include "Amr/AmrDefs.h"
7 
8 #include "Ippl.h"
9 
10 #include <MueLu.hpp>
11 #include <MueLu_Level.hpp>
12 #include <MueLu_Utilities.hpp>
13 
14 #include <MueLu_ParameterListInterpreter.hpp>
15 
16 template <class Level>
17 class MueLuBottomSolver : public BottomSolver<Teuchos::RCP<amr::matrix_t>,
18  Teuchos::RCP<amr::multivector_t>,
19  Level >
20 {
21 public:
30 
33 
34 // typedef amr::AmrGeometry_t AmrGeometry_t;
35 
36  typedef MueLu::Hierarchy<scalar_t, lo_t, go_t, node_t> hierarchy_t;
37  typedef MueLu::Level level_t;
38  typedef Xpetra::Matrix<scalar_t, lo_t, go_t, node_t> xmatrix_t;
39  typedef Xpetra::MultiVector<scalar_t, lo_t, go_t, node_t> xmv_t;
40  typedef MueLu::Utilities<scalar_t, lo_t, go_t, node_t> util_t;
41 
42  typedef MueLu::ParameterListInterpreter<scalar_t, lo_t, go_t, node_t> pListInterpreter_t;
43  typedef MueLu::HierarchyManager<scalar_t, lo_t, go_t, node_t> manager_t;
44 
45 public:
46 
47  MueLuBottomSolver(const bool& rebalance,
48  const std::string& reuse);
49 
50  void solve(const Teuchos::RCP<mv_t>& x,
51  const Teuchos::RCP<mv_t>& b);
52 
53  void setOperator(const Teuchos::RCP<matrix_t>& A,
54  Level* level_p = nullptr);
55 
56  std::size_t getNumIters();
57 
58  /*
59  * MueLu reuse option.
60  * Either none, RP, RAP or full
61  */
62  static std::string convertToMueLuReuseOption(const std::string& reuse);
63 
64 private:
65  void initMueLuList_m(const std::string& reuse);
66 
67 private:
68  Teuchos::RCP<hierarchy_t> hierarchy_mp;
69 
70  Teuchos::RCP<manager_t> factory_mp;
71 
72  Teuchos::RCP<xmatrix_t> A_mp;
73 
75 
76  Teuchos::ParameterList mueluList_m;
77 
78  bool rebalance_m;
79 
81 };
82 
83 #include "MueLuBottomSolver.hpp"
84 
85 #endif
void solve(const Teuchos::RCP< mv_t > &x, const Teuchos::RCP< mv_t > &b)
MueLu::HierarchyManager< scalar_t, lo_t, go_t, node_t > manager_t
void setOperator(const Teuchos::RCP< matrix_t > &A, Level *level_p=nullptr)
MueLu::Hierarchy< scalar_t, lo_t, go_t, node_t > hierarchy_t
amr::scalar_t scalar_t
MueLuBottomSolver(const bool &rebalance, const std::string &reuse)
amr::multivector_t mv_t
MueLu::Utilities< scalar_t, lo_t, go_t, node_t > util_t
int local_ordinal_t
amr::AmrIntVect_t AmrIntVect_t
Xpetra::MultiVector< scalar_t, lo_t, go_t, node_t > xmv_t
amrex::IntVect AmrIntVect_t
Definition: AmrDefs.h:28
amr::AmrBox_t AmrBox_t
KokkosClassic::DefaultNode::DefaultNodeType node_t
void initMueLuList_m(const std::string &reuse)
Teuchos::RCP< hierarchy_t > hierarchy_mp
manages the multigrid hierarchy
Abstract base class for all base level solvers.
Definition: BottomSolver.h:8
amr::vector_t vector_t
static std::string convertToMueLuReuseOption(const std::string &reuse)
std::size_t getNumIters()
amrex::Box AmrBox_t
Definition: AmrDefs.h:30
IpplTimings::TimerRef setupTimer_m
Teuchos::RCP< manager_t > factory_mp
sets up hierarchy
lo_t nSweeps_m
the number of multigrid iterations
Tpetra::CrsMatrix< scalar_t, local_ordinal_t, global_ordinal_t, node_t > matrix_t
amr::matrix_t matrix_t
MueLu::ParameterListInterpreter< scalar_t, lo_t, go_t, node_t > pListInterpreter_t
Tpetra::MultiVector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > multivector_t
Tpetra::Operator< scalar_t, local_ordinal_t, global_ordinal_t, node_t > operator_t
long global_ordinal_t
amr::global_ordinal_t go_t
Teuchos::ParameterList mueluList_m
MueLu::Level level_t
Timing::TimerRef TimerRef
Definition: IpplTimings.h:176
Tpetra::Vector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > vector_t
Teuchos::RCP< xmatrix_t > A_mp
MueLu requires Xpetra.
Xpetra::Matrix< scalar_t, lo_t, go_t, node_t > xmatrix_t
amr::operator_t op_t
amr::local_ordinal_t lo_t
double scalar_t
bool rebalance_m
use subcommunicators (less communication)