Belos
Version of the Day
|
A Belos::StatusTest class for specifying convergence of LSQR. The outer status tests passes if an inner status passes a user specified number of times consecutively. The inner status test depends on information specific to LSQR iteration. More...
#include <BelosLSQRStatusTest.hpp>
Public Types | |
typedef Teuchos::ScalarTraits < ScalarType > | SCT |
typedef SCT::magnitudeType | MagnitudeType |
typedef Belos::MultiVecTraits < ScalarType, MV > | MVT |
Public Member Functions | |
Constructor/Destructor. | |
LSQRStatusTest (MagnitudeType condMax=0.0, int term_iter_max=1, MagnitudeType rel_rhs_err=0.0, MagnitudeType rel_mat_err=0.0) | |
Constructor. More... | |
virtual | ~LSQRStatusTest () |
Destructor. More... | |
Status method | |
Belos::StatusType | checkStatus (Belos::Iteration< ScalarType, MV, OP > *iSolver) |
Check convergence status of the iterative solver: Unconverged, Converged, Failed. More... | |
Belos::StatusType | getStatus () const |
Return the result of the most recent CheckStatus call. More... | |
Reset methods | |
void | reset () |
Resets the status test to the initial internal state. More... | |
int | setCondLim (MagnitudeType condMax) |
Set the tolerances. More... | |
int | setTermIterMax (int term_iter_max) |
int | setRelRhsErr (MagnitudeType rel_rhs_err) |
int | setRelMatErr (MagnitudeType rel_mat_err) |
Accessor methods | |
MagnitudeType | getCondMaxLim () const |
Returns the value of the upper limit of the condition number of Abar set in the constructor. More... | |
int | getTermIterMax () const |
Returns the number of successful convergent iterations required set in the constructor. More... | |
MagnitudeType | getRelRhsErr () const |
Returns the value of the estimate of the relative error in the data defining b set in the constructor. More... | |
MagnitudeType | getMatErr () const |
Returns the value of the estimate of the relative error in the data defining A set in the constructor. More... | |
MagnitudeType | getMatCondNum () const |
Returns the value of the observed condition number of Abar. More... | |
MagnitudeType | getMatNorm () const |
Returns the value of the observed (Frobenius) norm of A. More... | |
int | getTermIter () const |
!Returns the current number of successful iterations from the most recent StatusTest call. More... | |
MagnitudeType | getResidNorm () const |
Returns the value of the observed norm of the residual r = b-Ax. More... | |
MagnitudeType | getLSResidNorm () const |
Returns the value of the observed norm of the Least Squares residual A^T r. 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, Belos::StatusType type) const |
Print message for each status specific to this stopping test. More... | |
Misc. | |
Belos::StatusType | firstCallCheckStatusSetup (Belos::Iteration< ScalarType, MV, OP > *iSolver) |
Called in checkStatus exactly once, on the first call to checkStatus. More... | |
Overridden from Teuchos::Describable | |
std::string | description () const |
Method to return description of the maximum iteration status test. 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 | |
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 |
DescribableStreamManipulatorState | describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default) |
std::ostream & | operator<< (std::ostream &os, const DescribableStreamManipulatorState &d) |
Additional Inherited Members | |
Static Public Attributes inherited from Teuchos::Describable | |
static const EVerbosityLevel | verbLevel_default |
A Belos::StatusTest class for specifying convergence of LSQR. The outer status tests passes if an inner status passes a user specified number of times consecutively. The inner status test depends on information specific to LSQR iteration.
Definition at line 33 of file BelosLSQRStatusTest.hpp.
typedef Teuchos::ScalarTraits<ScalarType> Belos::LSQRStatusTest< ScalarType, MV, OP >::SCT |
Definition at line 38 of file BelosLSQRStatusTest.hpp.
typedef SCT::magnitudeType Belos::LSQRStatusTest< ScalarType, MV, OP >::MagnitudeType |
Definition at line 39 of file BelosLSQRStatusTest.hpp.
typedef Belos::MultiVecTraits<ScalarType,MV> Belos::LSQRStatusTest< ScalarType, MV, OP >::MVT |
Definition at line 40 of file BelosLSQRStatusTest.hpp.
Belos::LSQRStatusTest< ScalarType, MV, OP >::LSQRStatusTest | ( | MagnitudeType | condMax = 0.0 , |
int | term_iter_max = 1 , |
||
MagnitudeType | rel_rhs_err = 0.0 , |
||
MagnitudeType | rel_mat_err = 0.0 |
||
) |
Constructor.
The constructor has four optional arguments, specifying the maximum observed condition number of Abar,
the number of successive convergent iterations, an estimate of the relative error in the data defining the right-hand side (b), and an estimate of the relative error in the data defining the coefficinet matrix (A). The default termIterMax is 1, and the other three parameters default to 0. The defaults specified in LSQRSolMgr are more realistic.
Definition at line 209 of file BelosLSQRStatusTest.hpp.
|
virtual |
Destructor.
Definition at line 227 of file BelosLSQRStatusTest.hpp.
|
virtual |
Check convergence status of the iterative solver: Unconverged, Converged, Failed.
This method checks if the convergence criteria are met using the current information from the iterative solver.
Implements Belos::StatusTest< ScalarType, MV, OP >.
Definition at line 237 of file BelosLSQRStatusTest.hpp.
|
inlinevirtual |
Return the result of the most recent CheckStatus call.
Implements Belos::StatusTest< ScalarType, MV, OP >.
Definition at line 69 of file BelosLSQRStatusTest.hpp.
|
virtual |
Resets the status test to the initial internal state.
Implements Belos::StatusTest< ScalarType, MV, OP >.
Definition at line 231 of file BelosLSQRStatusTest.hpp.
|
inline |
Set the tolerances.
Definition at line 80 of file BelosLSQRStatusTest.hpp.
|
inline |
Definition at line 85 of file BelosLSQRStatusTest.hpp.
|
inline |
Definition at line 91 of file BelosLSQRStatusTest.hpp.
|
inline |
Definition at line 95 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the value of the upper limit of the condition number of Abar set in the constructor.
Definition at line 105 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the number of successful convergent iterations required set in the constructor.
Definition at line 108 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the value of the estimate of the relative error in the data defining b set in the constructor.
Definition at line 111 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the value of the estimate of the relative error in the data defining A set in the constructor.
Definition at line 114 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the value of the observed condition number of Abar.
Definition at line 117 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the value of the observed (Frobenius) norm of A.
Definition at line 120 of file BelosLSQRStatusTest.hpp.
|
inline |
!Returns the current number of successful iterations from the most recent StatusTest call.
Definition at line 123 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the value of the observed norm of the residual r = b-Ax.
Definition at line 126 of file BelosLSQRStatusTest.hpp.
|
inline |
Returns the value of the observed norm of the Least Squares residual A^T r.
Definition at line 129 of file BelosLSQRStatusTest.hpp.
|
virtual |
Output formatted description of stopping test to output stream.
Implements Belos::StatusTest< ScalarType, MV, OP >.
Definition at line 363 of file BelosLSQRStatusTest.hpp.
|
virtual |
Print message for each status specific to this stopping test.
Reimplemented from Belos::StatusTest< ScalarType, MV, OP >.
Definition at line 373 of file BelosLSQRStatusTest.hpp.
Belos::StatusType Belos::LSQRStatusTest< ScalarType, MV, OP >::firstCallCheckStatusSetup | ( | Belos::Iteration< ScalarType, MV, OP > * | iSolver | ) |
Called in checkStatus exactly once, on the first call to checkStatus.
|
inlinevirtual |
Method to return description of the maximum iteration status test.
Reimplemented from Teuchos::Describable.
Definition at line 156 of file BelosLSQRStatusTest.hpp.