10 #ifndef THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
11 #define THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
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"
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
157 const bool allowInexactFwdSolve =
false,
159 const bool allowInexactAdjSolve =
false,
160 const double aztecSolverScalar = 1.0
187 bool *isExternalPrec = NULL,
190 bool *allowInexactFwdSolve = NULL,
192 bool *allowInexactAdjSolve = NULL,
193 double *aztecSolverScalar = NULL
257 bool isExternalPrec_;
260 bool allowInexactFwdSolve_;
262 bool allowInexactAdjSolve_;
263 double aztecSolverScalar_;
265 void assertInitialized()
const;
272 #endif // THYRA_AZTECOO_LINEAR_OP_WITH_SOLVE_HPP
virtual bool opSupportedImpl(EOpTransp M_trans) const
RCP< const LinearOpSourceBase< double > > extract_fwdOpSrc()
Extract the forward LinearOpBase<double> object so that it can be modified.
bool isExternalPrec() const
Determine if the preconditioner was external or not.
std::string description() const
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<double> 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)