42 #ifndef BELOS_STATUS_TEST_IMPRESNORM_H
43 #define BELOS_STATUS_TEST_IMPRESNORM_H
53 #include "Teuchos_as.hpp"
103 template <
class ScalarType,
class MV,
class OP>
107 typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType
MagnitudeType;
112 typedef Teuchos::ScalarTraits<ScalarType> STS;
113 typedef Teuchos::ScalarTraits<MagnitudeType> STM;
134 bool showMaxResNormOnly =
false);
180 tolerance_ = tolerance;
193 showMaxResNormOnly_ = showMaxResNormOnly;
226 void print(std::ostream& os,
int indent = 0)
const;
267 return currTolerance_;
271 const std::vector<MagnitudeType>*
getTestValue()
const {
return(&testvector_);};
301 std::ostringstream oss;
302 oss <<
"Belos::StatusTestImpResNorm<>: " << resFormStr();
303 oss <<
", tol = " << tolerance_;
315 std::string resFormStr()
const
317 std::ostringstream oss;
319 oss << ((resnormtype_==
OneNorm) ?
"1-Norm" : (resnormtype_==
TwoNorm) ?
"2-Norm" :
"Inf-Norm");
323 if (scaletype_!=
None)
330 oss <<
" (User Scale)";
333 oss << ((scalenormtype_==
OneNorm) ?
"1-Norm" : (resnormtype_==
TwoNorm) ?
"2-Norm" :
"Inf-Norm");
357 bool showMaxResNormOnly_;
372 std::vector<MagnitudeType> scalevector_;
375 std::vector<MagnitudeType> resvector_;
378 std::vector<MagnitudeType> testvector_;
381 Teuchos::RCP<MV> curSoln_;
384 std::vector<int> ind_;
396 std::vector<int> curLSIdx_;
405 bool firstcallCheckStatus_;
408 bool firstcallDefineResForm_;
411 bool firstcallDefineScaleForm_;
419 template <
class ScalarType,
class MV,
class OP>
422 : tolerance_(Tolerance),
423 currTolerance_(Tolerance),
425 showMaxResNormOnly_(showMaxResNormOnly),
435 firstcallCheckStatus_(true),
436 firstcallDefineResForm_(true),
437 firstcallDefineScaleForm_(true),
444 template <
class ScalarType,
class MV,
class OP>
448 template <
class ScalarType,
class MV,
class OP>
457 currTolerance_ = tolerance_;
458 firstcallCheckStatus_ =
true;
459 lossDetected_ =
false;
460 curSoln_ = Teuchos::null;
463 template <
class ScalarType,
class MV,
class OP>
466 TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineResForm_==
false,
StatusTestError,
467 "StatusTestResNorm::defineResForm(): The residual form has already been defined.");
468 firstcallDefineResForm_ =
false;
470 resnormtype_ = TypeOfNorm;
475 template <
class ScalarType,
class MV,
class OP>
479 TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineScaleForm_==
false,
StatusTestError,
480 "StatusTestResNorm::defineScaleForm(): The scaling type has already been defined.");
481 firstcallDefineScaleForm_ =
false;
483 scaletype_ = TypeOfScaling;
484 scalenormtype_ = TypeOfNorm;
485 scalevalue_ = ScaleValue;
490 template <
class ScalarType,
class MV,
class OP>
501 if (firstcallCheckStatus_) {
502 StatusType status = firstCallCheckStatusSetup (iSolver);
518 curBlksz_ = (int)curLSIdx_.size();
520 for (
int i=0; i<curBlksz_; ++i) {
521 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
524 curNumRHS_ = validLS;
525 curSoln_ = Teuchos::null;
554 std::vector<MagnitudeType> tmp_resvector( curBlksz_ );
556 if (! residMV.is_null ()) {
558 tmp_resvector.resize (MVT::GetNumberVecs (*residMV));
559 MVT::MvNorm (*residMV, tmp_resvector, resnormtype_);
560 typename std::vector<int>::iterator p = curLSIdx_.begin();
561 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
564 resvector_[*p] = tmp_resvector[i];
568 typename std::vector<int>::iterator p = curLSIdx_.begin();
569 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
572 resvector_[*p] = tmp_resvector[i];
579 if (scalevector_.size () > 0) {
581 typename std::vector<int>::iterator p = curLSIdx_.begin();
582 for (; p<curLSIdx_.end(); ++p) {
586 if ( scalevector_[ *p ] != zero ) {
588 testvector_[ *p ] = resvector_[ *p ] / scalevector_[ *p ] / scalevalue_;
590 testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
596 typename std::vector<int>::iterator p = curLSIdx_.begin();
597 for (; p<curLSIdx_.end(); ++p) {
600 testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
613 ind_.resize( curLSIdx_.size() );
614 std::vector<int> lclInd( curLSIdx_.size() );
615 typename std::vector<int>::iterator p = curLSIdx_.begin();
616 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
619 if (testvector_[ *p ] > currTolerance_) {
621 }
else if (testvector_[ *p ] <= currTolerance_) {
634 "StatusTestImpResNorm::checkStatus(): One or more of the current "
635 "implicit residual norms is NaN.");
650 RCP<MV> cur_res = MVT::Clone (*curSoln_, MVT::GetNumberVecs (*curSoln_));
652 tmp_resvector.resize (MVT::GetNumberVecs (*cur_res));
653 std::vector<MagnitudeType> tmp_testvector (have);
654 MVT::MvNorm (*cur_res, tmp_resvector, resnormtype_);
657 if ( scalevector_.size() > 0 ) {
658 for (
int i=0; i<have; ++i) {
660 if ( scalevector_[ ind_[i] ] != zero ) {
662 tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevector_[ ind_[i] ] / scalevalue_;
664 tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevalue_;
669 for (
int i=0; i<have; ++i) {
670 tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevalue_;
681 for (
int i = 0; i < have; ++i) {
691 if (tmp_testvector[i] <= tolerance_) {
692 ind_[have2] = ind_[i];
698 const MagnitudeType diff = STM::magnitude (testvector_[ind_[i]] - tmp_testvector[i]);
699 if (diff > currTolerance_) {
707 lossDetected_ =
true;
708 ind_[have2] = ind_[i];
728 const MagnitudeType onePointFive = as<MagnitudeType>(3) / as<MagnitudeType> (2);
729 const MagnitudeType oneTenth = STM::one () / as<MagnitudeType> (10);
731 currTolerance_ = currTolerance_ - onePointFive * diff;
732 while (currTolerance_ < STM::zero ()) {
733 currTolerance_ += oneTenth * diff;
744 int need = (quorum_ == -1) ? curNumRHS_: quorum_;
751 template <
class ScalarType,
class MV,
class OP>
754 for (
int j = 0; j < indent; j ++)
756 printStatus(os, status_);
759 os <<
", tol = " << tolerance_ << std::endl;
762 if(showMaxResNormOnly_ && curBlksz_ > 1) {
764 testvector_.begin()+curLSIdx_[0],testvector_.begin()+curLSIdx_[curBlksz_-1]
766 for (
int j = 0; j < indent + 13; j ++)
768 os <<
"max{residual["<<curLSIdx_[0]<<
"..."<<curLSIdx_[curBlksz_-1]<<
"]} = " << maxRelRes
769 << ( maxRelRes <= tolerance_ ?
" <= " :
" > " ) << tolerance_ << std::endl;
772 for (
int i=0; i<numrhs_; i++ ) {
773 for (
int j = 0; j < indent + 13; j ++)
775 os <<
"residual [ " << i <<
" ] = " << testvector_[ i ];
776 os << ((testvector_[i]<tolerance_) ?
" < " : (testvector_[i]==tolerance_) ?
" == " : (testvector_[i]>tolerance_) ?
" > " :
" " ) << tolerance_ << std::endl;
783 template <
class ScalarType,
class MV,
class OP>
786 os << std::left << std::setw(13) << std::setfill(
'.');
793 os <<
"Unconverged (LoA)";
802 os << std::left << std::setfill(
' ');
806 template <
class ScalarType,
class MV,
class OP>
815 if (firstcallCheckStatus_) {
819 firstcallCheckStatus_ =
false;
822 Teuchos::RCP<const MV> rhs = lp.
getRHS();
823 numrhs_ = MVT::GetNumberVecs( *rhs );
824 scalevector_.resize( numrhs_ );
825 MVT::MvNorm( *rhs, scalevector_, scalenormtype_ );
829 numrhs_ = MVT::GetNumberVecs( *init_res );
830 scalevector_.resize( numrhs_ );
831 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
835 numrhs_ = MVT::GetNumberVecs( *init_res );
836 scalevector_.resize( numrhs_ );
837 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
840 numrhs_ = MVT::GetNumberVecs( *(lp.
getRHS()) );
843 resvector_.resize( numrhs_ );
844 testvector_.resize( numrhs_ );
848 curBlksz_ = (int)curLSIdx_.size();
850 for (i=0; i<curBlksz_; ++i) {
851 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
854 curNumRHS_ = validLS;
857 for (i=0; i<numrhs_; i++) { testvector_[i] = one; }
860 if (scalevalue_ == zero) {