AztecOO  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
AztecOO_StatusTestResNorm Class Reference

AztecOO_StatusTestResNorm: An implementation of AztecOO_StatusTest using a family of residual norms. More...

#include <AztecOO_StatusTestResNorm.h>

Inheritance diagram for AztecOO_StatusTestResNorm:
Inheritance graph
[legend]
Collaboration diagram for AztecOO_StatusTestResNorm:
Collaboration graph
[legend]

Public Types

enum  ResType { Implicit, Explicit }
 Select how the residual vector is produced. More...
 
enum  NormType { OneNorm, TwoNorm, InfNorm }
 Select the norm type, where $ x$ is the norm argument and $ w $ 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, $ \tau $, set in the constructor.
 
double GetTestValue () const
 Returns the test value, $ \frac{\|r\|}{\sigma} $, computed in most recent call to CheckStatus.
 
double GetResNormValue () const
 Returns the residual norm value, $ \|r\| $, computed in most recent call to CheckStatus.
 
double GetScaledNormValue () const
 Returns the scaled norm value, $ \sigma $.
 
- 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)
 

Detailed Description

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

\[ \frac{\|r\|}{\sigma} \le \tau \]

where

Warning
Presently it is not valid to associate one status test instance with two different AztecOO objects.

Member Enumeration Documentation

Select the norm type, where $ x$ is the norm argument and $ w $ is an optionally defined weighting vector.

Enumerator
OneNorm 

Use the one-norm $\sum_{i=1}^{n}(|x_i w_i|)$.

TwoNorm 

Use the two-norm $\sqrt(\sum_{i=1}^{n}((x_i w_i)^2)$.

InfNorm 

Use the inf-norm $(\max_{i=1}^{n}\{|x_i w_i|\})$.

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.

Select the scale type.

Enumerator
NormOfRHS 

Use the norm of the right-hand-side.

NormOfInitRes 

Use the initial residual vector (exlicitly computed).

None 

Use unscaled residual.

UserProvided 

User provides an explicit value that the norm of the residual will be divided by.

Constructor & Destructor Documentation

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 ( $\tau$). If none of the form definition methods are called, we use $\|r\|_2/\|r^{(0)}\|_2 \le \tau$ as the stopping criterion, where $\|r\|_2$ uses the least costly form of the 2-norm of residual available from the iterative method and $\|r^{(0)}\|_2$ 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.

Parameters
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().

Member Function Documentation

AztecOO_StatusType AztecOO_StatusTestResNorm::CheckStatus ( int  CurrentIter,
Epetra_MultiVector CurrentResVector,
double  CurrentResNormEst,
bool  SolutionUpdated 
)
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.

Parameters
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.
Returns
AztecOO_StatusType: Unconverged, Converged or Failed.

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 $\|r\|$. We specify:

  • Whether the residual vector should be explicitly computed, or taken from the iterative method.
  • The norm to be used on the residual (this may be different than the norm used in DefineScaleForm()).
  • A weighting vector that will be multiplied, element by element, with the residual vector prior to computing the norm. This argument defaults to a zero vector and no weighting will be done.

References Explicit, and TwoNorm.

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:

  1. User-provided scaling value:

    • Set argument TypeOfScaling to UserProvided.
    • Set ScaleValue to a non-zero value that the residual norm will be divided by.
    • TypeOfNorm and Weights arguments will be ignored.
    • Sample use: Define ScaleValue = $\|A\|_{\infty}$ where $ A $ is the matrix of the linear problem.

  2. Use a supported Scaling Form:
    • Define TypeOfScaling to be the norm of the right hand side, the initial residual vector, or to none.
    • Define norm to be used on the scaling vector (this may be different than the norm used in DefineResForm()).
    • Define a weighting vector that will be multiplied, element by element, with the residual vector prior to computing the norm. This argument defaults to a zero vector and no weighting will be done.

References NormOfInitRes, and TwoNorm.

int AztecOO_StatusTestResNorm::GetMaxNumExtraIterations ( )
inline

Return the maximum number of extra iterations that are performed if the implicit residual succeeds while the true residual fails.

int AztecOO_StatusTestResNorm::ResetTolerance ( double  Tolerance)
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.

bool AztecOO_StatusTestResNorm::ResidualVectorRequired ( ) const
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.

int AztecOO_StatusTestResNorm::SetMaxNumExtraIterations ( int  maxNumExtraIterations)
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.

Parameters
maxNumExtraIterations(In) Maximum number of extra iterations that will be performed in case the implicit residual succeeds that the true residual fails.

The documentation for this class was generated from the following files: