|
Stratimikos
Version of the Day
|
Concrete LinearOpWithSolveBase subclass implemented using AztecOO.
More...
#include <Thyra_AztecOOLinearOpWithSolve.hpp>

Constructors/initializers/accessors | |
| AztecOOLinearOpWithSolve (const int fwdDefaultMaxIterations=400, const double fwdDefaultTol=1e-6, const int adjDefaultMaxIterations=400, const double adjDefaultTol=1e-6, const bool outputEveryRhs=false) | |
| STANDARD_MEMBER_COMPOSITION_MEMBERS (int, fwdDefaultMaxIterations) | |
| The default maximum number of iterations for forward solves. More... | |
| STANDARD_MEMBER_COMPOSITION_MEMBERS (double, fwdDefaultTol) | |
| The default solution tolerance on the residual for forward solves. More... | |
| STANDARD_MEMBER_COMPOSITION_MEMBERS (int, adjDefaultMaxIterations) | |
| The default maximum number of iterations for adjoint solves. More... | |
| STANDARD_MEMBER_COMPOSITION_MEMBERS (double, adjDefaultTol) | |
| The default solution tolerance on the residual for adjoint solves. More... | |
| STANDARD_MEMBER_COMPOSITION_MEMBERS (bool, outputEveryRhs) | |
| Determine if output for every RHS will be printed or not. More... | |
| 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. More... | |
| RCP< const LinearOpSourceBase < double > > | extract_fwdOpSrc () |
Extract the forward LinearOpBase<double> object so that it can be modified. More... | |
| RCP< const PreconditionerBase < double > > | extract_prec () |
| Extract the preconditioner. More... | |
| bool | isExternalPrec () const |
| Determine if the preconditioner was external or not. More... | |
| RCP< const LinearOpSourceBase < double > > | extract_approxFwdOpSrc () |
Extract the approximate forward LinearOpBase<double> object used to build the preconditioner. More... | |
| 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. More... | |
Overridden from LinearOpBase | |
| RCP< const VectorSpaceBase < double > > | range () const |
| RCP< const VectorSpaceBase < double > > | domain () const |
| RCP< const LinearOpBase< double > > | clone () const |
| virtual bool | opSupportedImpl (EOpTransp M_trans) const |
| virtual void | applyImpl (const EOpTransp M_trans, const MultiVectorBase< double > &X, const Ptr< MultiVectorBase< double > > &Y, const double alpha, const double beta) const |
Overridden from Teuchos::Describable | |
| std::string | description () const |
| void | describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const |
Overridden from LinearOpWithSolveBase. | |
| virtual bool | solveSupportsImpl (EOpTransp M_trans) const |
| virtual bool | solveSupportsSolveMeasureTypeImpl (EOpTransp M_trans, const SolveMeasureType &solveMeasureType) 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.
This subclass is designed to be very flexible and handle a number of different use cases. It supports forward and optionally adjoint (transpose) solves. I can support inexact solves based on a relative residual norm tolerance or just allow for a default (i.e. tight) linear solve tolerance.
This subclass is not designed to be used directly by users but instead by subclasses of LinearOpWithSolveFactoryBase. One standard implementation that is fairly flexible (and will be make more flexible in the future) is AztecOOLinearOpWithSolveFactory.
This subclass allows for user-defined preconditioners or for built-in aztec preconditioners.
ToDo: Finish documentation!
Definition at line 77 of file Thyra_AztecOOLinearOpWithSolve.hpp.
| Thyra::AztecOOLinearOpWithSolve::AztecOOLinearOpWithSolve | ( | const int | fwdDefaultMaxIterations = 400, |
| const double | fwdDefaultTol = 1e-6, |
||
| const int | adjDefaultMaxIterations = 400, |
||
| const double | adjDefaultTol = 1e-6, |
||
| const bool | outputEveryRhs = false |
||
| ) |
Construct uninitialized but with default option values.
Note, these defaults where taken from NOX::EpetraNew::LinearSystemAztecOO::applyJacobianInverse(...) on 2005/08/15.
Definition at line 196 of file Thyra_AztecOOLinearOpWithSolve.cpp.
| Thyra::AztecOOLinearOpWithSolve::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | int | , |
| fwdDefaultMaxIterations | |||
| ) |
The default maximum number of iterations for forward solves.
| Thyra::AztecOOLinearOpWithSolve::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | double | , |
| fwdDefaultTol | |||
| ) |
The default solution tolerance on the residual for forward solves.
| Thyra::AztecOOLinearOpWithSolve::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | int | , |
| adjDefaultMaxIterations | |||
| ) |
The default maximum number of iterations for adjoint solves.
| Thyra::AztecOOLinearOpWithSolve::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | double | , |
| adjDefaultTol | |||
| ) |
The default solution tolerance on the residual for adjoint solves.
| Thyra::AztecOOLinearOpWithSolve::STANDARD_MEMBER_COMPOSITION_MEMBERS | ( | bool | , |
| outputEveryRhs | |||
| ) |
Determine if output for every RHS will be printed or not.
| void Thyra::AztecOOLinearOpWithSolve::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.
| fwdOp | [in] The forward operator object that defines this objects LinearOpBase interface. interface. |
| fwdOpSrc | [in] The source for the forward operator object fwdOp. This also should be the exact same object that is passed in through a LinearOpWithSolveFactoryBase interface. |
| prec | [in] The original abstract preconditioner object that was passed through the LinearOpWithSolveFactoryBase interface. This object is not used for anything and can be set as prec==Teuchos::null. |
| isExternalPrec | [in] True if the precondition was created externally from the LinearOpWithSolveFactoryBase object, false otherwise. |
| approxFwdOpSrc | [in] The source for the original abstract approximate forward operator object that was passed through the LinearOpWithSolveFactoryBase interface. This object is not used for anything and can be set as approxFwdOpSrc==Teuchos::null. |
| aztecFwdSolver | [in] The AztecOO object used to perform forward solves. This object must be be ready to call aztecFwdSolver->SetRHS() and aztecFwdSolver->SetLHS() and then call aztecFwdSolver->Solve(). |
| allowInexactFwdSolve | [in] Determines if this->solveSupportsSolveTolType(NOTRANS,SOLVE_TOL_REL_RESIDUAL_NORM) returns true or not. With the current design, an inexact forward solve can not be supported if there is left scaling or a left preconditioner aggregated with *aztecFwdOp. |
| aztecAdjSolver | [in] The AztecOO object used to perform adjoint solves. This object must be be ready to call aztecAdjSolver->SetRHS() and aztecAdjSolver->SetLHS() and then call aztecAdjSolver->Solve(). |
| allowInexactAdjSolve | [in] Determines if this->solveSupportsSolveTolType(TRANS,SOLVE_TOL_REL_RESIDUAL_NORM) returns true or not. With the current design, an inexact forward solve can not be supported if there is left scaling or a left preconditioner aggregated with *aztecFwdOp. |
| linearSystemTransformer | [in] This is a transformation object that is called to pre-preprocess the linear problem before a forward and adjoint linear solver and post-process the linear problem after forward and adjoint linear solve. This abstract object is used to deal with scaling and aggregated preconditioners. It is what makes this implementation fairly flexible. |
Preconditions:
fwdOp.get()!=NULL fwdOpSrc.get()!=NULL fwdFwdSolver.get()!=NULL Postconditions:
this->range() == fwdOp->range() this->domain() == fwdOp->domain() this->opSupports(M_trans) == opSupports(*fwdOp,M_trans) this->solveSupportsTrans(M_trans) == (aztecAdjSolver.get()!=NULL) this->solveSupportsSolveTolType([NOTRANS,CONJ], SolveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MESURE_NORM_RHS)) == allowInexactFwdSolve this->solveSupportsSolveTolType([TRANS,CONJTRANS], SolveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MESURE_NORM_RHS)) == (aztecAdjSolver.get()!=NULL&&allowInexactAdjSolve) ToDo: Finish documentation!
Definition at line 215 of file Thyra_AztecOOLinearOpWithSolve.cpp.
| RCP< const LinearOpSourceBase< double > > Thyra::AztecOOLinearOpWithSolve::extract_fwdOpSrc | ( | ) |
Extract the forward LinearOpBase<double> object so that it can be modified.
Definition at line 250 of file Thyra_AztecOOLinearOpWithSolve.cpp.
| RCP< const PreconditionerBase< double > > Thyra::AztecOOLinearOpWithSolve::extract_prec | ( | ) |
Extract the preconditioner.
Definition at line 260 of file Thyra_AztecOOLinearOpWithSolve.cpp.
| bool Thyra::AztecOOLinearOpWithSolve::isExternalPrec | ( | ) | const |
Determine if the preconditioner was external or not.
Definition at line 269 of file Thyra_AztecOOLinearOpWithSolve.cpp.
| RCP< const LinearOpSourceBase< double > > Thyra::AztecOOLinearOpWithSolve::extract_approxFwdOpSrc | ( | ) |
Extract the approximate forward LinearOpBase<double> object used to build the preconditioner.
Definition at line 276 of file Thyra_AztecOOLinearOpWithSolve.cpp.
| void Thyra::AztecOOLinearOpWithSolve::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.
Definition at line 285 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 326 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
virtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 333 of file Thyra_AztecOOLinearOpWithSolve.cpp.
| RCP< const LinearOpBase< double > > Thyra::AztecOOLinearOpWithSolve::clone | ( | ) | const |
Definition at line 340 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 349 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
virtual |
Reimplemented from Teuchos::Describable.
Definition at line 362 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
protectedvirtual |
Implements Thyra::LinearOpBase< Scalar >.
Definition at line 440 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
protectedvirtual |
Definition at line 446 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
protectedvirtual |
Reimplemented from Thyra::LinearOpWithSolveBase< double >.
Definition at line 462 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
protectedvirtual |
Reimplemented from Thyra::LinearOpWithSolveBase< double >.
Definition at line 470 of file Thyra_AztecOOLinearOpWithSolve.cpp.
|
protectedvirtual |
Implements Thyra::LinearOpWithSolveBase< double >.
Definition at line 544 of file Thyra_AztecOOLinearOpWithSolve.cpp.
1.8.5