42 #ifndef BELOS_STATUS_TEST_GEN_RESNORM_H
43 #define BELOS_STATUS_TEST_GEN_RESNORM_H
78 template <
class ScalarType,
class MV,
class OP>
84 typedef Teuchos::ScalarTraits<ScalarType>
SCT;
167 int setQuorum(
int quorum) {quorum_ = quorum;
return(0);}
201 void print(std::ostream& os,
int indent = 0)
const;
213 Teuchos::RCP<MV>
getSolution() {
if (restype_==
Implicit) {
return Teuchos::null; }
else {
return curSoln_; } }
229 const std::vector<MagnitudeType>*
getTestValue()
const {
return(&testvector_);};
261 std::ostringstream oss;
262 oss <<
"Belos::StatusTestGenResNorm<>: " << resFormStr();
263 oss <<
", tol = " << tolerance_;
275 std::string resFormStr()
const
277 std::ostringstream oss;
279 oss << ((resnormtype_==
OneNorm) ?
"1-Norm" : (resnormtype_==
TwoNorm) ?
"2-Norm" :
"Inf-Norm");
280 oss << ((restype_==
Explicit) ?
" Exp" :
" Imp");
284 if (scaletype_!=
None)
291 oss <<
" (User Scale)";
294 oss << ((scalenormtype_==
OneNorm) ?
"1-Norm" : (resnormtype_==
TwoNorm) ?
"2-Norm" :
"Inf-Norm");
320 bool showMaxResNormOnly_;
338 std::vector<MagnitudeType> scalevector_;
341 std::vector<MagnitudeType> resvector_;
344 std::vector<MagnitudeType> testvector_;
347 std::vector<int> ind_;
350 Teuchos::RCP<MV> curSoln_;
362 std::vector<int> curLSIdx_;
371 bool firstcallCheckStatus_;
374 bool firstcallDefineResForm_;
377 bool firstcallDefineScaleForm_;
383 template <
class ScalarType,
class MV,
class OP>
386 : tolerance_(Tolerance),
388 showMaxResNormOnly_(showMaxResNormOnly),
399 firstcallCheckStatus_(true),
400 firstcallDefineResForm_(true),
401 firstcallDefineScaleForm_(true)
407 template <
class ScalarType,
class MV,
class OP>
411 template <
class ScalarType,
class MV,
class OP>
420 firstcallCheckStatus_ =
true;
421 curSoln_ = Teuchos::null;
424 template <
class ScalarType,
class MV,
class OP>
427 TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineResForm_==
false,
StatusTestError,
428 "StatusTestGenResNorm::defineResForm(): The residual form has already been defined.");
429 firstcallDefineResForm_ =
false;
431 restype_ = TypeOfResidual;
432 resnormtype_ = TypeOfNorm;
437 template <
class ScalarType,
class MV,
class OP>
441 TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineScaleForm_==
false,
StatusTestError,
442 "StatusTestGenResNorm::defineScaleForm(): The scaling type has already been defined.");
443 firstcallDefineScaleForm_ =
false;
445 scaletype_ = TypeOfScaling;
446 scalenormtype_ = TypeOfNorm;
447 scalevalue_ = ScaleValue;
452 template <
class ScalarType,
class MV,
class OP>
455 MagnitudeType zero = Teuchos::ScalarTraits<MagnitudeType>::zero();
458 if (firstcallCheckStatus_) {
459 StatusType status = firstCallCheckStatusSetup(iSolver);
474 curBlksz_ = (int)curLSIdx_.size();
476 for (
int i=0; i<curBlksz_; ++i) {
477 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
480 curNumRHS_ = validLS;
481 curSoln_ = Teuchos::null;
487 if (status_==
Passed) {
return status_; }
489 if (restype_==Implicit) {
495 std::vector<MagnitudeType> tmp_resvector( curBlksz_ );
497 if ( residMV != Teuchos::null ) {
498 tmp_resvector.resize( MVT::GetNumberVecs( *residMV ) );
499 MVT::MvNorm( *residMV, tmp_resvector, resnormtype_ );
500 typename std::vector<int>::iterator p = curLSIdx_.begin();
501 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
504 resvector_[*p] = tmp_resvector[i];
507 typename std::vector<int>::iterator p = curLSIdx_.begin();
508 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
511 resvector_[*p] = tmp_resvector[i];
515 else if (restype_==Explicit) {
521 Teuchos::RCP<MV> cur_res = MVT::Clone( *curSoln_, MVT::GetNumberVecs( *curSoln_ ) );
523 std::vector<MagnitudeType> tmp_resvector( MVT::GetNumberVecs( *cur_res ) );
524 MVT::MvNorm( *cur_res, tmp_resvector, resnormtype_ );
525 typename std::vector<int>::iterator p = curLSIdx_.begin();
526 for (
int i=0; p<curLSIdx_.end(); ++p, ++i) {
529 resvector_[*p] = tmp_resvector[i];
536 if ( scalevector_.size() > 0 ) {
537 typename std::vector<int>::iterator p = curLSIdx_.begin();
538 for (; p<curLSIdx_.end(); ++p) {
542 if ( scalevector_[ *p ] != zero ) {
544 testvector_[ *p ] = resvector_[ *p ] / scalevector_[ *p ] / scalevalue_;
546 testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
552 typename std::vector<int>::iterator p = curLSIdx_.begin();
553 for (; p<curLSIdx_.end(); ++p) {
556 testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
562 ind_.resize( curLSIdx_.size() );
563 typename std::vector<int>::iterator p = curLSIdx_.begin();
564 for (; p<curLSIdx_.end(); ++p) {
568 if (testvector_[ *p ] > tolerance_) {
570 }
else if (testvector_[ *p ] <= tolerance_) {
576 TEUCHOS_TEST_FOR_EXCEPTION(
true,
StatusTestError,
"StatusTestGenResNorm::checkStatus(): NaN has been detected.");
581 int need = (quorum_ == -1) ? curNumRHS_: quorum_;
588 template <
class ScalarType,
class MV,
class OP>
591 for (
int j = 0; j < indent; j ++)
593 printStatus(os, status_);
596 os <<
", tol = " << tolerance_ << std::endl;
599 if(showMaxResNormOnly_ && curBlksz_ > 1) {
601 testvector_.begin()+curLSIdx_[0],testvector_.begin()+curLSIdx_[curBlksz_-1]
603 for (
int j = 0; j < indent + 13; j ++)
605 os <<
"max{residual["<<curLSIdx_[0]<<
"..."<<curLSIdx_[curBlksz_-1]<<
"]} = " << maxRelRes
606 << ( maxRelRes <= tolerance_ ?
" <= " :
" > " ) << tolerance_ << std::endl;
609 for (
int i=0; i<numrhs_; i++ ) {
610 for (
int j = 0; j < indent + 13; j ++)
612 os <<
"residual [ " << i <<
" ] = " << testvector_[ i ];
613 os << ((testvector_[i]<tolerance_) ?
" < " : (testvector_[i]==tolerance_) ?
" == " : (testvector_[i]>tolerance_) ?
" > " :
" " ) << tolerance_ << std::endl;
620 template <
class ScalarType,
class MV,
class OP>
623 os << std::left << std::setw(13) << std::setfill(
'.');
636 os << std::left << std::setfill(
' ');
640 template <
class ScalarType,
class MV,
class OP>
644 MagnitudeType zero = Teuchos::ScalarTraits<MagnitudeType>::zero();
645 MagnitudeType one = Teuchos::ScalarTraits<MagnitudeType>::one();
648 if (firstcallCheckStatus_) {
652 firstcallCheckStatus_ =
false;
655 Teuchos::RCP<const MV> rhs = lp.
getRHS();
656 numrhs_ = MVT::GetNumberVecs( *rhs );
657 scalevector_.resize( numrhs_ );
658 MVT::MvNorm( *rhs, scalevector_, scalenormtype_ );
662 numrhs_ = MVT::GetNumberVecs( *init_res );
663 scalevector_.resize( numrhs_ );
664 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
668 numrhs_ = MVT::GetNumberVecs( *init_res );
669 scalevector_.resize( numrhs_ );
670 MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
673 numrhs_ = MVT::GetNumberVecs( *(lp.
getRHS()) );
676 resvector_.resize( numrhs_ );
677 testvector_.resize( numrhs_ );
681 curBlksz_ = (int)curLSIdx_.size();
683 for (i=0; i<curBlksz_; ++i) {
684 if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
687 curNumRHS_ = validLS;
690 for (i=0; i<numrhs_; i++) { testvector_[i] = one; }
693 if (scalevalue_ == zero) {