Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Thyra_AztecOOLinearOpWithSolve.hpp
1 // @HEADER
2 // *****************************************************************************
3 // Stratimikos: Thyra-based strategies for linear solvers
4 //
5 // Copyright 2006 NTESS and the Stratimikos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
11 #define THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
12 
13 #include "Thyra_LinearOpWithSolveBase.hpp"
14 #include "Thyra_LinearOpSourceBase.hpp"
15 #include "Thyra_EpetraLinearOp.hpp"
16 #include "Thyra_PreconditionerBase.hpp"
17 #include "Teuchos_StandardMemberCompositionMacros.hpp"
18 #include "AztecOO.h"
19 
20 
21 namespace Thyra {
22 
23 
45 class AztecOOLinearOpWithSolve : virtual public LinearOpWithSolveBase<double>
46 {
47 public:
48 
51 
59  const int fwdDefaultMaxIterations = 400,
60  const double fwdDefaultTol = 1e-6,
61  const int adjDefaultMaxIterations = 400,
62  const double adjDefaultTol = 1e-6,
63  const bool outputEveryRhs = false
64  );
65 
67  STANDARD_MEMBER_COMPOSITION_MEMBERS( int, fwdDefaultMaxIterations );
69  STANDARD_MEMBER_COMPOSITION_MEMBERS( double, fwdDefaultTol );
71  STANDARD_MEMBER_COMPOSITION_MEMBERS( int, adjDefaultMaxIterations );
73  STANDARD_MEMBER_COMPOSITION_MEMBERS( double, adjDefaultTol );
75  STANDARD_MEMBER_COMPOSITION_MEMBERS( bool, outputEveryRhs );
76 
150  void initialize(
151  const RCP<const LinearOpBase<double> > &fwdOp,
152  const RCP<const LinearOpSourceBase<double> > &fwdOpSrc,
153  const RCP<const PreconditionerBase<double> > &prec,
154  const bool isExternalPrec,
155  const RCP<const LinearOpSourceBase<double> > &approxFwdOpSrc,
156  const RCP<AztecOO> &aztecFwdSolver,
157  const bool allowInexactFwdSolve = false,
158  const RCP<AztecOO> &aztecAdjSolver = Teuchos::null,
159  const bool allowInexactAdjSolve = false,
160  const double aztecSolverScalar = 1.0
161  );
162 
167 
171 
174  bool isExternalPrec() const;
175 
176 
181 
183  void uninitialize(
184  RCP<const LinearOpBase<double> > *fwdOp = NULL,
185  RCP<const LinearOpSourceBase<double> > *fwdOpSrc = NULL,
186  RCP<const PreconditionerBase<double> > *prec = NULL,
187  bool *isExternalPrec = NULL,
188  RCP<const LinearOpSourceBase<double> > *approxFwdOpSrc = NULL,
189  RCP<AztecOO> *aztecFwdSolver = NULL,
190  bool *allowInexactFwdSolve = NULL,
191  RCP<AztecOO> *aztecAdjSolver = NULL,
192  bool *allowInexactAdjSolve = NULL,
193  double *aztecSolverScalar = NULL
194  );
195 
197 
207 
211  std::string description() const;
213  void describe(
215  const Teuchos::EVerbosityLevel verbLevel
216  ) const;
218 
219 protected:
220 
224  virtual bool opSupportedImpl(EOpTransp M_trans) const;
226  virtual void applyImpl(
227  const EOpTransp M_trans,
228  const MultiVectorBase<double> &X,
229  const Ptr<MultiVectorBase<double> > &Y,
230  const double alpha,
231  const double beta
232  ) const;
234 
238  virtual bool solveSupportsImpl(EOpTransp M_trans) const;
241  EOpTransp M_trans, const SolveMeasureType& solveMeasureType
242  ) const;
245  const EOpTransp M_trans,
246  const MultiVectorBase<double> &B,
247  const Ptr<MultiVectorBase<double> > &X,
248  const Ptr<const SolveCriteria<double> > solveCriteria
249  ) const;
251 
252 private:
253 
257  bool isExternalPrec_;
258  RCP<const LinearOpSourceBase<double> > approxFwdOpSrc_;
259  RCP<AztecOO> aztecFwdSolver_;
260  bool allowInexactFwdSolve_;
261  RCP<AztecOO> aztecAdjSolver_;
262  bool allowInexactAdjSolve_;
263  double aztecSolverScalar_;
264 
265  void assertInitialized() const;
266 
267 };
268 
269 
270 } // namespace Thyra
271 
272 #endif // THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
virtual bool opSupportedImpl(EOpTransp M_trans) const
EOpTransp
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
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 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)

Generated on Thu Nov 21 2024 09:22:16 for Stratimikos by doxygen 1.8.5