AztecOO
Development
|
AztecOO_StatusTestResNorm: An implementation of AztecOO_StatusTest using a family of residual norms. More...
#include <AztecOO_StatusTestResNorm.h>
Public Types | |
enum | ResType { Implicit, Explicit } |
Select how the residual vector is produced. More... | |
enum | NormType { OneNorm, TwoNorm, InfNorm } |
Select the norm type, where is the norm argument and is an optionally defined weighting vector. More... | |
enum | ScaleType { NormOfRHS, NormOfInitRes, None, UserProvided } |
Select the scale type. More... | |
Public Member Functions | |
AztecOO_StatusTestResNorm (const Epetra_Operator &Operator, const Epetra_Vector &LHS, const Epetra_Vector &RHS, double Tolerance) | |
Constructor. More... | |
virtual | ~AztecOO_StatusTestResNorm () |
Destructor. | |
int | DefineResForm (ResType TypeOfResidual, NormType TypeOfNorm, Epetra_Vector *Weights=0) |
Define form of the residual, its norm and optional weighting vector. More... | |
int | DefineScaleForm (ScaleType TypeOfScaling, NormType TypeOfNorm, Epetra_Vector *Weights=0, double ScaleValue=1.0) |
Define form of the scaling, its norm, its optional weighting vector, or, alternatively, define an explicit value. More... | |
int | ResetTolerance (double Tolerance) |
Reset the value of the tolerance. More... | |
int | SetMaxNumExtraIterations (int maxNumExtraIterations) |
Set the maximum number of extra iterations that will be performed in case the implicit residual succeeds that the true residual fails (default is 0). More... | |
int | GetMaxNumExtraIterations () |
bool | ResidualVectorRequired () const |
Indicates if residual vector is required by this convergence test. More... | |
AztecOO_StatusType | CheckStatus (int CurrentIter, Epetra_MultiVector *CurrentResVector, double CurrentResNormEst, bool SolutionUpdated) |
Check convergence status: Unconverged, Converged, Failed. More... | |
AztecOO_StatusType | GetStatus () const |
Return the result of the most recent checkStatus call. | |
std::ostream & | Print (std::ostream &stream, int indent=0) const |
Output formatted description of stopping test to output stream. | |
void | ResetStatus () |
Reset state of status test object. | |
double | GetTolerance () const |
Returns the value of the tolerance, , set in the constructor. | |
double | GetTestValue () const |
Returns the test value, , computed in most recent call to CheckStatus. | |
double | GetResNormValue () const |
Returns the residual norm value, , computed in most recent call to CheckStatus. | |
double | GetScaledNormValue () const |
Returns the scaled norm value, . | |
Public Member Functions inherited from AztecOO_StatusTest | |
virtual void | PrintStatus (std::ostream &os, AztecOO_StatusType type) const |
AztecOO_StatusTest () | |
Constructor. | |
virtual | ~AztecOO_StatusTest () |
Destructor. | |
Protected Member Functions | |
double | ComputeNorm (const Epetra_Vector &vec, NormType typeofnorm) |
AztecOO_StatusTestResNorm: An implementation of AztecOO_StatusTest using a family of residual norms.
AztecOO_StatusTestResNorm is an implementation of AztecOO_StatusTest that allows a user to construct one of a family of residual tests for use as a status/convergence test for AztecOO. The form of the test is
where
Select how the residual vector is produced.
Enumerator | |
---|---|
Implicit |
Use the residual vector produced by the iterative solver. |
Explicit |
Explicitly compute the residual vector r = b - A*x using the linear problem. NOTE: If if the ResType is set to Explicit, if the SolutionUpdated argument to CheckStatus() is false, we will not be able to compute an explicit result and we will proceed with the residual implicitly defined. |
AztecOO_StatusTestResNorm::AztecOO_StatusTestResNorm | ( | const Epetra_Operator & | Operator, |
const Epetra_Vector & | LHS, | ||
const Epetra_Vector & | RHS, | ||
double | Tolerance | ||
) |
Constructor.
The constructor takes a single argument specifying the tolerance ( ). If none of the form definition methods are called, we use as the stopping criterion, where uses the least costly form of the 2-norm of residual available from the iterative method and is the corresponding norm of the initial residual. The least costly form of the 2-norm depends on the chosen iterative method. Most Krylov methods produce the preconditioned residual vector in a form that would be exact in infinite precision arithmetic. This vector may be different from the true residual either because left scaling or preconditioning was used, or because round-off error has introduced significant error, or both.
Operator | (In) The original linear operator that was passed in to the AztecOO solver object. |
LHS | (In) The original left hand side vector that was passed in to the AztecOO solver object. NOTE: AztecOO accepts multivector objects, but AztecOO_StatusTestResNorm does not handle residual tests for multiple vectors. Most AztecOO users tend to have Epetra_Vectors, even though they are passing them in as Epetra_MultiVectors, so this should not be an issue. If you are truly using Epetra_MultiVector objects, remember that for a multivector object mv, the ith vector can be extracted as mv(i). |
RHS | (In) The original right hand side vector that was passed in to the AztecOO solver object. See note for LHS. |
Tolerance | (In) A double value that is used to test the convergence criterion. Can be reset using ResetTolerance(). |
|
virtual |
Check convergence status: Unconverged, Converged, Failed.
This method checks to see if the convergence criteria are met. Depending on how the residual test is constructed this method will return the appropriate status type.
CurrentIter | (In) Current iteration of iterative method. |
CurrentResVector | (In) The current residuals of the iterative process. |
CurrentResNormEst | (In) Estimate of the two-norm of the residual. The value will be set to -1.0 if no estimate is available. |
SolutionUpdated | (In) If this argument is true, then the solution vector that is part of the Epetra_LinearProblem object being solved is consistent with the residual. |
Implements AztecOO_StatusTest.
References Epetra_Operator::Apply(), Converged, Explicit, Failed, Implicit, NaN, NormOfInitRes, NormOfRHS, PartialFailed, TwoNorm, and Unconverged.
int AztecOO_StatusTestResNorm::DefineResForm | ( | ResType | TypeOfResidual, |
NormType | TypeOfNorm, | ||
Epetra_Vector * | Weights = 0 |
||
) |
Define form of the residual, its norm and optional weighting vector.
This method defines the form of . We specify:
int AztecOO_StatusTestResNorm::DefineScaleForm | ( | ScaleType | TypeOfScaling, |
NormType | TypeOfNorm, | ||
Epetra_Vector * | Weights = 0 , |
||
double | ScaleValue = 1.0 |
||
) |
Define form of the scaling, its norm, its optional weighting vector, or, alternatively, define an explicit value.
This method defines the form of how the residual is scaled (if at all). It operates in two modes:
User-provided scaling value:
References NormOfInitRes, and TwoNorm.
|
inline |
Return the maximum number of extra iterations that are performed if the implicit residual succeeds while the true residual fails.
|
inline |
Reset the value of the tolerance.
We allow the tolerance to be reset for cases where, in the process of testing the residual, we find that the initial tolerance was too tight or too lax.
|
virtual |
Indicates if residual vector is required by this convergence test.
The value returned by this method will depend on several factors. Once an AztecOO_StatusTestResNorm object is constructed and the DefineResForm and DefineScaleForm methods are optionally called, this method can tested. For most Krylov solvers, there is no extra cost to providing the residual vector. However, GMRES and Transpose-free QMR will need to explicitly compute this vector if ResidualVectorRequired() returns true, so this is an extra cost for these two iterative methods.
Implements AztecOO_StatusTest.
|
inline |
Set the maximum number of extra iterations that will be performed in case the implicit residual succeeds that the true residual fails (default is 0).
In some instance,especially with GMRES, the implicitly computed residual is an optimistic estimate of the true residual. In these cases, especially when the tolerance is set very small, the iterative solver can never satisfy the tolerance with the explicit residual, so we allow the user to limit the number of extra iterations that will be performed. If the implicit residual is satisfied, then the value set here will determine exactly how many extra iterations will be allowed to try and converge the explicit residual. This value has no impact status tests that are based on explicit residuals.
maxNumExtraIterations | (In) Maximum number of extra iterations that will be performed in case the implicit residual succeeds that the true residual fails. |