Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Thyra_BelosTpetrasSolverAdapter.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Stratimikos: Thyra-based strategies for linear solvers
4 //
5 // Copyright 2006 NTESS and the Stratimikos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef __Thyra_BelosTpetraSolverAdapter_hpp
11 #define __Thyra_BelosTpetraSolverAdapter_hpp
12 
13 #include "BelosConfigDefs.hpp"
14 #include "BelosLinearProblem.hpp"
15 // BelosThyraAdapter.hpp only includes this file if HAVE_BELOS_TSQR is
16 // defined. Thus, it's OK to include TSQR header files here.
17 
18 #include "Thyra_MultiVectorBase.hpp"
19 #include "Thyra_TpetraThyraWrappers.hpp"
20 
21 #include "Teuchos_ParameterList.hpp"
22 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
23 
24 #include "Tpetra_CrsMatrix.hpp"
25 #include "BelosTpetraAdapter.hpp"
26 #include "solvers/Belos_Tpetra_Krylov_parameters.hpp"
27 #include "solvers/Belos_Tpetra_Krylov.hpp"
28 #include "solvers/Belos_Tpetra_Gmres.hpp"
29 #include "solvers/Belos_Tpetra_GmresPipeline.hpp"
30 #include "solvers/Belos_Tpetra_GmresSingleReduce.hpp"
31 #include "solvers/Belos_Tpetra_GmresSstep.hpp"
32 
33 namespace Thyra {
34 
37  template<class SC, class MV, class OP>
38  class BelosTpetraKrylov : public Belos::SolverManager<SC, MV, OP> {
39  public:
40  using tpetra_base_solver_type = BelosTpetra::Impl::Krylov<SC>;
42 
44  BelosTpetraKrylov() = default;
45 
47  virtual Teuchos::RCP<Belos::SolverManager<SC, MV, OP> > clone () const override = 0;
48 
50  void setProblem( const Teuchos::RCP<Belos::LinearProblem<SC, MV, OP> > &problem ) override {
51  problem_ = problem;
52  }
53 
54  const Belos::LinearProblem<SC, MV, OP>& getProblem() const override {
55  return *problem_;
56  }
57 
59  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override {
60  params_ = params;
61  tpetra_solver->setParameters(*params);
62  }
63 
64  Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const override {
65  return params_;
66  }
67 
68  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const override {
69 
70  Teuchos::RCP<Teuchos::ParameterList> params = Teuchos::parameterList();
71 
72  bool defaultValues = true;
73  tpetra_solver->getParameters (*params, defaultValues);
74 
75  return params;
76  }
77 
79  int getNumIters() const override { return solver_output.numIters; }
80  bool isLOADetected() const override { return false; }
81  void reset( const Belos::ResetType type ) override {
82  if ((type & Belos::Problem) && !Teuchos::is_null(problem_))
83  problem_->setProblem();
84  }
85 
86  typename Teuchos::ScalarTraits<SC>::magnitudeType achievedTol() const override {
87  return tpetra_solver->achievedTol();
88  }
89 
91  Belos::ReturnType solve() override {
92  auto A = problem_->getOperator();
93  auto tpetraA = converter::getConstTpetraOperator(A);
94  tpetra_solver->setMatrix(tpetraA);
95 
96  auto lM = problem_->getLeftPrec();
97  if (lM != Teuchos::null) {
98  auto tpetraLM = converter::getConstTpetraOperator(lM);
99  tpetra_solver->setLeftPrec(tpetraLM);
100  }
101  auto rM = problem_->getRightPrec();
102  if (rM != Teuchos::null) {
103  auto tpetraRM = converter::getConstTpetraOperator(rM);
104  tpetra_solver->setRightPrec(tpetraRM);
105  }
106 
107  auto B = this->problem_->getRHS();
108  auto X = this->problem_->getLHS();
109  auto tpetraB = converter::getConstTpetraMultiVector(B);
110  auto tpetraX = converter::getTpetraMultiVector(X);
111 
112  solver_output = tpetra_solver->solve(*tpetraX, *tpetraB);
113 
114  // copy output
115  Belos::ReturnType belos_output = (solver_output.converged ? Belos::Converged : Belos::Unconverged);
116  return belos_output;
117  }
118 
119  protected:
120 
122  BelosTpetra::Impl::SolverOutput<SC> solver_output;
123 
126 
129  };
130 
131  //
132  // Gmres
133  template<class SC, class MV, class OP>
134  class BelosTpetraGmres : public BelosTpetraKrylov<SC, MV, OP> {
135  public:
136  using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
137  using tpetra_solver_type = BelosTpetra::Impl::Gmres<SC>;
138 
140  BelosTpetraGmres() :
141  krylov_base_solver_type ()
142  {
143  this->tpetra_solver = Teuchos::rcp(new tpetra_solver_type ());
144  }
145 
147  Teuchos::RCP<Belos::SolverManager<SC, MV, OP> > clone () const override {
148  return Teuchos::rcp(new BelosTpetraGmres<SC, MV, OP>);
149  }
150  };
151 
152  //
153  // Pipelined Gmres
154  template<class SC, class MV, class OP>
155  class BelosTpetraGmresPipeline : public BelosTpetraKrylov<SC, MV, OP> {
156  public:
157  using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
158  using tpetra_solver_type = BelosTpetra::Impl::GmresPipeline<SC>;
159 
161  BelosTpetraGmresPipeline() :
162  krylov_base_solver_type ()
163  {
164  this->tpetra_solver = Teuchos::rcp(new tpetra_solver_type ());
165  }
166 
168  Teuchos::RCP<Belos::SolverManager<SC, MV, OP> > clone () const override {
169  return Teuchos::rcp(new BelosTpetraGmresPipeline<SC, MV, OP>);
170  }
171  };
172 
173  //
174  // SingleReduce Gmres
175  template<class SC, class MV, class OP>
176  class BelosTpetraGmresSingleReduce : public BelosTpetraKrylov<SC, MV, OP> {
177  public:
178  using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
179  using tpetra_solver_type = BelosTpetra::Impl::GmresSingleReduce<SC>;
180 
182  BelosTpetraGmresSingleReduce() :
183  krylov_base_solver_type ()
184  {
185  this->tpetra_solver = Teuchos::rcp(new tpetra_solver_type ());
186  }
187 
189  Teuchos::RCP<Belos::SolverManager<SC, MV, OP> > clone () const override {
190  return Teuchos::rcp(new BelosTpetraGmresSingleReduce<SC, MV, OP>);
191  }
192  };
193 
194  //
195  // s-step Gmres
196  template<class SC, class MV, class OP>
197  class BelosTpetraGmresSstep : public BelosTpetraKrylov<SC, MV, OP> {
198  public:
199  using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
200  using tpetra_solver_type = BelosTpetra::Impl::GmresSstep<SC>;
201 
203  BelosTpetraGmresSstep() :
204  krylov_base_solver_type ()
205  {
206  this->tpetra_solver = Teuchos::rcp(new tpetra_solver_type ());
207  }
208 
210  Teuchos::RCP<Belos::SolverManager<SC, MV, OP> > clone () const override {
211  return Teuchos::rcp(new BelosTpetraGmresSstep<SC, MV, OP>);
212  }
213  };
214 
215 } // namespace Thyra
216 
217 #endif // __Thyra_TpetraSolverAdapter_hpp
218 
virtual void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)=0
void reset()
bool is_null(const std::shared_ptr< T > &p)
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraMultiVector(const RCP< MultiVectorBase< Scalar > > &mv)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
ParameterList & setParameters(const ParameterList &source)
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraMultiVector(const RCP< const MultiVectorBase< Scalar > > &mv)
static RCP< const Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraOperator(const RCP< const LinearOpBase< Scalar > > &op)

Generated on Thu Nov 21 2024 09:22:16 for Stratimikos by doxygen 1.8.5