10#ifndef BELOS_STATUS_TEST_IMPRESNORM_H
11#define BELOS_STATUS_TEST_IMPRESNORM_H
21#include "Teuchos_as.hpp"
71template <
class ScalarType,
class MV,
class OP>
75 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
MagnitudeType;
80 typedef Teuchos::ScalarTraits<ScalarType> STS;
81 typedef Teuchos::ScalarTraits<MagnitudeType> STM;
102 bool showMaxResNormOnly =
false);
148 tolerance_ = tolerance;
161 showMaxResNormOnly_ = showMaxResNormOnly;
194 void print(std::ostream& os,
int indent = 0)
const;
235 return currTolerance_;
239 const std::vector<MagnitudeType>*
getTestValue()
const {
return(&testvector_);};
269 std::ostringstream oss;
270 oss <<
"Belos::StatusTestImpResNorm<>: " << resFormStr();
271 oss <<
", tol = " << tolerance_;
283 std::string resFormStr()
const
285 std::ostringstream oss;
287 oss << ((resnormtype_==
OneNorm) ?
"1-Norm" : (resnormtype_==
TwoNorm) ?
"2-Norm" :
"Inf-Norm");
291 if (scaletype_!=
None)
298 oss <<
" (User Scale)";
301 oss << ((scalenormtype_==
OneNorm) ?
"1-Norm" : (resnormtype_==
TwoNorm) ?
"2-Norm" :
"Inf-Norm");
325 bool showMaxResNormOnly_;
340 std::vector<MagnitudeType> scalevector_;
343 std::vector<MagnitudeType> resvector_;
346 std::vector<MagnitudeType> testvector_;
349 Teuchos::RCP<MV> curSoln_;
352 std::vector<int> ind_;
364 std::vector<int> curLSIdx_;
373 bool firstcallCheckStatus_;
376 bool firstcallDefineResForm_;
379 bool firstcallDefineScaleForm_;
387template <
class ScalarType,
class MV,
class OP>
390 : tolerance_(Tolerance),
391 currTolerance_(Tolerance),
393 showMaxResNormOnly_(showMaxResNormOnly),
403 firstcallCheckStatus_(true),
404 firstcallDefineResForm_(true),
405 firstcallDefineScaleForm_(true),
412template <
class ScalarType,
class MV,
class OP>
416template <
class ScalarType,
class MV,
class OP>
425 currTolerance_ = tolerance_;
426 firstcallCheckStatus_ =
true;
427 lossDetected_ =
false;
428 curSoln_ = Teuchos::null;
431template <
class ScalarType,
class MV,
class OP>
434 TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineResForm_==
false,
StatusTestError,
435 "StatusTestResNorm::defineResForm(): The residual form has already been defined.");
436 firstcallDefineResForm_ =
false;
438 resnormtype_ = TypeOfNorm;
443template <
class ScalarType,
class MV,
class OP>
447 TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineScaleForm_==
false,
StatusTestError,
448 "StatusTestResNorm::defineScaleForm(): The scaling type has already been defined.");
449 firstcallDefineScaleForm_ =
false;
451 scaletype_ = TypeOfScaling;
452 scalenormtype_ = TypeOfNorm;
453 scalevalue_ = ScaleValue;
458template <
class ScalarType,
class MV,
class OP>
469 if (firstcallCheckStatus_) {
486 curBlksz_ = (int)curLSIdx_.size();
488 for (
int i=0; i<curBlksz_; ++i) {
489 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
492 curNumRHS_ = validLS;
493 curSoln_ = Teuchos::null;
522 std::vector<MagnitudeType> tmp_resvector( curBlksz_ );
524 if (! residMV.is_null ()) {
527 MVT::MvNorm (*residMV, tmp_resvector, resnormtype_);
528 typename std::vector<int>::iterator p = curLSIdx_.begin();
529 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
532 resvector_[*p] = tmp_resvector[i];
536 typename std::vector<int>::iterator p = curLSIdx_.begin();
537 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
540 resvector_[*p] = tmp_resvector[i];
547 if (scalevector_.size () > 0) {
549 typename std::vector<int>::iterator p = curLSIdx_.begin();
550 for (; p<curLSIdx_.end(); ++p) {
554 if ( scalevector_[ *p ] != zero ) {
556 testvector_[ *p ] = resvector_[ *p ] / scalevector_[ *p ] / scalevalue_;
558 testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
564 typename std::vector<int>::iterator p = curLSIdx_.begin();
565 for (; p<curLSIdx_.end(); ++p) {
568 testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
581 ind_.resize( curLSIdx_.size() );
582 std::vector<int> lclInd( curLSIdx_.size() );
583 typename std::vector<int>::iterator p = curLSIdx_.begin();
584 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
587 if (testvector_[ *p ] > currTolerance_) {
589 }
else if (testvector_[ *p ] <= currTolerance_) {
602 "StatusTestImpResNorm::checkStatus(): One or more of the current "
603 "implicit residual norms is NaN.");
621 std::vector<MagnitudeType> tmp_testvector (have);
622 MVT::MvNorm (*cur_res, tmp_resvector, resnormtype_);
625 if ( scalevector_.size() > 0 ) {
626 for (
int i=0; i<have; ++i) {
628 if ( scalevector_[ ind_[i] ] != zero ) {
630 tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevector_[ ind_[i] ] / scalevalue_;
632 tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevalue_;
637 for (
int i=0; i<have; ++i) {
638 tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevalue_;
649 for (
int i = 0; i < have; ++i) {
659 if (tmp_testvector[i] <= tolerance_) {
660 ind_[have2] = ind_[i];
666 const MagnitudeType diff = STM::magnitude (testvector_[ind_[i]] - tmp_testvector[i]);
667 if (diff > currTolerance_) {
675 lossDetected_ =
true;
676 ind_[have2] = ind_[i];
696 const MagnitudeType onePointFive = as<MagnitudeType>(3) / as<MagnitudeType> (2);
697 const MagnitudeType oneTenth = STM::one () / as<MagnitudeType> (10);
699 currTolerance_ = currTolerance_ - onePointFive * diff;
700 while (currTolerance_ < STM::zero ()) {
701 currTolerance_ += oneTenth * diff;
712 int need = (quorum_ == -1) ? curNumRHS_: quorum_;
719template <
class ScalarType,
class MV,
class OP>
722 for (
int j = 0; j < indent; j ++)
727 os <<
", tol = " << tolerance_ << std::endl;
730 if(showMaxResNormOnly_ && curBlksz_ > 1) {
732 testvector_.begin()+curLSIdx_[0],testvector_.begin()+curLSIdx_[curBlksz_-1]
734 for (
int j = 0; j < indent + 13; j ++)
736 os <<
"max{residual["<<curLSIdx_[0]<<
"..."<<curLSIdx_[curBlksz_-1]<<
"]} = " << maxRelRes
737 << ( maxRelRes <= tolerance_ ?
" <= " :
" > " ) << tolerance_ << std::endl;
740 for (
int i=0; i<numrhs_; i++ ) {
741 for (
int j = 0; j < indent + 13; j ++)
743 os <<
"residual [ " << i <<
" ] = " << testvector_[ i ];
744 os << ((testvector_[i]<tolerance_) ?
" < " : (testvector_[i]==tolerance_) ?
" == " : (testvector_[i]>tolerance_) ?
" > " :
" " ) << tolerance_ << std::endl;
751template <
class ScalarType,
class MV,
class OP>
754 os << std::left << std::setw(13) << std::setfill(
'.');
761 os <<
"Unconverged (LoA)";
770 os << std::left << std::setfill(
' ');
774template <
class ScalarType,
class MV,
class OP>
783 if (firstcallCheckStatus_) {
787 firstcallCheckStatus_ =
false;
790 Teuchos::RCP<const MV> rhs = lp.
getRHS();
792 scalevector_.resize( numrhs_ );
798 scalevector_.resize( numrhs_ );
799 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
804 scalevector_.resize( numrhs_ );
805 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
811 resvector_.resize( numrhs_ );
812 testvector_.resize( numrhs_ );
816 curBlksz_ = (int)curLSIdx_.size();
818 for (i=0; i<curBlksz_; ++i) {
819 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
822 curNumRHS_ = validLS;
825 for (i=0; i<numrhs_; i++) { testvector_[i] = one; }
828 if (scalevalue_ == zero) {
Class which describes the linear problem to be solved by the iterative solver.
Declaration of basic traits for the multivector type.
Belos::StatusTest abstract class for specifying a residual norm stopping criteria.
virtual Teuchos::RCP< const MV > getNativeResiduals(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > *norms) const =0
virtual Teuchos::RCP< MV > getCurrentUpdate() const =0
Get the current update to the linear system.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const =0
Get a constant reference to the linear problem.
int getLSNumber() const
The number of linear systems that have been set.
Teuchos::RCP< const MV > getInitPrecResVec() const
A pointer to the preconditioned initial residual vector.
const std::vector< int > & getLSIndex() const
(Zero-based) indices of the linear system(s) currently being solved.
Teuchos::RCP< const MV > getRHS() const
A pointer to the right-hand side B.
Teuchos::RCP< const MV > getInitResVec() const
A pointer to the initial unpreconditioned residual 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.
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.
Traits class which defines basic operations on multivectors.
static void MvNorm(const MV &mv, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &normvec, NormType type=TwoNorm)
Compute the norm of each individual vector of mv. Upon return, normvec[i] holds the value of ,...
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
StatusTestImpResNorm(MagnitudeType Tolerance, int quorum=-1, bool showMaxResNormOnly=false)
Constructor.
const std::vector< MagnitudeType > * getScaledNormValue() const
Returns the scaled norm value, .
virtual ~StatusTestImpResNorm()
Destructor (virtual for memory safety).
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
The type of the magnitude (absolute value) of a ScalarType.
MagnitudeType getTolerance() const
"Original" convergence tolerance as set by user.
MagnitudeType getCurrTolerance() const
Current convergence tolerance; may be changed to prevent loss of accuracy.
std::string description() const
Method to return description of the maximum iteration status test.
std::vector< int > convIndices()
Returns the vector containing the indices of the residuals that passed the test.
void printStatus(std::ostream &os, StatusType type) const
Print message for each status specific to this stopping test.
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,...
const std::vector< MagnitudeType > * getResNormValue() const
Returns the residual norm value, , computed in most recent call to CheckStatus.
Teuchos::RCP< MV > getSolution()
Returns the current solution estimate that was computed for the most recent residual test.
void print(std::ostream &os, int indent=0) const
Output formatted description of stopping test to output stream.
void reset()
Resets the internal configuration to the initial state.
bool getShowMaxResNormOnly()
Returns whether the only maximum residual norm is displayed when the print() method is called.
bool getLOADetected() const
Returns a boolean indicating a loss of accuracy has been detected in computing the residual.
StatusType getStatus() const
Return the result of the most recent CheckStatus call.
StatusType checkStatus(Iteration< ScalarType, MV, OP > *iSolver)
Check convergence status: Passed, Failed, or Undefined.
int setTolerance(MagnitudeType tolerance)
Set the value of the tolerance.
const std::vector< MagnitudeType > * getTestValue() const
Returns the test value, , computed in most recent call to CheckStatus.
int defineResForm(NormType TypeOfNorm)
Define form of the residual, its norm and optional weighting vector.
int setShowMaxResNormOnly(bool showMaxResNormOnly)
Set whether the only maximum residual norm is displayed when the print() method is called.
StatusType firstCallCheckStatusSetup(Iteration< ScalarType, MV, OP > *iSolver)
Call to setup initial scaling vector.
int setQuorum(int quorum)
An abstract class of StatusTest for stopping criteria using residual norms.
NormType
The type of vector norm to compute.
StatusType
Whether the StatusTest wants iteration to stop.
ScaleType
The type of scaling to use on the residual norm value.