Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | List of all members
Belos::StatusTestImpResNorm< ScalarType, MV, OP > Class Template Reference

Convergence test using the implicit residual norm(s), with an explicit residual norm(s) check for loss of accuracy if necessary. More...

#include <BelosStatusTestImpResNorm.hpp>

Inheritance diagram for Belos::StatusTestImpResNorm< ScalarType, MV, OP >:
Inheritance graph
[legend]

Public Types

typedef Teuchos::ScalarTraits
< ScalarType >::magnitudeType 
MagnitudeType
 The type of the magnitude (absolute value) of a ScalarType. More...
 
- Public Types inherited from Belos::StatusTestResNorm< ScalarType, MV, OP >
typedef Teuchos::ScalarTraits
< ScalarType > 
SCT
 
typedef SCT::magnitudeType MagnitudeType
 
typedef MultiVecTraits
< ScalarType, MV > 
MVT
 

Abbreviations for method implementations

typedef Teuchos::ScalarTraits
< ScalarType > 
STS
 
typedef Teuchos::ScalarTraits
< MagnitudeType
STM
 
typedef MultiVecTraits
< ScalarType, MV > 
MVT
 

Constructors and destructor.

 StatusTestImpResNorm (MagnitudeType Tolerance, int quorum=-1, bool showMaxResNormOnly=false)
 Constructor. More...
 
virtual ~StatusTestImpResNorm ()
 Destructor (virtual for memory safety). More...
 

Form and parameter definition methods.

int defineResForm (NormType TypeOfNorm)
 Define form of the residual, its norm and optional weighting vector. More...
 
int defineScaleForm (ScaleType TypeOfScaling, NormType TypeOfNorm, MagnitudeType ScaleValue=Teuchos::ScalarTraits< MagnitudeType >::one())
 Define form of the scaling, its norm, its optional weighting vector, or, alternatively, define an explicit value. More...
 
int setTolerance (MagnitudeType tolerance)
 Set the value of the tolerance. More...
 
int setQuorum (int quorum)
 Sets the number of residuals that must pass the convergence test before Passed is returned. More...
 
int setShowMaxResNormOnly (bool showMaxResNormOnly)
 Set whether the only maximum residual norm is displayed when the print() method is called. More...
 

Status methods

StatusType checkStatus (Iteration< ScalarType, MV, OP > *iSolver)
 Check convergence status: Passed, Failed, or Undefined. More...
 
StatusType getStatus () const
 Return the result of the most recent CheckStatus call. More...
 

Reset methods

void reset ()
 Resets the internal configuration to the initial state. More...
 

Print methods

void print (std::ostream &os, int indent=0) const
 Output formatted description of stopping test to output stream. More...
 
void printStatus (std::ostream &os, StatusType type) const
 Print message for each status specific to this stopping test. More...
 

Methods to access data members.

Teuchos::RCP< MV > getSolution ()
 Returns the current solution estimate that was computed for the most recent residual test. More...
 
int getQuorum () const
 Returns the number of residuals that must pass the convergence test before Passed is returned. More...
 
bool getShowMaxResNormOnly ()
 Returns whether the only maximum residual norm is displayed when the print() method is called. More...
 
std::vector< int > convIndices ()
 Returns the vector containing the indices of the residuals that passed the test. More...
 
MagnitudeType getTolerance () const
 "Original" convergence tolerance $\tau$ as set by user. More...
 
MagnitudeType getCurrTolerance () const
 Current convergence tolerance; may be changed to prevent loss of accuracy. More...
 
const std::vector
< MagnitudeType > * 
getTestValue () const
 Returns the test value, $ \frac{\|r\|}{\sigma} $, computed in most recent call to CheckStatus. More...
 
const std::vector
< MagnitudeType > * 
getResNormValue () const
 Returns the residual norm value, $ \|r\| $, computed in most recent call to CheckStatus. More...
 
