Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_BelosLinearOpWithSolve_decl.hpp
Go to the documentation of this file.
1 /*
2 // @HEADER
3 // ***********************************************************************
4 //
5 // Stratimikos: Thyra-based strategies for linear solvers
6 // Copyright (2006) Sandia Corporation
7 //
8 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 // license for use of this work by or on behalf of the U.S. Government.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Roscoe A. Bartlett (rabartl@sandia.gov)
39 //
40 // ***********************************************************************
41 // @HEADER
42 */
43 
44 
45 #ifndef THYRA_BELOS_LINEAR_OP_WITH_SOLVE_DECL_HPP
46 #define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_DECL_HPP
47 
48 #include "Thyra_LinearOpWithSolveBase.hpp"
49 #include "Thyra_LinearOpSourceBase.hpp"
50 #include "BelosSolverManager.hpp"
51 #include "BelosThyraAdapter.hpp"
53 
54 
55 namespace Thyra {
56 
57 
65 template<class Scalar>
66 class BelosLinearOpWithSolve : virtual public LinearOpWithSolveBase<Scalar>
67 {
68 public:
69 
72 
74  typedef MultiVectorBase<Scalar> MV_t;
76  typedef LinearOpBase<Scalar> LO_t;
77 
79 
82 
85 
125  void initialize(
126  const RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > &lp,
127  const RCP<Teuchos::ParameterList> &solverPL,
128  const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver,
129  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
130  const RCP<const PreconditionerBase<Scalar> > &prec,
131  const bool isExternalPrec,
132  const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
133  const ESupportSolveUse &supportSolveUse,
134  const int convergenceTestFrequency
135  );
136 
139 
142 
144  bool isExternalPrec() const;
145 
148 
150  ESupportSolveUse supportSolveUse() const;
151 
156  void uninitialize(
157  RCP<Belos::LinearProblem<Scalar,MV_t,LO_t> > *lp = NULL,
158  RCP<Teuchos::ParameterList> *solverPL = NULL,
159  RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver = NULL,
160  RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc = NULL,
161  RCP<const PreconditionerBase<Scalar> > *prec = NULL,
162  bool *isExternalPrec = NULL,
163  RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc = NULL,
164  ESupportSolveUse *supportSolveUse = NULL
165  );
166 
168 
178 
182  std::string description() const;
184  void describe(
186  const Teuchos::EVerbosityLevel verbLevel
187  ) const;
189 
192 
194  void setParameterList(RCP<Teuchos::ParameterList> const& paramList);
203 
205 
206 protected:
207 
211  virtual bool opSupportedImpl(EOpTransp M_trans) const;
213  virtual void applyImpl(
214  const EOpTransp M_trans,
215  const MultiVectorBase<Scalar> &X,
216  const Ptr<MultiVectorBase<Scalar> > &Y,
217  const Scalar alpha,
218  const Scalar beta
219  ) const;
221 
225  virtual bool solveSupportsImpl(EOpTransp M_trans) const;
227  virtual bool solveSupportsNewImpl(EOpTransp transp,
228  const Ptr<const SolveCriteria<Scalar> > solveCriteria) const;
231  EOpTransp M_trans, const SolveMeasureType& solveMeasureType
232  ) const;
234  virtual SolveStatus<Scalar> solveImpl(
235  const EOpTransp transp,
236  const MultiVectorBase<Scalar> &B,
237  const Ptr<MultiVectorBase<Scalar> > &X,
238  const Ptr<const SolveCriteria<Scalar> > solveCriteria
239  ) const;
241 
242 private:
243 
244  // ///////////////////////////////
245  // Private data members
246 
247 
252 
257  ESupportSolveUse supportSolveUse_;
258 
260 
261  void assertInitialized() const;
262 
264  mutable int counter_;
265 
266 };
267 
268 
269 } // namespace Thyra
270 
271 
272 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_DECL_HPP
virtual bool solveSupportsImpl(EOpTransp M_trans) const
RCP< const VectorSpaceBase< Scalar > > domain() const
RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > lp_
RCP< const Teuchos::ParameterList > getValidParameters() const
virtual SolveStatus< Scalar > solveImpl(const EOpTransp transp, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
void initialize(const RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > &lp, const RCP< Teuchos::ParameterList > &solverPL, const RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > &iterativeSolver, const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const RCP< const PreconditionerBase< Scalar > > &prec, const bool isExternalPrec, const RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, const ESupportSolveUse &supportSolveUse, const int convergenceTestFrequency)
Initializes given precreated solver objects.
Concrete LinearOpWithSolveBase subclass in terms of Belos.
RCP< const PreconditionerBase< Scalar > > extract_prec()
virtual bool solveSupportsNewImpl(EOpTransp transp, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
RCP< const PreconditionerBase< Scalar > > prec_
RCP< const VectorSpaceBase< Scalar > > range() const
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
void setParameterList(RCP< Teuchos::ParameterList > const &paramList)
Thyra specializations of MultiVecTraits and OperatorTraits.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()
Teuchos::ScalarTraits< Scalar >::magnitudeType defaultTol_
BelosLinearOpWithSolve()
Construct to unintialize.
void uninitialize(RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > *lp=NULL, RCP< Teuchos::ParameterList > *solverPL=NULL, RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > *iterativeSolver=NULL, RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc=NULL, RCP< const PreconditionerBase< Scalar > > *prec=NULL, bool *isExternalPrec=NULL, RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc=NULL, ESupportSolveUse *supportSolveUse=NULL)
Uninitializes and returns stored quantities.
RCP< const LinearOpSourceBase< Scalar > > approxFwdOpSrc_
RCP< const LinearOpBase< Scalar > > clone() const
RCP< const Teuchos::ParameterList > getParameterList() const
RCP< Teuchos::ParameterList > unsetParameterList()
virtual bool opSupportedImpl(EOpTransp M_trans) const
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const
RCP< const LinearOpSourceBase< Scalar > > fwdOpSrc_
RCP< Teuchos::ParameterList > getNonconstParameterList()
RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > iterativeSolver_