10 #ifndef __Thyra_BelosTpetraSolverAdapter_hpp
11 #define __Thyra_BelosTpetraSolverAdapter_hpp
18 #include "Thyra_MultiVectorBase.hpp"
19 #include "Thyra_TpetraThyraWrappers.hpp"
21 #include "Teuchos_ParameterList.hpp"
22 #include "Teuchos_ParameterListAcceptorDefaultBase.hpp"
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"
37 template<
class SC,
class MV,
class OP>
40 using tpetra_base_solver_type = BelosTpetra::Impl::Krylov<SC>;
44 BelosTpetraKrylov() =
default;
72 bool defaultValues =
true;
73 tpetra_solver->getParameters (*params, defaultValues);
79 int getNumIters()
const override {
return solver_output.numIters; }
80 bool isLOADetected()
const override {
return false; }
83 problem_->setProblem();
87 return tpetra_solver->achievedTol();
92 auto A = problem_->getOperator();
94 tpetra_solver->setMatrix(tpetraA);
96 auto lM = problem_->getLeftPrec();
97 if (lM != Teuchos::null) {
99 tpetra_solver->setLeftPrec(tpetraLM);
101 auto rM = problem_->getRightPrec();
102 if (rM != Teuchos::null) {
104 tpetra_solver->setRightPrec(tpetraRM);
107 auto B = this->problem_->getRHS();
108 auto X = this->problem_->getLHS();
112 solver_output = tpetra_solver->solve(*tpetraX, *tpetraB);
122 BelosTpetra::Impl::SolverOutput<SC> solver_output;
133 template<
class SC,
class MV,
class OP>
134 class BelosTpetraGmres :
public BelosTpetraKrylov<SC, MV, OP> {
136 using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
137 using tpetra_solver_type = BelosTpetra::Impl::Gmres<SC>;
141 krylov_base_solver_type ()
143 this->tpetra_solver =
Teuchos::rcp(
new tpetra_solver_type ());
154 template<
class SC,
class MV,
class OP>
155 class BelosTpetraGmresPipeline :
public BelosTpetraKrylov<SC, MV, OP> {
157 using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
158 using tpetra_solver_type = BelosTpetra::Impl::GmresPipeline<SC>;
161 BelosTpetraGmresPipeline() :
162 krylov_base_solver_type ()
164 this->tpetra_solver =
Teuchos::rcp(
new tpetra_solver_type ());
169 return Teuchos::rcp(
new BelosTpetraGmresPipeline<SC, MV, OP>);
175 template<
class SC,
class MV,
class OP>
176 class BelosTpetraGmresSingleReduce :
public BelosTpetraKrylov<SC, MV, OP> {
178 using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
179 using tpetra_solver_type = BelosTpetra::Impl::GmresSingleReduce<SC>;
182 BelosTpetraGmresSingleReduce() :
183 krylov_base_solver_type ()
185 this->tpetra_solver =
Teuchos::rcp(
new tpetra_solver_type ());
190 return Teuchos::rcp(
new BelosTpetraGmresSingleReduce<SC, MV, OP>);
196 template<
class SC,
class MV,
class OP>
197 class BelosTpetraGmresSstep :
public BelosTpetraKrylov<SC, MV, OP> {
199 using krylov_base_solver_type = BelosTpetraKrylov<SC, MV, OP>;
200 using tpetra_solver_type = BelosTpetra::Impl::GmresSstep<SC>;
203 BelosTpetraGmresSstep() :
204 krylov_base_solver_type ()
206 this->tpetra_solver =
Teuchos::rcp(
new tpetra_solver_type ());
211 return Teuchos::rcp(
new BelosTpetraGmresSstep<SC, MV, OP>);
217 #endif // __Thyra_TpetraSolverAdapter_hpp
virtual void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)=0
bool is_null(const std::shared_ptr< T > &p)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
ParameterList & setParameters(const ParameterList &source)