Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_AmesosLinearOpWithSolve.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 #ifndef THYRA_AMESOS_LINEAR_OP_WITH_SOLVE_HPP
45 #define THYRA_AMESOS_LINEAR_OP_WITH_SOLVE_HPP
46 
47 #include "Thyra_LinearOpWithSolveBase.hpp"
48 #include "Thyra_LinearOpSourceBase.hpp"
50 #include "Epetra_LinearProblem.h"
51 #include "Amesos_BaseSolver.h"
52 
53 
54 namespace Thyra {
55 
56 
70 class AmesosLinearOpWithSolve : virtual public LinearOpWithSolveBase<double>
71 {
72 public:
73 
76 
79 
82  const Teuchos::RCP<const LinearOpBase<double> > &fwdOp,
83  const Teuchos::RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
84  const Teuchos::RCP<Epetra_LinearProblem> &epetraLP,
85  const Teuchos::RCP<Amesos_BaseSolver> &amesosSolver,
86  const EOpTransp amesosSolverTransp,
87  const double amesosSolverScalar
88  );
89 
130  void initialize(
131  const Teuchos::RCP<const LinearOpBase<double> > &fwdOp,
132  const Teuchos::RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
133  const Teuchos::RCP<Epetra_LinearProblem> &epetraLP,
134  const Teuchos::RCP<Amesos_BaseSolver> &amesosSolver,
135  const EOpTransp amesosSolverTransp,
136  const double amesosSolverScalar
137  );
138 
147 
150 
153 
156 
159 
161  EOpTransp get_amesosSolverTransp() const;
162 
164  double get_amesosSolverScalar() const;
165 
168  void uninitialize(
169  Teuchos::RCP<const LinearOpBase<double> > *fwdOp = NULL,
170  Teuchos::RCP<const LinearOpSourceBase<double> > *fwdOpSrc = NULL,
171  Teuchos::RCP<Epetra_LinearProblem> *epetraLP = NULL,
172  Teuchos::RCP<Amesos_BaseSolver> *amesosSolver = NULL,
173  EOpTransp *amesosSolverTransp = NULL,
174  double *amesosSolverScalar = NULL
175  );
176 
178 
188 
192  std::string description() const;
194  void describe(
196  const Teuchos::EVerbosityLevel verbLevel
197  ) const;
199 
200 protected:
201 
205  virtual bool opSupportedImpl(EOpTransp M_trans) const;
207  virtual void applyImpl(
208  const EOpTransp M_trans,
209  const MultiVectorBase<double> &X,
210  const Ptr<MultiVectorBase<double> > &Y,
211  const double alpha,
212  const double beta
213  ) const;
215 
219  virtual bool solveSupportsImpl(EOpTransp M_trans) const;
222  EOpTransp M_trans, const SolveMeasureType& solveMeasureType
223  ) const;
225  SolveStatus<double> solveImpl(
226  const EOpTransp M_trans,
227  const MultiVectorBase<double> &B,
228  const Ptr<MultiVectorBase<double> > &X,
229  const Ptr<const SolveCriteria<double> > solveCriteria
230  ) const;
232 
233 private:
234 
241 
242  void assertInitialized() const;
243 
244 };
245 
246 // ///////////////////////////
247 // Inline members
248 
249 inline
252 {
253  return fwdOp_;
254 }
255 
256 inline
259 {
260  return fwdOpSrc_;
261 }
262 
263 inline
266 {
267  return epetraLP_;
268 }
269 
270 inline
273 {
274  return amesosSolver_;
275 }
276 
277 inline
279 {
280  return amesosSolverTransp_;
281 }
282 
283 inline
285 {
286  return amesosSolverScalar_;
287 }
288 
289 } // namespace Thyra
290 
291 #endif // THYRA_AMESOS_LINEAR_OP_WITH_SOLVE_HPP
Teuchos::RCP< Epetra_LinearProblem > epetraLP_
Teuchos::RCP< const LinearOpSourceBase< double > > extract_fwdOpSrc()
Extract the LinearOpSourceBase&lt;double&gt; object so that it can be modified.
Teuchos::RCP< Amesos_BaseSolver > amesosSolver_
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual bool opSupportedImpl(EOpTransp M_trans) const
Teuchos::RCP< Amesos_BaseSolver > get_amesosSolver() const
virtual bool solveSupportsImpl(EOpTransp M_trans) const
SolveStatus< double > solveImpl(const EOpTransp M_trans, const MultiVectorBase< double > &B, const Ptr< MultiVectorBase< double > > &X, const Ptr< const SolveCriteria< double > > solveCriteria) const
Teuchos::RCP< const VectorSpaceBase< double > > range() const
Teuchos::RCP< const LinearOpSourceBase< double > > fwdOpSrc_
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const
void initialize(const Teuchos::RCP< const LinearOpBase< double > > &fwdOp, const Teuchos::RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const Teuchos::RCP< Epetra_LinearProblem > &epetraLP, const Teuchos::RCP< Amesos_BaseSolver > &amesosSolver, const EOpTransp amesosSolverTransp, const double amesosSolverScalar)
First initialization.
Concrete LinearOpWithSolveBase subclass that adapts any Amesos_BaseSolver object. ...
Teuchos::RCP< const LinearOpBase< double > > get_fwdOp() const
void uninitialize(Teuchos::RCP< const LinearOpBase< double > > *fwdOp=NULL, Teuchos::RCP< const LinearOpSourceBase< double > > *fwdOpSrc=NULL, Teuchos::RCP< Epetra_LinearProblem > *epetraLP=NULL, Teuchos::RCP< Amesos_BaseSolver > *amesosSolver=NULL, EOpTransp *amesosSolverTransp=NULL, double *amesosSolverScalar=NULL)
Uninitialize.
Teuchos::RCP< const LinearOpSourceBase< double > > get_fwdOpSrc() const
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha, const double beta) const
Teuchos::RCP< const LinearOpBase< double > > fwdOp_
Teuchos::RCP< const LinearOpBase< double > > clone() const
Teuchos::RCP< Epetra_LinearProblem > get_epetraLP() const
Teuchos::RCP< const VectorSpaceBase< double > > domain() const
AmesosLinearOpWithSolve()
Construct to uninitialized.