42 #ifndef BELOS_STATUS_TEST_GEN_RESNORM_H
43 #define BELOS_STATUS_TEST_GEN_RESNORM_H
78 template <
class ScalarType,
class MV,
class OP>
203 void print(std::ostream& os,
int indent = 0)
const;
263 std::ostringstream oss;
264 oss <<
"Belos::StatusTestGenResNorm<>: " <<
resFormStr();
279 std::ostringstream oss;
293 oss <<
" (User Scale)";
385 template <
class ScalarType,
class MV,
class OP>
388 : tolerance_(Tolerance),
390 showMaxResNormOnly_(showMaxResNormOnly),
401 firstcallCheckStatus_(true),
402 firstcallDefineResForm_(true),
403 firstcallDefineScaleForm_(true)
409 template <
class ScalarType,
class MV,
class OP>
413 template <
class ScalarType,
class MV,
class OP>
422 firstcallCheckStatus_ =
true;
426 template <
class ScalarType,
class MV,
class OP>
430 "StatusTestGenResNorm::defineResForm(): The residual form has already been defined.");
431 firstcallDefineResForm_ =
false;
433 restype_ = TypeOfResidual;
434 resnormtype_ = TypeOfNorm;
439 template <
class ScalarType,
class MV,
class OP>
444 "StatusTestGenResNorm::defineScaleForm(): The scaling type has already been defined.");
445 firstcallDefineScaleForm_ =
false;
447 scaletype_ = TypeOfScaling;
448 scalenormtype_ = TypeOfNorm;
449 scalevalue_ = ScaleValue;
454 template <
class ScalarType,
class MV,
class OP>
460 if (firstcallCheckStatus_) {
461 StatusType status = firstCallCheckStatusSetup(iSolver);
476 curBlksz_ = (int)curLSIdx_.size();
478 for (
int i=0; i<curBlksz_; ++i) {
479 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
482 curNumRHS_ = validLS;
489 if (status_==
Passed) {
return status_; }
491 if (restype_==Implicit) {
497 std::vector<MagnitudeType> tmp_resvector( curBlksz_ );
500 tmp_resvector.resize( MVT::GetNumberVecs( *residMV ) );
501 MVT::MvNorm( *residMV, tmp_resvector, resnormtype_ );
502 typename std::vector<int>::iterator p = curLSIdx_.begin();
503 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
506 resvector_[*p] = tmp_resvector[i];
509 typename std::vector<int>::iterator p = curLSIdx_.begin();
510 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
513 resvector_[*p] = tmp_resvector[i];
517 else if (restype_==Explicit) {
523 Teuchos::RCP<MV> cur_res = MVT::Clone( *curSoln_, MVT::GetNumberVecs( *curSoln_ ) );
525 std::vector<MagnitudeType> tmp_resvector( MVT::GetNumberVecs( *cur_res ) );
526 MVT::MvNorm( *cur_res, tmp_resvector, resnormtype_ );
527 typename std::vector<int>::iterator p = curLSIdx_.begin();
528 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
531 resvector_[*p] = tmp_resvector[i];
538 if ( scalevector_.size() > 0 ) {
539 typename std::vector<int>::iterator p = curLSIdx_.begin();
540 for (; p<curLSIdx_.end(); ++p) {
545 scalevector_[ *p ] != zero? resvector_[ *p ] / (scalevector_[ *p ] * scalevalue_) : resvector_[ *p ] / scalevalue_;
550 typename std::vector<int>::iterator p = curLSIdx_.begin();
551 for (; p<curLSIdx_.end(); ++p) {
554 testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
560 ind_.resize( curLSIdx_.size() );
561 typename std::vector<int>::iterator p = curLSIdx_.begin();
562 for (; p<curLSIdx_.end(); ++p) {
566 if (testvector_[ *p ] > tolerance_) {
568 }
else if (testvector_[ *p ] <= tolerance_) {
579 int need = (quorum_ == -1) ? curNumRHS_: quorum_;
586 template <
class ScalarType,
class MV,
class OP>
589 for (
int j = 0; j < indent; j ++)
591 printStatus(os, status_);
594 os <<
", tol = " << tolerance_ << std::endl;
597 if(showMaxResNormOnly_ && curBlksz_ > 1) {
599 testvector_.begin()+curLSIdx_[0],testvector_.begin()+curLSIdx_[curBlksz_-1]
601 for (
int j = 0; j < indent + 13; j ++)
603 os <<
"max{residual["<<curLSIdx_[0]<<
"..."<<curLSIdx_[curBlksz_-1]<<
"]} = " << maxRelRes
604 << ( maxRelRes <= tolerance_ ?
" <= " :
" > " ) << tolerance_ << std::endl;
607 for (
int i=0; i<numrhs_; i++ ) {
608 for (
int j = 0; j < indent + 13; j ++)
610 os <<
"residual [ " << i <<
" ] = " << testvector_[ i ];
611 os << ((testvector_[i]<tolerance_) ?
" < " : (testvector_[i]==tolerance_) ?
" == " : (testvector_[i]>tolerance_) ?
" > " :
" " ) << tolerance_ << std::endl;
618 template <
class ScalarType,
class MV,
class OP>
621 os << std::left << std::setw(13) << std::setfill(
'.');
634 os << std::left << std::setfill(
' ');
638 template <
class ScalarType,
class MV,
class OP>
646 if (firstcallCheckStatus_) {
650 firstcallCheckStatus_ =
false;
654 numrhs_ = MVT::GetNumberVecs( *rhs );
655 scalevector_.resize( numrhs_ );
656 MVT::MvNorm( *rhs, scalevector_, scalenormtype_ );
660 numrhs_ = MVT::GetNumberVecs( *init_res );
661 scalevector_.resize( numrhs_ );
662 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
666 numrhs_ = MVT::GetNumberVecs( *init_res );
667 scalevector_.resize( numrhs_ );
668 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
671 numrhs_ = MVT::GetNumberVecs( *(lp.
getRHS()) );
674 resvector_.resize( numrhs_ );
675 testvector_.resize( numrhs_ );
679 curBlksz_ = (int)curLSIdx_.size();
681 for (i=0; i<curBlksz_; ++i) {
682 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
685 curNumRHS_ = validLS;
688 for (i=0; i<numrhs_; i++) { testvector_[i] = one; }
691 if (scalevalue_ == zero) {
StatusType status_
Status.
virtual Teuchos::RCP< const MV > getNativeResiduals(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > *norms) const =0
Get the residuals native to the solver.
bool firstcallDefineResForm_
Is this the first time DefineResForm is called?
bool showMaxResNormOnly_
Determines if the entries for all of the residuals are shown or just the max.
std::string resFormStr() const
Description of current residual form.
int defineScaleForm(ScaleType TypeOfScaling, NormType TypeOfNorm, MagnitudeType ScaleValue=Teuchos::ScalarTraits< MagnitudeType >::one())
Define form of the scaling, its norm, its optional weighting std::vector, or, alternatively, define an explicit value.
MagnitudeType getTolerance() const
Returns the value of the tolerance, , set in the constructor.
const std::vector< MagnitudeType > * getTestValue() const
Returns the test value, , computed in most recent call to CheckStatus.
int curLSNum_
The current number of linear systems that have been loaded into the linear problem.
ScaleType
The type of scaling to use on the residual norm value.
virtual ~StatusTestGenResNorm()
Destructor.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
int setShowMaxResNormOnly(bool showMaxResNormOnly)
Set whether the only maximum residual norm is displayed when the print() method is called...
std::vector< int > curLSIdx_
The indices of the current number of right-hand sides being solved for.
StatusType firstCallCheckStatusSetup(Iteration< ScalarType, MV, OP > *iSolver)
Call to setup initial scaling std::vector.
virtual Teuchos::RCP< MV > updateSolution(const Teuchos::RCP< MV > &update=Teuchos::null, bool updateLP=false, ScalarType scale=Teuchos::ScalarTraits< ScalarType >::one())
Compute the new solution to the linear system using the given update vector.
ScaleType scaletype_
Type of scaling to use (Norm of RHS, Norm of Initial Residual, None or User provided) ...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual Teuchos::RCP< MV > getCurrentUpdate() const =0
Get the current update to the linear system.
Declaration of basic traits for the multivector type.
std::string description() const
Method to return description of the maximum iteration status test.
std::vector< MagnitudeType > testvector_
Test std::vector = resvector_ / scalevector_.
std::vector< int > ind_
Vector containing the indices for the vectors that passed the test.
Teuchos::RCP< const MV > getInitPrecResVec() const
A pointer to the preconditioned initial residual vector.
An abstract class of StatusTest for stopping criteria using residual norms.
void printStatus(std::ostream &os, StatusType type) const
Print message for each status specific to this stopping test.
An implementation of StatusTestResNorm using a family of residual norms.
int curBlksz_
The current blocksize of the linear system being solved.
Teuchos::RCP< const MV > getRHS() const
A pointer to the right-hand side B.
std::vector< MagnitudeType > resvector_
Residual norm std::vector.
StatusType
Whether the StatusTest wants iteration to stop.
const std::vector< MagnitudeType > * getResNormValue() const
Returns the residual norm value, , computed in most recent call to CheckStatus.
int numrhs_
The total number of right-hand sides being solved for.
Traits class which defines basic operations on multivectors.
SCT::magnitudeType MagnitudeType
NormType getResNormType()
const std::vector< int > & getLSIndex() const
(Zero-based) indices of the linear system(s) currently being solved.
MultiVecTraits< ScalarType, MV > MVT
NormType resnormtype_
Type of norm to use on residual (OneNorm, TwoNorm, or InfNorm).
bool getShowMaxResNormOnly()
Returns whether the only maximum residual norm is displayed when the print() method is called...
int defineResForm(ResType TypeOfResidual, NormType TypeOfNorm)
Define form of the residual, its norm and optional weighting std::vector.
bool firstcallCheckStatus_
Is this the first time CheckStatus is called?
Teuchos::ScalarTraits< ScalarType > SCT
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
int getLSNumber() const
The number of linear systems that have been set.
NormType scalenormtype_
Type of norm to use on the scaling (OneNorm, TwoNorm, or InfNorm)
int setQuorum(int quorum)
Sets the number of residuals that must pass the convergence test before Passed is returned...
std::vector< int > convIndices()
Returns the std::vector containing the indices of the residuals that passed the test.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const =0
Get a constant reference to the linear problem.
Teuchos::RCP< const MV > getInitResVec() const
A pointer to the initial unpreconditioned residual vector.
Belos::StatusTest abstract class for specifying a residual norm stopping criteria.
StatusType getStatus() const
Return the result of the most recent CheckStatus call.
int curNumRHS_
The current number of right-hand sides being solved for.
int setTolerance(MagnitudeType tolerance)
Set the value of the tolerance.
NormType
The type of vector norm to compute.
StatusType checkStatus(Iteration< ScalarType, MV, OP > *iSolver)
Check convergence status: Passed, Failed, or Undefined.
std::vector< MagnitudeType > scalevector_
Scaling std::vector.
virtual void computeCurrResVec(MV *R, const MV *X=0, const MV *B=0) const
Compute a residual R for this operator given a solution X, and right-hand side B. ...
ResType restype_
Type of residual to use (explicit or implicit)
void print(std::ostream &os, int indent=0) const
Output formatted description of stopping test to output stream.
ResType
Select how the residual std::vector is produced.
int quorum_
Number of residuals that must pass the convergence test before Passed is returned.
void reset()
Resets the internal configuration to the initial state.
bool firstcallDefineScaleForm_
Is this the first time DefineScaleForm is called?
int getQuorum() const
Returns the number of residuals that must pass the convergence test before Passed is returned...
MagnitudeType tolerance_
Tolerance used to determine convergence.
StatusTestGenResNorm(MagnitudeType Tolerance, int quorum=-1, bool showMaxResNormOnly=false)
Constructor.
MagnitudeType scalevalue_
Scaling value.
Teuchos::RCP< MV > curSoln_
Most recent solution vector used by this status test.
const std::vector< MagnitudeType > * getScaledNormValue() const
Returns the scaled norm value, .
Teuchos::RCP< MV > getSolution()
Returns the current solution estimate that was computed for the most recent residual test...
bool getLOADetected() const
Returns a boolean indicating a loss of accuracy has been detected in computing the residual...