#include <AztecOO_StatusTestResNorm.h>
Inheritance diagram for AztecOO_StatusTestResNorm:
Public Types | |
Enums. | |
| 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 | |
Constructors/destructors. | |
| AztecOO_StatusTestResNorm (const Epetra_Operator &Operator, const Epetra_Vector &LHS, const Epetra_Vector &RHS, double Tolerance) | |
| Constructor. | |
| virtual | ~AztecOO_StatusTestResNorm () |
| Destructor. | |
Form and parameter definition methods. | |
| int | DefineResForm (ResType TypeOfResidual, NormType TypeOfNorm, Epetra_Vector *Weights=0) |
| Define form of the residual, its norm and optional weighting vector. | |
| 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. | |
| int | ResetTolerance (double Tolerance) |
| Reset the value of the tolerance. | |
Methods that implement the AztecOO_StatusTest interface. | |
| bool | ResidualVectorRequired () const |
| Indicates if residual vector is required by this convergence test. | |
| AztecOO_StatusType | CheckStatus (int CurrentIter, Epetra_MultiVector *CurrentResVector, double CurrentResNormEst, bool SolutionUpdated) |
| Check convergence status: Unconverged, Converged, Failed. | |
| AztecOO_StatusType | GetStatus () const |
| Return the result of the most recent checkStatus call. | |
| ostream & | Print (ostream &stream, int indent=0) const |
| Output formatted description of stopping test to output stream. | |
Methods to access data members. | |
| 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, . | |
Protected Member Functions | |
Internal Methods. | |
| double | ComputeNorm (const Epetra_Vector &vec, NormType typeofnorm) |
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
is the residual vector, implicitly or explicitly computed (determined by enum ResType),
is the residual norm determined by the enum NormType (1-norm, 2-norm or inf-norm),
is the scale factor that can be passed in as a precomputed double precision number, or can be selected from by the enum ScaleType (norm of RHS, norm of initial residual).
is the tolerance that is passed in as a double precision number to the constructor. The value of
can be reset using the ResetTolerance() method.
|
|
Select the norm type, where
|
|
|
Select how the residual vector is produced.
|
|
|
Select the scale type.
|
|
||||||||||||||||||||
|
Constructor.
The constructor takes a single argument specifying the tolerance (
|
|
||||||||||||||||||||
|
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.
Implements AztecOO_StatusTest. |
|
||||||||||||||||
|
Define form of the residual, its norm and optional weighting vector.
This method defines the form of
|
|
||||||||||||||||||||
|
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:
|
|
|
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. |
|
|
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. |
1.3.9.1