OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
BelosBottomSolver.h
Go to the documentation of this file.
1 #ifndef BELOS_SOLVER_H
2 #define BELOS_SOLVER_H
3 
4 #include "BottomSolver.h"
5 
6 #include <BelosLinearProblem.hpp>
7 #include <BelosTpetraAdapter.hpp>
8 #include <BelosSolverFactory.hpp>
9 // #include <BelosSolverFactory_Tpetra.hpp>
10 
11 #include "AmrPreconditioner.h"
12 
13 #include <string>
14 
16 template <class Level>
17 class BelosBottomSolver : public BottomSolver<Teuchos::RCP<amr::matrix_t>,
18  Teuchos::RCP<amr::multivector_t>,
19  Level>
20 {
21 public:
30 
31  typedef Belos::SolverManager<scalar_t, mv_t, op_t> solver_t;
32  typedef Belos::LinearProblem<scalar_t, mv_t, op_t> problem_t;
33 
35 
36 public:
37 
42  BelosBottomSolver(std::string solvertype = "Pseudoblock CG",
43  const std::shared_ptr<prec_t>& prec_p = nullptr);
44 
45  void solve(const Teuchos::RCP<mv_t>& x,
46  const Teuchos::RCP<mv_t>& b);
47 
48  void setOperator(const Teuchos::RCP<matrix_t>& A,
49  Level* level_p = nullptr);
50 
51  std::size_t getNumIters();
52 
53 private:
58  void initSolver_m(std::string solvertype);
59 
60 private:
61  Teuchos::RCP<problem_t> problem_mp;
62  Teuchos::RCP<Teuchos::ParameterList> params_mp;
63  Teuchos::RCP<solver_t> solver_mp;
64  std::shared_ptr<prec_t> prec_mp;
65  Teuchos::RCP<matrix_t> A_mp;
66 
68 
70  int maxiter_m;
71 };
72 
73 #include "BelosBottomSolver.hpp"
74 
75 #endif
Teuchos::RCP< matrix_t > A_mp
copy of matrix (has to be positive definite)
std::shared_ptr< prec_t > prec_mp
preconditioner
Belos::LinearProblem< scalar_t, mv_t, op_t > problem_t
amr::operator_t op_t
Bottom solver preconditioners.
void setOperator(const Teuchos::RCP< matrix_t > &A, Level *level_p=nullptr)
int local_ordinal_t
Teuchos::RCP< solver_t > solver_mp
solver instance
std::size_t getNumIters()
KokkosClassic::DefaultNode::DefaultNodeType node_t
amr::local_ordinal_t lo_t
amr::multivector_t mv_t
amr::vector_t vector_t
Abstract base class for all base level solvers.
Definition: BottomSolver.h:8
amr::global_ordinal_t go_t
void solve(const Teuchos::RCP< mv_t > &x, const Teuchos::RCP< mv_t > &b)
Teuchos::RCP< problem_t > problem_mp
represents linear problem
amr::matrix_t matrix_t
void initSolver_m(std::string solvertype)
BelosBottomSolver(std::string solvertype="Pseudoblock CG", const std::shared_ptr< prec_t > &prec_p=nullptr)
AmrPreconditioner< matrix_t, Level > prec_t
Belos::SolverManager< scalar_t, mv_t, op_t > solver_t
amr::scalar_t scalar_t
Teuchos::RCP< Teuchos::ParameterList > params_mp
parameter list of solver
Interface to Belos solvers of the Trilinos package.
Tpetra::CrsMatrix< scalar_t, local_ordinal_t, global_ordinal_t, node_t > matrix_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
Tpetra::Vector< scalar_t, local_ordinal_t, global_ordinal_t, node_t > vector_t
int maxiter_m
allowed number of steps for iterative solvers
double scalar_t
scalar_t reltol_m
relative tolerance