Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_BelosTpetrasSolverAdapter.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Stratimikos: Thyra-based strategies for linear solvers
5 // Copyright (2006) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
38 //
39 // ***********************************************************************
40 // @HEADER
41 
42 #ifndef __Thyra_BelosTpetraSolverAdapter_hpp
43 #define __Thyra_BelosTpetraSolverAdapter_hpp
44 
45 #include "BelosConfigDefs.hpp"
46 #include "BelosLinearProblem.hpp"
47 // BelosThyraAdapter.hpp only includes this file if HAVE_BELOS_TSQR is
48 // defined. Thus, it's OK to include TSQR header files here.
49 
50 #include "Thyra_MultiVectorBase.hpp"
52 
55 
56 #include "Tpetra_CrsMatrix.hpp"
57 #include "BelosTpetraAdapter.hpp"
58 #include "solvers/Belos_Tpetra_Krylov_parameters.hpp"
59 #include "solvers/Belos_Tpetra_Krylov.hpp"
60 #include "solvers/Belos_Tpetra_Gmres.hpp"
61 #include "solvers/Belos_Tpetra_GmresPipeline.hpp"
62 #include "solvers/Belos_Tpetra_GmresSingleReduce.hpp"
63 #include "solvers/Belos_Tpetra_GmresSstep.hpp"
64 
65 namespace Thyra {
66 
69  template<class SC, class MV, class OP>
70  class BelosTpetraKrylov : public Belos::SolverManager<SC, MV, OP> {
71  public:
72  using tpetra_base_solver_type = BelosTpetra::Impl::Krylov<SC>;
74 
76  BelosTpetraKrylov() = default;
77 
79  virtual Teuchos::RCP<Belos::SolverManager<SC, MV, OP> > clone () const override = 0;
80 
82  void setProblem( const Teuchos::RCP<Belos::LinearProblem<SC, MV, OP> > &problem ) override {
83  problem_ = problem;
84  }
85 
86  const Belos::LinearProblem<SC, MV, OP>& getProblem() const override {
87  return *problem_;
88  }
89 
91  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override {
92  params_ = params;
93  tpetra_solver->setParameters(*params);
94  }
95 
97  return params_;
98  }
99 
101 
102  Teuchos::RCP<Teuchos::ParameterList> params = Teuchos::parameterList();
103 
104  bool defaultValues = true;
105  tpetra_solver->getParameters (*params, defaultValues);
106 
107  return params;
108  }
109 
111  int getNumIters() const override { return solver_output.numIters; }
112  bool isLOADetected() const override { return false; }
113  void reset( const Belos::ResetType type ) override {
114  if ((type & Belos::Problem) && !Teuchos::is_null(problem_))
115  problem_->setProblem();
116  }
117 
119  return tpetra_solver->achievedTol();
120  }
121 
123  Belos::ReturnType solve() override {
124  auto A = problem_->getOperator();
125  auto tpetraA = converter::getConstTpetraOperator(A);
126  tpetra_solver->setMatrix(tpetraA);
127 
128  auto lM = problem_->getLeftPrec();
129  if (lM != Teuchos::null) {
130  auto tpetraLM = converter::getConstTpetraOperator(lM);
131  tpetra_solver->setLeftPrec(tpetraLM);
132  }
133  auto rM = problem_->getRightPrec();
134  if (rM != Teuchos::null) {
135  auto tpetraRM = converter::getConstTpetraOperator(rM);
136  tpetra_solver->setRightPrec(tpetraRM);
137  }
138 
139  auto B = this->problem_->getRHS();
140  auto X = this->problem_->getLHS();
141  auto tpetraB = converter::getConstTpetraMultiVector(B);
142  auto tpetraX = converter::getTpetraMultiVector(X);
143 
144  solver_output = tpetra_solver->solve(*tpetraX, *tpetraB);
145 
146  // copy output
147  Belos::ReturnType belos_output = (solver_output.converged ? Belos::Converged : Belos::Unconverged);
148  return belos_output;
149  }
150 
151  protected:
152 
154  BelosTpetra::Impl::SolverOutput<SC> solver_output;
155 
158 
161  };
162 
163  //
164  // Gmres
165  template<class SC, class MV, class OP>
166  class BelosTpetraGmres : public BelosTpetraKrylov<SC, MV, OP> {
167  public:
169  using tpetra_solver_type = BelosTpetra::Impl::Gmres<SC>;
170 
174  {
176  }
177 
181  }
182  };
183 
184  //
185  // Pipelined Gmres
186  template<class SC, class MV, class OP>
187  class BelosTpetraGmresPipeline : public BelosTpetraKrylov<SC, MV, OP> {
188  public:
190  using tpetra_solver_type = BelosTpetra::Impl::GmresPipeline<SC>;
191 
195  {
197  }
198 
202  }
203  };
204 
205  //
206  // SingleReduce Gmres
207  template<class SC, class MV, class OP>
208  class BelosTpetraGmresSingleReduce : public BelosTpetraKrylov<SC, MV, OP> {
209  public:
211  using tpetra_solver_type = BelosTpetra::Impl::GmresSingleReduce<SC>;
212 
216  {
218  }
219 
223  }
224  };
225 
226  //
227  // s-step Gmres
228  template<class SC, class MV, class OP>
229  class BelosTpetraGmresSstep : public BelosTpetraKrylov<SC, MV, OP> {
230  public:
232  using tpetra_solver_type = BelosTpetra::Impl::GmresSstep<SC>;
233 
237  {
239  }
240 
244  }
245  };
246 
247 } // namespace Thyra
248 
249 #endif // __Thyra_TpetraSolverAdapter_hpp
250 
Teuchos::RCP< Teuchos::ParameterList > params_
Current parameter list.
void reset(const Belos::ResetType type) override
BelosTpetra::Impl::GmresSingleReduce< SC > tpetra_solver_type
bool is_null(const std::shared_ptr< T > &p)
Teuchos::RCP< tpetra_base_solver_type > tpetra_solver
BelosTpetra::Impl::Gmres< SC > tpetra_solver_type
BelosTpetraKrylov()=default
constructor
Teuchos::RCP< Belos::SolverManager< SC, MV, OP > > clone() const override
clone for Inverted Injection (DII)
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraMultiVector(const RCP< MultiVectorBase< Scalar > > &mv)
virtual Teuchos::RCP< Belos::SolverManager< SC, MV, OP > > clone() const override=0
clone for Inverted Injection (DII)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
const Belos::LinearProblem< SC, MV, OP > & getProblem() const override
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Belos::ReturnType solve() override
solve
BelosTpetra::Impl::Krylov< SC > tpetra_base_solver_type
Teuchos::RCP< Belos::SolverManager< SC, MV, OP > > clone() const override
clone for Inverted Injection (DII)
Teuchos::RCP< Belos::LinearProblem< SC, MV, OP > > problem_
The linear problem to solve.
int getNumIters() const override
Get the iteration count for the most recent call to solve().
void setProblem(const Teuchos::RCP< Belos::LinearProblem< SC, MV, OP > > &problem) override
set/get problem
BelosTpetra::Impl::SolverOutput< SC > solver_output
Teuchos::RCP< Belos::SolverManager< SC, MV, OP > > clone() const override
clone for Inverted Injection (DII)
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraMultiVector(const RCP< const MultiVectorBase< Scalar > > &mv)
BelosTpetra::Impl::GmresSstep< SC > tpetra_solver_type
static RCP< const Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraOperator(const RCP< const LinearOpBase< Scalar > > &op)
BelosTpetra::Impl::GmresPipeline< SC > tpetra_solver_type
Teuchos::ScalarTraits< SC >::magnitudeType achievedTol() const override
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
set/get parameters
Teuchos::RCP< Belos::SolverManager< SC, MV, OP > > clone() const override
clone for Inverted Injection (DII)