Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosStatusTestGenResSubNorm.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef BELOS_STATUS_TEST_GEN_RESSUBNORM_H
11 #define BELOS_STATUS_TEST_GEN_RESSUBNORM_H
12 
19 #include "BelosLinearProblem.hpp"
20 #include "BelosMultiVecTraits.hpp"
21 #include "BelosOperatorTraits.hpp"
22 
23 #ifdef HAVE_BELOS_THYRA
24 #include <Thyra_MultiVectorBase.hpp>
25 #include <Thyra_MultiVectorStdOps.hpp>
26 #include <Thyra_ProductMultiVectorBase.hpp>
27 #endif
28 
37 namespace Belos {
38 
39 template <class ScalarType, class MV, class OP>
40 class StatusTestGenResSubNorm: public StatusTestResNorm<ScalarType,MV,OP> {
41 
42  public:
43  // Convenience typedefs
47 
49 
50 
64  StatusTestGenResSubNorm( MagnitudeType /* Tolerance */, size_t /* subIdx */, int /* quorum */ = -1, bool /* showMaxResNormOnly */ = false ) {
66  "StatusTestGenResSubNorm::StatusTestGenResSubNorm(): StatusTestGenResSubNorm only available for blocked operators (e.g., Thyra).");
67  }
68 
70  virtual ~StatusTestGenResSubNorm() { };
72 
74 
75 
77 
83  int defineResForm( NormType /* TypeOfNorm */) {
85  "StatusTestGenResSubNorm::defineResForm(): StatusTestGenResSubNorm only available for blocked operators (e.g., Thyra).");
87  }
88 
90 
110  int defineScaleForm( ScaleType /* TypeOfScaling */, NormType /* TypeOfNorm */, MagnitudeType /* ScaleValue */ = Teuchos::ScalarTraits<MagnitudeType>::one()) {
112  "StatusTestGenResSubNorm::defineScaleForm(): StatusTestGenResSubNorm only available for blocked operators (e.g., Thyra).");
114  }
115 
117 
120  int setTolerance(MagnitudeType /* tolerance */) { return 0; }
121 
123 
125  int setSubIdx ( size_t subIdx ) { return 0;}
126 
129  int setQuorum(int /* quorum */) { return 0; }
130 
132  int setShowMaxResNormOnly(bool /* showMaxResNormOnly */) { return 0; }
133 
135 
137 
138 
146 
148  StatusType getStatus() const {return Undefined;}
150 
152 
153 
155  void reset() { }
156 
158 
160 
161 
163  void print(std::ostream& /* os */, int /* indent */ = 0) const { }
164 
166  void printStatus(std::ostream& /* os */, StatusType /* type */) const { }
168 
170 
171 
175  Teuchos::RCP<MV> getSolution() { return Teuchos::null; }
176 
179  int getQuorum() const { return -1; }
180 
182  size_t getSubIdx() const { return 0; }
183 
185  bool getShowMaxResNormOnly() { return false; }
186 
188  std::vector<int> convIndices() { return std::vector<int>(0); }
189 
192 
194  const std::vector<MagnitudeType>* getTestValue() const {return NULL;};
195 
197  const std::vector<MagnitudeType>* getResNormValue() const {return NULL;};
198 
200  const std::vector<MagnitudeType>* getScaledNormValue() const {return NULL;};
201 
204  bool getLOADetected() const { return false; }
205 
207 
208 
211 
218  return Undefined;
219  }
221 
224 
226  std::string description() const
227  { return std::string(""); }
229 };
230 
231 #ifdef HAVE_BELOS_THYRA
232 
233 // specialization for Thyra
234 template <class ScalarType>
235 class StatusTestGenResSubNorm<ScalarType,Thyra::MultiVectorBase<ScalarType>,Thyra::LinearOpBase<ScalarType> >
236  : public StatusTestResNorm<ScalarType,Thyra::MultiVectorBase<ScalarType>,Thyra::LinearOpBase<ScalarType> > {
237 
238  public:
239  // Convenience typedefs
240  typedef Thyra::MultiVectorBase<ScalarType> MV;
241  typedef Thyra::LinearOpBase<ScalarType> OP;
242 
244  typedef typename SCT::magnitudeType MagnitudeType;
245  typedef MultiVecTraits<ScalarType,MV> MVT;
246  typedef OperatorTraits<ScalarType,MV,OP> OT;
247 
249 
250 
263  StatusTestGenResSubNorm( MagnitudeType Tolerance, size_t subIdx, int quorum = -1, bool showMaxResNormOnly = false )
264  : tolerance_(Tolerance),
265  subIdx_(subIdx),
266  quorum_(quorum),
267  showMaxResNormOnly_(showMaxResNormOnly),
268  resnormtype_(TwoNorm),
269  scaletype_(NormOfInitRes),
270  scalenormtype_(TwoNorm),
271  scalevalue_(Teuchos::ScalarTraits<MagnitudeType>::one ()),
272  status_(Undefined),
273  curBlksz_(0),
274  curNumRHS_(0),
275  curLSNum_(0),
276  numrhs_(0),
277  firstcallCheckStatus_(true),
278  firstcallDefineResForm_(true),
279  firstcallDefineScaleForm_(true) { }
280 
282  virtual ~StatusTestGenResSubNorm() { };
284 
286 
287 
289 
295  int defineResForm(NormType TypeOfNorm) {
296  TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineResForm_==false,StatusTestError,
297  "StatusTestGenResSubNorm::defineResForm(): The residual form has already been defined.");
298  firstcallDefineResForm_ = false;
299 
300  resnormtype_ = TypeOfNorm;
301 
302  return(0);
303  }
304 
306 
326  int defineScaleForm( ScaleType TypeOfScaling, NormType TypeOfNorm, MagnitudeType ScaleValue = Teuchos::ScalarTraits<MagnitudeType>::one()) {
327  TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineScaleForm_==false,StatusTestError,
328  "StatusTestGenResSubNorm::defineScaleForm(): The scaling type has already been defined.");
329  firstcallDefineScaleForm_ = false;
330 
331  scaletype_ = TypeOfScaling;
332  scalenormtype_ = TypeOfNorm;
333  scalevalue_ = ScaleValue;
334 
335  return(0);
336  }
337 
339 
342  int setTolerance(MagnitudeType tolerance) {tolerance_ = tolerance; return(0);}
343 
345 
347  int setSubIdx ( size_t subIdx ) { subIdx_ = subIdx; return(0);}
348 
351  int setQuorum(int quorum) {quorum_ = quorum; return(0);}
352 
354  int setShowMaxResNormOnly(bool showMaxResNormOnly) {showMaxResNormOnly_ = showMaxResNormOnly; return(0);}
355 
357 
359 
360 
367  StatusType checkStatus(Iteration<ScalarType,MV,OP>* iSolver) {
369  const LinearProblem<ScalarType,MV,OP>& lp = iSolver->getProblem();
370  // Compute scaling term (done once for each block that's being solved)
371  if (firstcallCheckStatus_) {
372  StatusType status = firstCallCheckStatusSetup(iSolver);
373  if(status==Failed) {
374  status_ = Failed;
375  return(status_);
376  }
377  }
378 
379  //
380  // This section computes the norm of the residual std::vector
381  //
382  if ( curLSNum_ != lp.getLSNumber() ) {
383  //
384  // We have moved on to the next rhs block
385  //
386  curLSNum_ = lp.getLSNumber();
387  curLSIdx_ = lp.getLSIndex();
388  curBlksz_ = (int)curLSIdx_.size();
389  int validLS = 0;
390  for (int i=0; i<curBlksz_; ++i) {
391  if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
392  validLS++;
393  }
394  curNumRHS_ = validLS;
395  curSoln_ = Teuchos::null;
396  //
397  } else {
398  //
399  // We are in the same rhs block, return if we are converged
400  //
401  if (status_==Passed) { return status_; }
402  }
403 
404  //
405  // Request the true residual for this block of right-hand sides.
406  //
407  Teuchos::RCP<MV> cur_update = iSolver->getCurrentUpdate();
408  curSoln_ = lp.updateSolution( cur_update );
409  Teuchos::RCP<MV> cur_res = MVT::Clone( *curSoln_, MVT::GetNumberVecs( *curSoln_ ) );
410  lp.computeCurrResVec( &*cur_res, &*curSoln_ );
411  std::vector<MagnitudeType> tmp_resvector( MVT::GetNumberVecs( *cur_res ) );
412  MvSubNorm( *cur_res, subIdx_, tmp_resvector, resnormtype_ );
413 
414  typename std::vector<int>::iterator pp = curLSIdx_.begin();
415  for (int i=0; pp<curLSIdx_.end(); ++pp, ++i) {
416  // Check if this index is valid
417  if (*pp != -1)
418  resvector_[*pp] = tmp_resvector[i];
419  }
420 
421  //
422  // Compute the new linear system residuals for testing.
423  // (if any of them don't meet the tolerance or are NaN, then we exit with that status)
424  //
425  if ( scalevector_.size() > 0 ) {
426  typename std::vector<int>::iterator p = curLSIdx_.begin();
427  for (; p<curLSIdx_.end(); ++p) {
428  // Check if this index is valid
429  if (*p != -1) {
430  // Scale the std::vector accordingly
431  if ( scalevector_[ *p ] != zero ) {
432  // Don't intentionally divide by zero.
433  testvector_[ *p ] = resvector_[ *p ] / scalevector_[ *p ] / scalevalue_;
434  } else {
435  testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
436  }
437  }
438  }
439  }
440  else {
441  typename std::vector<int>::iterator ppp = curLSIdx_.begin();
442  for (; ppp<curLSIdx_.end(); ++ppp) {
443  // Check if this index is valid
444  if (*ppp != -1)
445  testvector_[ *ppp ] = resvector_[ *ppp ] / scalevalue_;
446  }
447  }
448  // Check status of new linear system residuals and see if we have the quorum.
449  int have = 0;
450  ind_.resize( curLSIdx_.size() );
451  typename std::vector<int>::iterator p2 = curLSIdx_.begin();
452  for (; p2<curLSIdx_.end(); ++p2) {
453  // Check if this index is valid
454  if (*p2 != -1) {
455  // Check if any of the residuals are larger than the tolerance.
456  if (testvector_[ *p2 ] > tolerance_) {
457  // do nothing.
458  } else if (testvector_[ *p2 ] <= tolerance_) {
459  ind_[have] = *p2;
460  have++;
461  } else {
462  // Throw an std::exception if a NaN is found.
463  status_ = Failed;
464  TEUCHOS_TEST_FOR_EXCEPTION(true,StatusTestNaNError,"StatusTestGenResSubNorm::checkStatus(): NaN has been detected.");
465  }
466  }
467  }
468  ind_.resize(have);
469  int need = (quorum_ == -1) ? curNumRHS_: quorum_;
470  status_ = (have >= need) ? Passed : Failed;
471  // Return the current status
472  return status_;
473  }
474 
476  StatusType getStatus() const {return(status_);};
478 
480 
481 
483  void reset() {
484  status_ = Undefined;
485  curBlksz_ = 0;
486  curLSNum_ = 0;
487  curLSIdx_.resize(0);
488  numrhs_ = 0;
489  ind_.resize(0);
490  firstcallCheckStatus_ = true;
491  curSoln_ = Teuchos::null;
492  }
493 
495 
497 
498 
500  void print(std::ostream& os, int indent = 0) const {
501  os.setf(std::ios_base::scientific);
502  for (int j = 0; j < indent; j ++)
503  os << ' ';
504  printStatus(os, status_);
505  os << resFormStr();
506  if (status_==Undefined)
507  os << ", tol = " << tolerance_ << std::endl;
508  else {
509  os << std::endl;
510  if(showMaxResNormOnly_ && curBlksz_ > 1) {
511  const MagnitudeType maxRelRes = *std::max_element(
512  testvector_.begin()+curLSIdx_[0],testvector_.begin()+curLSIdx_[curBlksz_-1]
513  );
514  for (int j = 0; j < indent + 13; j ++)
515  os << ' ';
516  os << "max{residual["<<curLSIdx_[0]<<"..."<<curLSIdx_[curBlksz_-1]<<"]} = " << maxRelRes
517  << ( maxRelRes <= tolerance_ ? " <= " : " > " ) << tolerance_ << std::endl;
518  }
519  else {
520  for ( int i=0; i<numrhs_; i++ ) {
521  for (int j = 0; j < indent + 13; j ++)
522  os << ' ';
523  os << "residual [ " << i << " ] = " << testvector_[ i ];
524  os << ((testvector_[i]<tolerance_) ? " < " : (testvector_[i]==tolerance_) ? " == " : (testvector_[i]>tolerance_) ? " > " : " " ) << tolerance_ << std::endl;
525  }
526  }
527  }
528  os << std::endl;
529  }
530 
532  void printStatus(std::ostream& os, StatusType type) const {
533  os << std::left << std::setw(13) << std::setfill('.');
534  switch (type) {
535  case Passed:
536  os << "Converged";
537  break;
538  case Failed:
539  os << "Unconverged";
540  break;
541  case Undefined:
542  default:
543  os << "**";
544  break;
545  }
546  os << std::left << std::setfill(' ');
547  return;
548  }
550 
552 
553 
555  Teuchos::RCP<MV> getSolution() { return curSoln_; }
556 
559  int getQuorum() const { return quorum_; }
560 
562  size_t getSubIdx() const { return subIdx_; }
563 
565  bool getShowMaxResNormOnly() { return showMaxResNormOnly_; }
566 
568  std::vector<int> convIndices() { return ind_; }
569 
571  MagnitudeType getTolerance() const {return(tolerance_);};
572 
574  const std::vector<MagnitudeType>* getTestValue() const {return(&testvector_);};
575 
577  const std::vector<MagnitudeType>* getResNormValue() const {return(&resvector_);};
578 
580  const std::vector<MagnitudeType>* getScaledNormValue() const {return(&scalevector_);};
581 
584  bool getLOADetected() const { return false; }
585 
587 
588 
591 
597  StatusType firstCallCheckStatusSetup(Iteration<ScalarType,MV,OP>* iSolver) {
598  int i;
601  const LinearProblem<ScalarType,MV,OP>& lp = iSolver->getProblem();
602  // Compute scaling term (done once for each block that's being solved)
603  if (firstcallCheckStatus_) {
604  //
605  // Get some current solver information.
606  //
607  firstcallCheckStatus_ = false;
608 
609  if (scaletype_== NormOfRHS) {
610  Teuchos::RCP<const MV> rhs = lp.getRHS();
611  numrhs_ = MVT::GetNumberVecs( *rhs );
612  scalevector_.resize( numrhs_ );
613  MvSubNorm( *rhs, subIdx_, scalevector_, scalenormtype_ );
614  }
615  else if (scaletype_==NormOfInitRes) {
616  Teuchos::RCP<const MV> init_res = lp.getInitResVec();
617  numrhs_ = MVT::GetNumberVecs( *init_res );
618  scalevector_.resize( numrhs_ );
619  MvSubNorm( *init_res, subIdx_, scalevector_, scalenormtype_ );
620  }
621  else if (scaletype_==NormOfPrecInitRes) {
622  Teuchos::RCP<const MV> init_res = lp.getInitPrecResVec();
623  numrhs_ = MVT::GetNumberVecs( *init_res );
624  scalevector_.resize( numrhs_ );
625  MvSubNorm( *init_res, subIdx_, scalevector_, scalenormtype_ );
626  }
627  else if (scaletype_==NormOfFullInitRes) {
628  Teuchos::RCP<const MV> init_res = lp.getInitResVec();
629  numrhs_ = MVT::GetNumberVecs( *init_res );
630  scalevector_.resize( numrhs_ );
631  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
633  }
634  else if (scaletype_==NormOfFullPrecInitRes) {
635  Teuchos::RCP<const MV> init_res = lp.getInitPrecResVec();
636  numrhs_ = MVT::GetNumberVecs( *init_res );
637  scalevector_.resize( numrhs_ );
638  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
640  }
641  else if (scaletype_==NormOfFullScaledInitRes) {
642  Teuchos::RCP<const MV> init_res = lp.getInitResVec();
643  numrhs_ = MVT::GetNumberVecs( *init_res );
644  scalevector_.resize( numrhs_ );
645  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
646  MvScalingRatio( *init_res, subIdx_, scalevalue_ );
647  }
648  else if (scaletype_==NormOfFullScaledPrecInitRes) {
649  Teuchos::RCP<const MV> init_res = lp.getInitPrecResVec();
650  numrhs_ = MVT::GetNumberVecs( *init_res );
651  scalevector_.resize( numrhs_ );
652  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
653  MvScalingRatio( *init_res, subIdx_, scalevalue_ );
654  }
655  else {
656  numrhs_ = MVT::GetNumberVecs( *(lp.getRHS()) );
657  }
658 
659  resvector_.resize( numrhs_ );
660  testvector_.resize( numrhs_ );
661 
662  curLSNum_ = lp.getLSNumber();
663  curLSIdx_ = lp.getLSIndex();
664  curBlksz_ = (int)curLSIdx_.size();
665  int validLS = 0;
666  for (i=0; i<curBlksz_; ++i) {
667  if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
668  validLS++;
669  }
670  curNumRHS_ = validLS;
671  //
672  // Initialize the testvector.
673  for (i=0; i<numrhs_; i++) { testvector_[i] = one; }
674 
675  // Return an error if the scaling is zero.
676  if (scalevalue_ == zero) {
677  return Failed;
678  }
679  }
680  return Undefined;
681  }
683 
686 
688  std::string description() const
689  {
690  std::ostringstream oss;
691  oss << "Belos::StatusTestGenResSubNorm<>: " << resFormStr();
692  oss << ", tol = " << tolerance_;
693  return oss.str();
694  }
696 
697  protected:
698 
699  private:
700 
702 
703 
704  std::string resFormStr() const
705  {
706  std::ostringstream oss;
707  oss << "(";
708  oss << ((resnormtype_==OneNorm) ? "1-Norm" : (resnormtype_==TwoNorm) ? "2-Norm" : "Inf-Norm");
709  oss << " Exp";
710  oss << " Res Vec [" << subIdx_ << "]) ";
711 
712  // If there is no residual scaling, return current string.
713  if (scaletype_!=None)
714  {
715  // Insert division sign.
716  oss << "/ ";
717 
718  // Determine output string for scaling, if there is any.
719  if (scaletype_==UserProvided)
720  oss << " (User Scale)";
721  else {
722  oss << "(";
723  oss << ((scalenormtype_==OneNorm) ? "1-Norm" : (resnormtype_==TwoNorm) ? "2-Norm" : "Inf-Norm");
724  if (scaletype_==NormOfInitRes)
725  oss << " Res0 [" << subIdx_ << "]";
726  else if (scaletype_==NormOfPrecInitRes)
727  oss << " Prec Res0 [" << subIdx_ << "]";
728  else if (scaletype_==NormOfFullInitRes)
729  oss << " Full Res0 [" << subIdx_ << "]";
730  else if (scaletype_==NormOfFullPrecInitRes)
731  oss << " Full Prec Res0 [" << subIdx_ << "]";
732  else if (scaletype_==NormOfFullScaledInitRes)
733  oss << " scaled Full Res0 [" << subIdx_ << "]";
734  else if (scaletype_==NormOfFullScaledPrecInitRes)
735  oss << " scaled Full Prec Res0 [" << subIdx_ << "]";
736  else
737  oss << " RHS [" << subIdx_ << "]";
738  oss << ")";
739  }
740  }
741 
742  // TODO add a tagging name
743 
744  return oss.str();
745  }
746 
748 
750 
751 
752  // calculate norm of partial multivector
753  void MvSubNorm( const MV& mv, size_t block, std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>& normVec, NormType type = TwoNorm) {
754  Teuchos::RCP<const MV> input = Teuchos::rcpFromRef(mv);
755 
756  typedef typename Thyra::ProductMultiVectorBase<ScalarType> TPMVB;
757  Teuchos::RCP<const TPMVB> thyProdVec = Teuchos::rcp_dynamic_cast<const TPMVB>(input);
758 
759  TEUCHOS_TEST_FOR_EXCEPTION(thyProdVec == Teuchos::null, std::invalid_argument,
760  "Belos::StatusTestGenResSubNorm::MvSubNorm (Thyra specialization): "
761  "mv must be a Thyra::ProductMultiVector, but is of type " << thyProdVec);
762 
763  Teuchos::RCP<const MV> thySubVec = thyProdVec->getMultiVectorBlock(block);
764 
765  MVT::MvNorm(*thySubVec,normVec,type);
766  }
767 
768  // calculate ration of sub-vector length to full vector length (for scalevalue_)
769  void MvScalingRatio( const MV& mv, size_t block, MagnitudeType& lengthRatio) {
770  Teuchos::RCP<const MV> input = Teuchos::rcpFromRef(mv);
771 
772  typedef typename Thyra::ProductMultiVectorBase<ScalarType> TPMVB;
773  Teuchos::RCP<const TPMVB> thyProdVec = Teuchos::rcp_dynamic_cast<const TPMVB>(input);
774 
775  TEUCHOS_TEST_FOR_EXCEPTION(thyProdVec == Teuchos::null, std::invalid_argument,
776  "Belos::StatusTestGenResSubNorm::MvScalingRatio (Thyra specialization): "
777  "mv must be a Thyra::ProductMultiVector, but is of type " << thyProdVec);
778 
779  Teuchos::RCP<const MV> thySubVec = thyProdVec->getMultiVectorBlock(block);
780 
781  lengthRatio = Teuchos::as<MagnitudeType>(thySubVec->range()->dim()) / Teuchos::as<MagnitudeType>(thyProdVec->range()->dim());
782  }
783 
785 
787 
788 
790  MagnitudeType tolerance_;
791 
793  size_t subIdx_;
794 
796  int quorum_;
797 
799  bool showMaxResNormOnly_;
800 
802  NormType resnormtype_;
803 
805  ScaleType scaletype_;
806 
808  NormType scalenormtype_;
809 
811  MagnitudeType scalevalue_;
812 
814  std::vector<MagnitudeType> scalevector_;
815 
817  std::vector<MagnitudeType> resvector_;
818 
820  std::vector<MagnitudeType> testvector_;
821 
823  std::vector<int> ind_;
824 
826  Teuchos::RCP<MV> curSoln_;
827 
829  StatusType status_;
830 
832  int curBlksz_;
833 
835  int curNumRHS_;
836 
838  std::vector<int> curLSIdx_;
839 
841  int curLSNum_;
842 
844  int numrhs_;
845 
847  bool firstcallCheckStatus_;
848 
850  bool firstcallDefineResForm_;
851 
853  bool firstcallDefineScaleForm_;
854 
856 
857 };
858 
859 #endif // HAVE_BELOS_THYRA
860 
861 } // end namespace Belos
862 
863 #endif /* BELOS_STATUS_TEST_RESSUBNORM_H */
Teuchos::RCP< MV > getSolution()
Returns the current solution estimate that was computed for the most recent residual test...
int setQuorum(int)
Sets the number of residuals that must pass the convergence test before Passed is returned...
const std::vector< MagnitudeType > * getResNormValue() const
Returns the residual norm value, , computed in most recent call to CheckStatus.
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:88
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
An implementation of StatusTestResNorm using a family of norms of subvectors of the residual vectors...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Declaration of basic traits for the multivector type.
void print(std::ostream &, int=0) const
Output formatted description of stopping test to output stream.
An abstract class of StatusTest for stopping criteria using residual norms.
int defineScaleForm(ScaleType, NormType, MagnitudeType=Teuchos::ScalarTraits< MagnitudeType >::one())
Define form of the scaling, its norm, its optional weighting std::vector, or, alternatively, define an explicit value.
void reset()
Resets the internal configuration to the initial state.
Class which defines basic traits for the operator type.
StatusType
Whether the StatusTest wants iteration to stop.
Definition: BelosTypes.hpp:157
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
MagnitudeType getTolerance() const
Returns the value of the tolerance, , set in the constructor.
const std::vector< MagnitudeType > * getScaledNormValue() const
Returns the scaled norm value, .
StatusType getStatus() const
Return the result of the most recent CheckStatus call.
MultiVecTraits< ScalarType, MV > MVT
Traits class which defines basic operations on multivectors.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
StatusTestGenResSubNorm(MagnitudeType, size_t, int=-1, bool=false)
Constructor.
void printStatus(std::ostream &, StatusType) const
Print message for each status specific to this stopping test.
std::string description() const
Method to return description of the maximum iteration status test.
StatusType firstCallCheckStatusSetup(Iteration< ScalarType, MV, OP > *iSolver)
Call to setup initial scaling std::vector.
Class which describes the linear problem to be solved by the iterative solver.
int getQuorum() const
Returns the number of residuals that must pass the convergence test before Passed is returned...
int setSubIdx(size_t subIdx)
Set the block index of which we want to check the norm of the sub-residuals.
size_t getSubIdx() const
Returns the index of the block row the norms are calculated for.
const std::vector< MagnitudeType > * getTestValue() const
Returns the test value, , computed in most recent call to CheckStatus.
int setTolerance(MagnitudeType)
Set the value of the tolerance.
Belos::StatusTest abstract class for specifying a residual norm stopping criteria.
Teuchos::ScalarTraits< ScalarType > SCT
static magnitudeType magnitude(ScalarTypea)
NormType
The type of vector norm to compute.
Definition: BelosTypes.hpp:65
bool getLOADetected() const
Returns a boolean indicating a loss of accuracy has been detected in computing the residual...
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
StatusType checkStatus(Iteration< ScalarType, MV, OP > *)
Check convergence status: Passed, Failed, or Undefined.
int setShowMaxResNormOnly(bool)
Set whether the only maximum residual norm is displayed when the print() method is called...
std::vector< int > convIndices()
Returns the std::vector containing the indices of the residuals that passed the test.
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 ...
bool getShowMaxResNormOnly()
Returns whether the only maximum residual norm is displayed when the print() method is called...
int defineResForm(NormType)
Define norm of the residual.

Generated on Fri Dec 20 2024 09:24:49 for Belos by doxygen 1.8.5