Stratimikos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Thyra_AztecOOLinearOpWithSolve.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_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
43 #define THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
44 
45 #include "Thyra_LinearOpWithSolveBase.hpp"
46 #include "Thyra_LinearOpSourceBase.hpp"
47 #include "Thyra_EpetraLinearOp.hpp"
48 #include "Thyra_PreconditionerBase.hpp"
50 #include "AztecOO.h"
51 
52 
53 namespace Thyra {
54 
55 
77 class AztecOOLinearOpWithSolve : virtual public LinearOpWithSolveBase<double>
78 {
79 public:
80 
83 
91  const int fwdDefaultMaxIterations = 400,
92  const double fwdDefaultTol = 1e-6,
93  const int adjDefaultMaxIterations = 400,
94  const double adjDefaultTol = 1e-6,
95  const bool outputEveryRhs = false
96  );
97 
99  STANDARD_MEMBER_COMPOSITION_MEMBERS( int, fwdDefaultMaxIterations );
101  STANDARD_MEMBER_COMPOSITION_MEMBERS( double, fwdDefaultTol );
103  STANDARD_MEMBER_COMPOSITION_MEMBERS( int, adjDefaultMaxIterations );
105  STANDARD_MEMBER_COMPOSITION_MEMBERS( double, adjDefaultTol );
107  STANDARD_MEMBER_COMPOSITION_MEMBERS( bool, outputEveryRhs );
108 
182  void initialize(
183  const RCP<const LinearOpBase<double> > &fwdOp,
184  const RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
185  const RCP<const PreconditionerBase<double> > &prec,
186  const bool isExternalPrec,
187  const RCP<const LinearOpSourceBase<double> > &approxFwdOpSrc,
188  const RCP<AztecOO> &aztecFwdSolver,
189  const bool allowInexactFwdSolve = false,
190  const RCP<AztecOO> &aztecAdjSolver = Teuchos::null,
191  const bool allowInexactAdjSolve = false,
192  const double aztecSolverScalar = 1.0
193  );
194 
199 
203 
206  bool isExternalPrec() const;
207 
208 
213 
215  void uninitialize(
216  RCP<const LinearOpBase<double> > *fwdOp = NULL,
217  RCP<const LinearOpSourceBase<double> > *fwdOpSrc = NULL,
218  RCP<const PreconditionerBase<double> > *prec = NULL,
219  bool *isExternalPrec = NULL,
220  RCP<const LinearOpSourceBase<double> > *approxFwdOpSrc = NULL,
221  RCP<AztecOO> *aztecFwdSolver = NULL,
222  bool *allowInexactFwdSolve = NULL,
223  RCP<AztecOO> *aztecAdjSolver = NULL,
224  bool *allowInexactAdjSolve = NULL,
225  double *aztecSolverScalar = NULL
226  );
227 
229 
239 
243  std::string description() const;
245  void describe(
247  const Teuchos::EVerbosityLevel verbLevel
248  ) const;
250 
251 protected:
252 
256  virtual bool opSupportedImpl(EOpTransp M_trans) const;
258  virtual void applyImpl(
259  const EOpTransp M_trans,
260  const MultiVectorBase<double> &X,
261  const Ptr<MultiVectorBase<double> > &Y,
262  const double alpha,
263  const double beta
264  ) const;
266 
270  virtual bool solveSupportsImpl(EOpTransp M_trans) const;
273  EOpTransp M_trans, const SolveMeasureType& solveMeasureType
274  ) const;
276  SolveStatus<double> solveImpl(
277  const EOpTransp M_trans,
278  const MultiVectorBase<double> &B,
279  const Ptr<MultiVectorBase<double> > &X,
280  const Ptr<const SolveCriteria<double> > solveCriteria
281  ) const;
283 
284 private:
285 
296 
297  void assertInitialized() const;
298 
299 };
300 
301 
302 } // namespace Thyra
303 
304 #endif // THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
virtual bool opSupportedImpl(EOpTransp M_trans) const
RCP< const PreconditionerBase< double > > prec_
RCP< const LinearOpSourceBase< double > > extract_fwdOpSrc()
Extract the forward LinearOpBase&lt;double&gt; object so that it can be modified.
bool isExternalPrec() const
Determine if the preconditioner was external or not.
STANDARD_MEMBER_COMPOSITION_MEMBERS(int, fwdDefaultMaxIterations)
The default maximum number of iterations for forward solves.
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha, const double beta) const
virtual bool solveSupportsImpl(EOpTransp M_trans) const
RCP< const VectorSpaceBase< double > > domain() const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
RCP< const LinearOpBase< double > > fwdOp_
RCP< const LinearOpSourceBase< double > > fwdOpSrc_
void uninitialize(RCP< const LinearOpBase< double > > *fwdOp=NULL, RCP< const LinearOpSourceBase< double > > *fwdOpSrc=NULL, RCP< const PreconditionerBase< double > > *prec=NULL, bool *isExternalPrec=NULL, RCP< const LinearOpSourceBase< double > > *approxFwdOpSrc=NULL, RCP< AztecOO > *aztecFwdSolver=NULL, bool *allowInexactFwdSolve=NULL, RCP< AztecOO > *aztecAdjSolver=NULL, bool *allowInexactAdjSolve=NULL, double *aztecSolverScalar=NULL)
Uninitialize.
RCP< const PreconditionerBase< double > > extract_prec()
Extract the preconditioner.
void initialize(const RCP< const LinearOpBase< double > > &fwdOp, const RCP< const LinearOpSourceBase< double > > &fwdOpSrc, const RCP< const PreconditionerBase< double > > &prec, const bool isExternalPrec, const RCP< const LinearOpSourceBase< double > > &approxFwdOpSrc, const RCP< AztecOO > &aztecFwdSolver, const bool allowInexactFwdSolve=false, const RCP< AztecOO > &aztecAdjSolver=Teuchos::null, const bool allowInexactAdjSolve=false, const double aztecSolverScalar=1.0)
Sets up this object.
RCP< const LinearOpBase< double > > clone() const
RCP< const LinearOpSourceBase< double > > approxFwdOpSrc_
RCP< const VectorSpaceBase< double > > range() const
SolveStatus< double > solveImpl(const EOpTransp M_trans, const MultiVectorBase< double > &B, const Ptr< MultiVectorBase< double > > &X, const Ptr< const SolveCriteria< double > > solveCriteria) const
Concrete LinearOpWithSolveBase subclass implemented using AztecOO.
RCP< const LinearOpSourceBase< double > > extract_approxFwdOpSrc()
Extract the approximate forward LinearOpBase&lt;double&gt; object used to build the preconditioner.
AztecOOLinearOpWithSolve(const int fwdDefaultMaxIterations=400, const double fwdDefaultTol=1e-6, const int adjDefaultMaxIterations=400, const double adjDefaultTol=1e-6, const bool outputEveryRhs=false)