const std::vector
< MagnitudeType > * 
getScaledNormValue () const
 Returns the scaled norm value, $ \sigma $. More...
 
bool getLOADetected () const
 Returns a boolean indicating a loss of accuracy has been detected in computing the residual. More...
 

Misc.

StatusType firstCallCheckStatusSetup (Iteration< ScalarType, MV, OP > *iSolver)
 Call to setup initial scaling vector. More...
 

Overridden from Teuchos::Describable

std::string description () const
 Method to return description of the maximum iteration status test. More...
 

Private methods.

std::string resFormStr () const
 Description of current residual form. More...
 

Private data members.

MagnitudeType tolerance_
 Current tolerance used to determine convergence and the default tolerance set by user. More...
 
MagnitudeType currTolerance_
 
int quorum_
 Number of residuals that must pass the convergence test before Passed is returned. More...
 
bool showMaxResNormOnly_
 Determines if the entries for all of the residuals are shown or just the max. More...
 
NormType resnormtype_
 Type of norm to use on residual (OneNorm, TwoNorm, or InfNorm). More...
 
ScaleType scaletype_
 Type of scaling to use (Norm of RHS, Norm of Initial Residual, None or User provided) More...
 
NormType scalenormtype_
 Type of norm to use on the scaling (OneNorm, TwoNorm, or InfNorm) More...
 
MagnitudeType scalevalue_
 Scaling value. More...
 
std::vector< MagnitudeTypescalevector_
 Scaling vector. More...
 
std::vector< MagnitudeTyperesvector_
 Residual norm vector. More...
 
std::vector< MagnitudeTypetestvector_
 Test vector = resvector_ / scalevector_. More...
 
Teuchos::RCP< MV > curSoln_
 Current solution vector. More...
 
std::vector< int > ind_
 Vector containing the indices for the vectors that passed the test. More...
 
StatusType status_
 Status. More...
 
int curBlksz_
 The current blocksize of the linear system being solved. More...
 
int curNumRHS_
 The current number of right-hand sides being solved for. More...
 
std::vector< int > curLSIdx_
 The indices of the current number of right-hand sides being solved for. More...
 
int curLSNum_
 The current number of linear systems that have been loaded into the linear problem. More...
 
int numrhs_
 The total number of right-hand sides being solved for. More...
 
bool firstcallCheckStatus_
 Is this the first time CheckStatus is called? More...
 
bool firstcallDefineResForm_
 Is this the first time DefineResForm is called? More...
 
bool firstcallDefineScaleForm_
 Is this the first time DefineScaleForm is called? More...
 
bool lossDetected_
 Has a loss in accuracy been detected? More...
 

Additional Inherited Members

- Public Member Functions inherited from Belos::StatusTestResNorm< ScalarType, MV, OP >
virtual int setTolerance (MagnitudeType tolerance)=0
 Set the value of the tolerance. More...
 
virtual int defineScaleForm (ScaleType TypeOfScaling, NormType TypeOfNorm, MagnitudeType ScaleValue=Teuchos::ScalarTraits< MagnitudeType >::one())=0
 Define the form of the scaling for the residual. More...
 
- Public Member Functions inherited from Belos::StatusTest< ScalarType, MV, OP >
 StatusTest ()
 Constructor. More...
 
virtual ~StatusTest ()
 Destructor. More...
 
- Public Member Functions inherited from Teuchos::Describable
DescribableStreamManipulatorState describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default)
 
std::ostream & operator<< (std::ostream &os, const DescribableStreamManipulatorState &d)
 
virtual void describe (FancyOStream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
void describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
virtual ~Describable ()
 
 LabeledObject ()
 
virtual ~LabeledObject ()
 
virtual void setObjectLabel (const std::string &objectLabel)
 
virtual std::string getObjectLabel () const
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 

Detailed Description

template<class ScalarType, class MV, class OP>
class Belos::StatusTestImpResNorm< ScalarType, MV, OP >

Convergence test using the implicit residual norm(s), with an explicit residual norm(s) check for loss of accuracy if necessary.

This test passes when a quorum of residual vectors (for all right-hand sides) passes a residual norm test. For residual vector $r_i$, the form of the test is $\|r_i\| / \sigma_i \le \tau$, where $\tau$ is the convergence tolerance set in the constructor or by setTolerance(). The default quorum consists of all the vectors, but you can set the number of vectors that must pass before the whole test passes, either by the constructor argument or by calling setQuorum().

The default residual vector norm (the norm used in the numerator of the test, and applied to the current residual vectors) is the 2-norm, but you can change the norm type (1-norm, 2-norm, or infinity norm) using defineResForm(). The default scaling factor $\sigma_i$ for residual vector i is the 2-norm of the initial residual vector, but you can change this using defineScaleForm(). The norm used for the scaling factors may be different than the norm used for the current residual vectors. Furthermore, each residual vector may use a different scaling factor $\sigma_i$.

This test starts by using the "implicit" residual norm, otherwise called the "recursive" or "native" residual norm. It comes from the iterative method's projected problem, unlike the "explicit" residual norm, which comes from explicitly computing $R = B - A * X$. Once the implicit residual norm reaches the convergence tolerance, this test then checks the explicit residual norm to make sure that it has reached the convergence tolerance as well.

We say that there is a "potential loss of accuracy" when we first detect that the implicit residual norm(s) have met the desired original convergence tolerance, but the explicit residual norm(s) have not. We don't actually call this a "loss of accuracy" (in the sense of getLOADetected() returning true) unless we tried to remedy this situation, but couldn't fix it. Upon detecting a potential loss of accuracy, this test tells the solver to "iterate a few more steps" by making the "current" residual tolerance (the value of getCurrTolerance()) smaller, and forcing the implicit residual norm(s) to meet the current tolerance before moving on to test the explicit norm(s). If that doesn't work, this test declares a loss of accuracy.

Definition at line 104 of file BelosStatusTestImpResNorm.hpp.

Member Typedef Documentation

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<ScalarType>::magnitudeType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::MagnitudeType

The type of the magnitude (absolute value) of a ScalarType.

Definition at line 107 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<ScalarType> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::STS
private

Definition at line 112 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
typedef Teuchos::ScalarTraits<MagnitudeType> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::STM
private

Definition at line 113 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
typedef MultiVecTraits<ScalarType,MV> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::MVT
private

Definition at line 114 of file BelosStatusTestImpResNorm.hpp.

Constructor & Destructor Documentation

template<class ScalarType , class MV , class OP >
Belos::StatusTestImpResNorm< ScalarType, MV, OP >::StatusTestImpResNorm ( MagnitudeType  Tolerance,
int  quorum = -1,
bool  showMaxResNormOnly = false 
)

Constructor.

Parameters
Tolerance[in] Convergence tolerance $\tau$.
quorum[in] The number of vectors in the problem that must pass the convergence criteria before this StatusTest passes. -1 means that all residuals must pass.
showMaxResNormOnly[in] Whether only the maximum residual norm (of all vectors) is displayed when the print() method is called.

Definition at line 421 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
Belos::StatusTestImpResNorm< ScalarType, MV, OP >::~StatusTestImpResNorm ( )
virtual

Destructor (virtual for memory safety).

Definition at line 445 of file BelosStatusTestImpResNorm.hpp.

Member Function Documentation

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::defineResForm ( NormType  TypeOfNorm)

Define form of the residual, its norm and optional weighting vector.

This method defines the form of $\|r\|$. We specify:

  • The norm to be used on the residual (this may be different than the norm used in defineScaleForm()).

Definition at line 464 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::defineScaleForm ( ScaleType  TypeOfScaling,
NormType  TypeOfNorm,
MagnitudeType  ScaleValue = Teuchos::ScalarTraits<MagnitudeType>::one() 
)

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

Definition at line 476 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::setTolerance ( MagnitudeType  tolerance)
inline

Set 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.

Definition at line 179 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::setQuorum ( int  quorum)
inlinevirtual

Sets the number of residuals that must pass the convergence test before Passed is returned.

Note
If quorum=-1 then all residuals must pass the convergence test before Passed is returned.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 186 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::setShowMaxResNormOnly ( bool  showMaxResNormOnly)
inlinevirtual

Set whether the only maximum residual norm is displayed when the print() method is called.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 192 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
StatusType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::checkStatus ( Iteration< ScalarType, MV, OP > *  iSolver)
virtual

Check convergence status: Passed, Failed, or Undefined.

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.

Returns
StatusType: Passed, Failed, or Undefined.

Implements Belos::StatusTest< ScalarType, MV, OP >.

Definition at line 492 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
StatusType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getStatus ( ) const
inlinevirtual

Return the result of the most recent CheckStatus call.

Implements Belos::StatusTest< ScalarType, MV, OP >.

Definition at line 211 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
void Belos::StatusTestImpResNorm< ScalarType, MV, OP >::reset ( )
virtual

Resets the internal configuration to the initial state.

Implements Belos::StatusTest< ScalarType, MV, OP >.

Definition at line 449 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
void Belos::StatusTestImpResNorm< ScalarType, MV, OP >::print ( std::ostream &  os,
int  indent = 0 
) const
virtual

Output formatted description of stopping test to output stream.

Implements Belos::StatusTest< ScalarType, MV, OP >.

Definition at line 752 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
void Belos::StatusTestImpResNorm< ScalarType, MV, OP >::printStatus ( std::ostream &  os,
StatusType  type 
) const
virtual

Print message for each status specific to this stopping test.

Reimplemented from Belos::StatusTest< ScalarType, MV, OP >.

Definition at line 784 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getSolution ( )
inlinevirtual

Returns the current solution estimate that was computed for the most recent residual test.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 236 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getQuorum ( ) const
inlinevirtual

Returns the number of residuals that must pass the convergence test before Passed is returned.

Note
If quorum=-1 then all residuals must pass the convergence test before Passed is returned.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 240 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getShowMaxResNormOnly ( )
inlinevirtual

Returns whether the only maximum residual norm is displayed when the print() method is called.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 243 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::vector<int> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::convIndices ( )
inlinevirtual

Returns the vector containing the indices of the residuals that passed the test.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 246 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getTolerance ( ) const
inlinevirtual

"Original" convergence tolerance $\tau$ as set by user.

This value is the convergence tolerance as set by the user in the constructor, or by the setTolerance() method. See this class' main documentation and the documentation of getCurrTolerance() for an explanation of the difference between the "original" and "current" tolerances.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 255 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getCurrTolerance ( ) const
inline

Current convergence tolerance; may be changed to prevent loss of accuracy.

The difference between "original" tolerance (the value of getTolerance()) and "current" tolerance (the value of this method) relates to the idea of "loss of accuracy." See this class' main documentation for details. We do not allow users to set the "current" tolerance.

Definition at line 266 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
const std::vector<MagnitudeType>* Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getTestValue ( ) const
inlinevirtual

Returns the test value, $ \frac{\|r\|}{\sigma} $, computed in most recent call to CheckStatus.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 271 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
const std::vector<MagnitudeType>* Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getResNormValue ( ) const
inline

Returns the residual norm value, $ \|r\| $, computed in most recent call to CheckStatus.

Definition at line 274 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
const std::vector<MagnitudeType>* Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getScaledNormValue ( ) const
inline

Returns the scaled norm value, $ \sigma $.

Definition at line 277 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::StatusTestImpResNorm< ScalarType, MV, OP >::getLOADetected ( ) const
inlinevirtual

Returns a boolean indicating a loss of accuracy has been detected in computing the residual.

Implements Belos::StatusTestResNorm< ScalarType, MV, OP >.

Definition at line 280 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
StatusType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::firstCallCheckStatusSetup ( Iteration< ScalarType, MV, OP > *  iSolver)

Call to setup initial scaling vector.

After this function is called getScaledNormValue() can be called to get the scaling vector.

Definition at line 808 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::StatusTestImpResNorm< ScalarType, MV, OP >::description ( ) const
inlinevirtual

Method to return description of the maximum iteration status test.

Reimplemented from Teuchos::Describable.

Definition at line 299 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::string Belos::StatusTestImpResNorm< ScalarType, MV, OP >::resFormStr ( ) const
inlineprivate

Description of current residual form.

Definition at line 315 of file BelosStatusTestImpResNorm.hpp.

Member Data Documentation

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::tolerance_
private

Current tolerance used to determine convergence and the default tolerance set by user.

Definition at line 351 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::currTolerance_
private

Definition at line 351 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::quorum_
private

Number of residuals that must pass the convergence test before Passed is returned.

Definition at line 354 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::StatusTestImpResNorm< ScalarType, MV, OP >::showMaxResNormOnly_
private

Determines if the entries for all of the residuals are shown or just the max.

Definition at line 357 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
NormType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::resnormtype_
private

Type of norm to use on residual (OneNorm, TwoNorm, or InfNorm).

Definition at line 360 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
ScaleType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::scaletype_
private

Type of scaling to use (Norm of RHS, Norm of Initial Residual, None or User provided)

Definition at line 363 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
NormType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::scalenormtype_
private

Type of norm to use on the scaling (OneNorm, TwoNorm, or InfNorm)

Definition at line 366 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
MagnitudeType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::scalevalue_
private

Scaling value.

Definition at line 369 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::vector<MagnitudeType> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::scalevector_
private

Scaling vector.

Definition at line 372 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::vector<MagnitudeType> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::resvector_
private

Residual norm vector.

Definition at line 375 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::vector<MagnitudeType> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::testvector_
private

Test vector = resvector_ / scalevector_.

Definition at line 378 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
Teuchos::RCP<MV> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::curSoln_
private

Current solution vector.

Definition at line 381 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::vector<int> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::ind_
private

Vector containing the indices for the vectors that passed the test.

Definition at line 384 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
StatusType Belos::StatusTestImpResNorm< ScalarType, MV, OP >::status_
private

Status.

Definition at line 387 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::curBlksz_
private

The current blocksize of the linear system being solved.

Definition at line 390 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::curNumRHS_
private

The current number of right-hand sides being solved for.

Definition at line 393 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
std::vector<int> Belos::StatusTestImpResNorm< ScalarType, MV, OP >::curLSIdx_
private

The indices of the current number of right-hand sides being solved for.

Definition at line 396 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::curLSNum_
private

The current number of linear systems that have been loaded into the linear problem.

Definition at line 399 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
int Belos::StatusTestImpResNorm< ScalarType, MV, OP >::numrhs_
private

The total number of right-hand sides being solved for.

Definition at line 402 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::StatusTestImpResNorm< ScalarType, MV, OP >::firstcallCheckStatus_
private

Is this the first time CheckStatus is called?

Definition at line 405 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::StatusTestImpResNorm< ScalarType, MV, OP >::firstcallDefineResForm_
private

Is this the first time DefineResForm is called?

Definition at line 408 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::StatusTestImpResNorm< ScalarType, MV, OP >::firstcallDefineScaleForm_
private

Is this the first time DefineScaleForm is called?

Definition at line 411 of file BelosStatusTestImpResNorm.hpp.

template<class ScalarType , class MV , class OP >
bool Belos::StatusTestImpResNorm< ScalarType, MV, OP >::lossDetected_
private

Has a loss in accuracy been detected?

Definition at line 414 of file BelosStatusTestImpResNorm.hpp.


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