Belos  Version of the Day
BelosStatusTestGenResNorm.hpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Belos: Block Linear Solvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 //@HEADER
41 
42 #ifndef BELOS_STATUS_TEST_GEN_RESNORM_H
43 #define BELOS_STATUS_TEST_GEN_RESNORM_H
44 
51 #include "BelosLinearProblem.hpp"
52 #include "BelosMultiVecTraits.hpp"
53 
76 namespace Belos {
77 
78 template <class ScalarType, class MV, class OP>
79 class StatusTestGenResNorm: public StatusTestResNorm<ScalarType,MV,OP> {
80 
81  public:
82 
83  // Convenience typedefs
84  typedef Teuchos::ScalarTraits<ScalarType> SCT;
85  typedef typename SCT::magnitudeType MagnitudeType;
87 
89 
90 
94  enum ResType {Implicit,
97  };
98 
100 
102 
103 
117  StatusTestGenResNorm( MagnitudeType Tolerance, int quorum = -1, bool showMaxResNormOnly = false );
118 
120  virtual ~StatusTestGenResNorm();
122 
124 
125 
127 
134  int defineResForm( ResType TypeOfResidual, NormType TypeOfNorm);
135 
137 
157  int defineScaleForm( ScaleType TypeOfScaling, NormType TypeOfNorm, MagnitudeType ScaleValue = Teuchos::ScalarTraits<MagnitudeType>::one());
158 
160 
163  int setTolerance(MagnitudeType tolerance) {tolerance_ = tolerance; return(0);}
164 
167  int setQuorum(int quorum) {quorum_ = quorum; return(0);}
168 
170  int setShowMaxResNormOnly(bool showMaxResNormOnly) {showMaxResNormOnly_ = showMaxResNormOnly; return(0);}
171 
173 
175 
176 
184 
186  StatusType getStatus() const {return(status_);};
188 
190 
191 
193  void reset();
194 
196 
198 
199 
201  void print(std::ostream& os, int indent = 0) const;
202 
204  void printStatus(std::ostream& os, StatusType type) const;
206 
208 
209 
213  Teuchos::RCP<MV> getSolution() { if (restype_==Implicit) { return Teuchos::null; } else { return curSoln_; } }
214 
217  int getQuorum() const { return quorum_; }
218 
220  bool getShowMaxResNormOnly() { return showMaxResNormOnly_; }
221 
223  std::vector<int> convIndices() { return ind_; }
224 
226  MagnitudeType getTolerance() const {return(tolerance_);};
227 
229  const std::vector<MagnitudeType>* getTestValue() const {return(&testvector_);};
230 
232  const std::vector<MagnitudeType>* getResNormValue() const {return(&resvector_);};
233 
235  const std::vector<MagnitudeType>* getScaledNormValue() const {return(&scalevector_);};
236 
239  bool getLOADetected() const { return false; }
240 
242 
243 
246 
254 
257 
259  std::string description() const
260  {
261  std::ostringstream oss;
262  oss << "Belos::StatusTestGenResNorm<>: " << resFormStr();
263  oss << ", tol = " << tolerance_;
264  return oss.str();
265  }
267 
268  protected:
269 
270  private:
271 
273 
274 
275  std::string resFormStr() const
276  {
277  std::ostringstream oss;
278  oss << "(";
279  oss << ((resnormtype_==OneNorm) ? "1-Norm" : (resnormtype_==TwoNorm) ? "2-Norm" : "Inf-Norm");
280  oss << ((restype_==Explicit) ? " Exp" : " Imp");
281  oss << " Res Vec) ";
282 
283  // If there is no residual scaling, return current string.
284  if (scaletype_!=None)
285  {
286  // Insert division sign.
287  oss << "/ ";
288 
289  // Determine output string for scaling, if there is any.
290  if (scaletype_==UserProvided)
291  oss << " (User Scale)";
292  else {
293  oss << "(";
294  oss << ((scalenormtype_==OneNorm) ? "1-Norm" : (resnormtype_==TwoNorm) ? "2-Norm" : "Inf-Norm");
295  if (scaletype_==NormOfInitRes || scaletype_==NormOfFullInitRes || scaletype_==NormOfFullScaledInitRes)
296  oss << " Res0";
297  else if (scaletype_==NormOfPrecInitRes || scaletype_==NormOfFullPrecInitRes || scaletype_==NormOfFullScaledPrecInitRes)
298  oss << " Prec Res0";
299  else
300  oss << " RHS ";
301  oss << ")";
302  }
303  }
304 
305  return oss.str();
306  }
307 
309 
311 
312 
314  MagnitudeType tolerance_;
315 
317  int quorum_;
318 
320  bool showMaxResNormOnly_;
321 
323  ResType restype_;
324 
326  NormType resnormtype_;
327 
329  ScaleType scaletype_;
330 
332  NormType scalenormtype_;
333 
335  MagnitudeType scalevalue_;
336 
338  std::vector<MagnitudeType> scalevector_;
339 
341  std::vector<MagnitudeType> resvector_;
342 
344  std::vector<MagnitudeType> testvector_;
345 
347  std::vector<int> ind_;
348 
350  Teuchos::RCP<MV> curSoln_;
351 
353  StatusType status_;
354 
356  int curBlksz_;
357 
359  int curNumRHS_;
360 
362  std::vector<int> curLSIdx_;
363 
365  int curLSNum_;
366 
368  int numrhs_;
369 
371  bool firstcallCheckStatus_;
372 
374  bool firstcallDefineResForm_;
375 
377  bool firstcallDefineScaleForm_;
378 
380 
381 };
382 
383 template <class ScalarType, class MV, class OP>
385 StatusTestGenResNorm (MagnitudeType Tolerance, int quorum, bool showMaxResNormOnly)
386  : tolerance_(Tolerance),
387  quorum_(quorum),
388  showMaxResNormOnly_(showMaxResNormOnly),
389  restype_(Implicit),
390  resnormtype_(TwoNorm),
391  scaletype_(NormOfInitRes),
392  scalenormtype_(TwoNorm),
393  scalevalue_(Teuchos::ScalarTraits<MagnitudeType>::one ()),
394  status_(Undefined),
395  curBlksz_(0),
396  curNumRHS_(0),
397  curLSNum_(0),
398  numrhs_(0),
399  firstcallCheckStatus_(true),
400  firstcallDefineResForm_(true),
401  firstcallDefineScaleForm_(true)
402 {
403  // This constructor will compute the residual ||r_i||/||r0_i|| <= tolerance using the 2-norm of
404  // the implicit residual std::vector.
405 }
406 
407 template <class ScalarType, class MV, class OP>
409 {}
410 
411 template <class ScalarType, class MV, class OP>
413 {
414  status_ = Undefined;
415  curBlksz_ = 0;
416  curLSNum_ = 0;
417  curLSIdx_.resize(0);
418  numrhs_ = 0;
419  ind_.resize(0);
420  firstcallCheckStatus_ = true;
421  curSoln_ = Teuchos::null;
422 }
423 
424 template <class ScalarType, class MV, class OP>
426 {
427  TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineResForm_==false,StatusTestError,
428  "StatusTestGenResNorm::defineResForm(): The residual form has already been defined.");
429  firstcallDefineResForm_ = false;
430 
431  restype_ = TypeOfResidual;
432  resnormtype_ = TypeOfNorm;
433 
434  return(0);
435 }
436 
437 template <class ScalarType, class MV, class OP>
439  MagnitudeType ScaleValue )
440 {
441  TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineScaleForm_==false,StatusTestError,
442  "StatusTestGenResNorm::defineScaleForm(): The scaling type has already been defined.");
443  firstcallDefineScaleForm_ = false;
444 
445  scaletype_ = TypeOfScaling;
446  scalenormtype_ = TypeOfNorm;
447  scalevalue_ = ScaleValue;
448 
449  return(0);
450 }
451 
452 template <class ScalarType, class MV, class OP>
454 {
455  MagnitudeType zero = Teuchos::ScalarTraits<MagnitudeType>::zero();
456  const LinearProblem<ScalarType,MV,OP>& lp = iSolver->getProblem();
457  // Compute scaling term (done once for each block that's being solved)
458  if (firstcallCheckStatus_) {
459  StatusType status = firstCallCheckStatusSetup(iSolver);
460  if(status==Failed) {
461  status_ = Failed;
462  return(status_);
463  }
464  }
465  //
466  // This section computes the norm of the residual std::vector
467  //
468  if ( curLSNum_ != lp.getLSNumber() ) {
469  //
470  // We have moved on to the next rhs block
471  //
472  curLSNum_ = lp.getLSNumber();
473  curLSIdx_ = lp.getLSIndex();
474  curBlksz_ = (int)curLSIdx_.size();
475  int validLS = 0;
476  for (int i=0; i<curBlksz_; ++i) {
477  if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
478  validLS++;
479  }
480  curNumRHS_ = validLS;
481  curSoln_ = Teuchos::null;
482  //
483  } else {
484  //
485  // We are in the same rhs block, return if we are converged
486  //
487  if (status_==Passed) { return status_; }
488  }
489  if (restype_==Implicit) {
490  //
491  // get the native residual norms from the solver for this block of right-hand sides.
492  // If the residual is returned in multivector form, use the resnormtype to compute the residual norms.
493  // Otherwise the native residual is assumed to be stored in the resvector_.
494  //
495  std::vector<MagnitudeType> tmp_resvector( curBlksz_ );
496  Teuchos::RCP<const MV> residMV = iSolver->getNativeResiduals( &tmp_resvector );
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) {
502  // Check if this index is valid
503  if (*p != -1)
504  resvector_[*p] = tmp_resvector[i];
505  }
506  } else {
507  typename std::vector<int>::iterator p = curLSIdx_.begin();
508  for (int i=0; p<curLSIdx_.end(); ++p, ++i) {
509  // Check if this index is valid
510  if (*p != -1)
511  resvector_[*p] = tmp_resvector[i];
512  }
513  }
514  }
515  else if (restype_==Explicit) {
516  //
517  // Request the true residual for this block of right-hand sides.
518  //
519  Teuchos::RCP<MV> cur_update = iSolver->getCurrentUpdate();
520  curSoln_ = lp.updateSolution( cur_update );
521  Teuchos::RCP<MV> cur_res = MVT::Clone( *curSoln_, MVT::GetNumberVecs( *curSoln_ ) );
522  lp.computeCurrResVec( &*cur_res, &*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) {
527  // Check if this index is valid
528  if (*p != -1)
529  resvector_[*p] = tmp_resvector[i];
530  }
531  }
532  //
533  // Compute the new linear system residuals for testing.
534  // (if any of them don't meet the tolerance or are NaN, then we exit with that status)
535  //
536  if ( scalevector_.size() > 0 ) {
537  typename std::vector<int>::iterator p = curLSIdx_.begin();
538  for (; p<curLSIdx_.end(); ++p) {
539  // Check if this index is valid
540  if (*p != -1) {
541  // Scale the std::vector accordingly
542  if ( scalevector_[ *p ] != zero ) {
543  // Don't intentionally divide by zero.
544  testvector_[ *p ] = resvector_[ *p ] / scalevector_[ *p ] / scalevalue_;
545  } else {
546  testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
547  }
548  }
549  }
550  }
551  else {
552  typename std::vector<int>::iterator p = curLSIdx_.begin();
553  for (; p<curLSIdx_.end(); ++p) {
554  // Check if this index is valid
555  if (*p != -1)
556  testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
557  }
558  }
559 
560  // Check status of new linear system residuals and see if we have the quorum.
561  int have = 0;
562  ind_.resize( curLSIdx_.size() );
563  typename std::vector<int>::iterator p = curLSIdx_.begin();
564  for (; p<curLSIdx_.end(); ++p) {
565  // Check if this index is valid
566  if (*p != -1) {
567  // Check if any of the residuals are larger than the tolerance.
568  if (testvector_[ *p ] > tolerance_) {
569  // do nothing.
570  } else if (testvector_[ *p ] <= tolerance_) {
571  ind_[have] = *p;
572  have++;
573  } else {
574  // Throw an std::exception if a NaN is found.
575  status_ = Failed;
576  TEUCHOS_TEST_FOR_EXCEPTION(true,StatusTestError,"StatusTestGenResNorm::checkStatus(): NaN has been detected.");
577  }
578  }
579  }
580  ind_.resize(have);
581  int need = (quorum_ == -1) ? curNumRHS_: quorum_;
582  status_ = (have >= need) ? Passed : Failed;
583 
584  // Return the current status
585  return status_;
586 }
587 
588 template <class ScalarType, class MV, class OP>
589 void StatusTestGenResNorm<ScalarType,MV,OP>::print(std::ostream& os, int indent) const
590 {
591  for (int j = 0; j < indent; j ++)
592  os << ' ';
593  printStatus(os, status_);
594  os << resFormStr();
595  if (status_==Undefined)
596  os << ", tol = " << tolerance_ << std::endl;
597  else {
598  os << std::endl;
599  if(showMaxResNormOnly_ && curBlksz_ > 1) {
600  const MagnitudeType maxRelRes = *std::max_element(
601  testvector_.begin()+curLSIdx_[0],testvector_.begin()+curLSIdx_[curBlksz_-1]
602  );
603  for (int j = 0; j < indent + 13; j ++)
604  os << ' ';
605  os << "max{residual["<<curLSIdx_[0]<<"..."<<curLSIdx_[curBlksz_-1]<<"]} = " << maxRelRes
606  << ( maxRelRes <= tolerance_ ? " <= " : " > " ) << tolerance_ << std::endl;
607  }
608  else {
609  for ( int i=0; i<numrhs_; i++ ) {
610  for (int j = 0; j < indent + 13; j ++)
611  os << ' ';
612  os << "residual [ " << i << " ] = " << testvector_[ i ];
613  os << ((testvector_[i]<tolerance_) ? " < " : (testvector_[i]==tolerance_) ? " == " : (testvector_[i]>tolerance_) ? " > " : " " ) << tolerance_ << std::endl;
614  }
615  }
616  }
617  os << std::endl;
618 }
619 
620 template <class ScalarType, class MV, class OP>
622 {
623  os << std::left << std::setw(13) << std::setfill('.');
624  switch (type) {
625  case Passed:
626  os << "Converged";
627  break;
628  case Failed:
629  os << "Unconverged";
630  break;
631  case Undefined:
632  default:
633  os << "**";
634  break;
635  }
636  os << std::left << std::setfill(' ');
637  return;
638 }
639 
640 template <class ScalarType, class MV, class OP>
642 {
643  int i;
644  MagnitudeType zero = Teuchos::ScalarTraits<MagnitudeType>::zero();
645  MagnitudeType one = Teuchos::ScalarTraits<MagnitudeType>::one();
646  const LinearProblem<ScalarType,MV,OP>& lp = iSolver->getProblem();
647  // Compute scaling term (done once for each block that's being solved)
648  if (firstcallCheckStatus_) {
649  //
650  // Get some current solver information.
651  //
652  firstcallCheckStatus_ = false;
653 
654  if (scaletype_== NormOfRHS) {
655  Teuchos::RCP<const MV> rhs = lp.getRHS();
656  numrhs_ = MVT::GetNumberVecs( *rhs );
657  scalevector_.resize( numrhs_ );
658  MVT::MvNorm( *rhs, scalevector_, scalenormtype_ );
659  }
660  else if (scaletype_==NormOfInitRes || scaletype_==NormOfFullInitRes || scaletype_==NormOfFullScaledInitRes) {
661  Teuchos::RCP<const MV> init_res = lp.getInitResVec();
662  numrhs_ = MVT::GetNumberVecs( *init_res );
663  scalevector_.resize( numrhs_ );
664  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
665  }
666  else if (scaletype_==NormOfPrecInitRes || scaletype_==NormOfFullPrecInitRes || scaletype_==NormOfFullScaledPrecInitRes) {
667  Teuchos::RCP<const MV> init_res = lp.getInitPrecResVec();
668  numrhs_ = MVT::GetNumberVecs( *init_res );
669  scalevector_.resize( numrhs_ );
670  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
671  }
672  else {
673  numrhs_ = MVT::GetNumberVecs( *(lp.getRHS()) );
674  }
675 
676  resvector_.resize( numrhs_ );
677  testvector_.resize( numrhs_ );
678 
679  curLSNum_ = lp.getLSNumber();
680  curLSIdx_ = lp.getLSIndex();
681  curBlksz_ = (int)curLSIdx_.size();
682  int validLS = 0;
683  for (i=0; i<curBlksz_; ++i) {
684  if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
685  validLS++;
686  }
687  curNumRHS_ = validLS;
688  //
689  // Initialize the testvector.
690  for (i=0; i<numrhs_; i++) { testvector_[i] = one; }
691 
692  // Return an error if the scaling is zero.
693  if (scalevalue_ == zero) {
694  return Failed;
695  }
696  }
697  return Undefined;
698 }
699 
700 } // end namespace Belos
701 
702 #endif /* BELOS_STATUS_TEST_RESNORM_H */
Belos::ScaleType
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:119
Belos::StatusTestGenResNorm::~StatusTestGenResNorm
virtual ~StatusTestGenResNorm()
Destructor.
Definition: BelosStatusTestGenResNorm.hpp:408
Belos::LinearProblem::getLSNumber
int getLSNumber() const
The number of linear systems that have been set.
Definition: BelosLinearProblem.hpp:404
Belos::StatusTestGenResNorm::getSolution
Teuchos::RCP< MV > getSolution()
Returns the current solution estimate that was computed for the most recent residual test.
Definition: BelosStatusTestGenResNorm.hpp:213
Belos::StatusTestGenResNorm::defineScaleForm
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,...
Definition: BelosStatusTestGenResNorm.hpp:438
Belos::StatusTestGenResNorm::firstCallCheckStatusSetup
StatusType firstCallCheckStatusSetup(Iteration< ScalarType, MV, OP > *iSolver)
Call to setup initial scaling std::vector.
Definition: BelosStatusTestGenResNorm.hpp:641
BelosLinearProblem.hpp
Class which describes the linear problem to be solved by the iterative solver.
Belos::NormOfFullPrecInitRes
Definition: BelosTypes.hpp:125
Belos::LinearProblem::computeCurrResVec
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.
Definition: BelosLinearProblem.hpp:1196
Belos::NormOfFullScaledPrecInitRes
Definition: BelosTypes.hpp:127
Belos::Iteration
Definition: BelosIteration.hpp:73
Belos::StatusTestError
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
Definition: BelosStatusTest.hpp:73
Belos::None
Definition: BelosTypes.hpp:122
Belos::StatusTestGenResNorm::description
std::string description() const
Method to return description of the maximum iteration status test
Definition: BelosStatusTestGenResNorm.hpp:259
BelosStatusTestResNorm.hpp
Belos::StatusTest abstract class for specifying a residual norm stopping criteria.
Belos::StatusTestGenResNorm::print
void print(std::ostream &os, int indent=0) const
Output formatted description of stopping test to output stream.
Definition: BelosStatusTestGenResNorm.hpp:589
Belos::StatusTestGenResNorm::printStatus
void printStatus(std::ostream &os, StatusType type) const
Print message for each status specific to this stopping test.
Definition: BelosStatusTestGenResNorm.hpp:621
Belos::Passed
Definition: BelosTypes.hpp:188
Belos::TwoNorm
Definition: BelosTypes.hpp:98
Belos::LinearProblem
A linear system to solve, and its associated information.
Definition: BelosIteration.hpp:61
Belos::LinearProblem::getInitPrecResVec
Teuchos::RCP< const MV > getInitPrecResVec() const
A pointer to the preconditioned initial residual vector.
Definition: BelosLinearProblem.hpp:958
Belos::OneNorm
Definition: BelosTypes.hpp:97
Belos::StatusTestGenResNorm::SCT
Teuchos::ScalarTraits< ScalarType > SCT
Definition: BelosStatusTestGenResNorm.hpp:84
Belos::NormOfFullScaledInitRes
Definition: BelosTypes.hpp:126
Belos
Definition: Belos_Details_EBelosSolverType.cpp:45
Belos::StatusType
StatusType
Whether the StatusTest wants iteration to stop.
Definition: BelosTypes.hpp:188
Belos::StatusTestGenResNorm::setQuorum
int setQuorum(int quorum)
Sets the number of residuals that must pass the convergence test before Passed is returned.
Definition: BelosStatusTestGenResNorm.hpp:167
Belos::Iteration::getCurrentUpdate
virtual Teuchos::RCP< MV > getCurrentUpdate() const =0
Get the current update to the linear system.
Belos::Undefined
Definition: BelosTypes.hpp:190
Belos::StatusTestGenResNorm::MagnitudeType
SCT::magnitudeType MagnitudeType
Definition: BelosStatusTestGenResNorm.hpp:85
Belos::StatusTestGenResNorm::getShowMaxResNormOnly
bool getShowMaxResNormOnly()
Returns whether the only maximum residual norm is displayed when the print() method is called.
Definition: BelosStatusTestGenResNorm.hpp:220
Belos::LinearProblem::getLSIndex
const std::vector< int > getLSIndex() const
(Zero-based) indices of the linear system(s) currently being solved.
Definition: BelosLinearProblem.hpp:396
Belos::UserProvided
Definition: BelosTypes.hpp:123
Belos::StatusTestGenResNorm::MVT
MultiVecTraits< ScalarType, MV > MVT
Definition: BelosStatusTestGenResNorm.hpp:86
Belos::StatusTestGenResNorm::getResNormValue
const std::vector< MagnitudeType > * getResNormValue() const
Returns the residual norm value, , computed in most recent call to CheckStatus.
Definition: BelosStatusTestGenResNorm.hpp:232
Belos::NormOfRHS
Definition: BelosTypes.hpp:119
Belos::StatusTestGenResNorm::setTolerance
int setTolerance(MagnitudeType tolerance)
Set the value of the tolerance.
Definition: BelosStatusTestGenResNorm.hpp:163
Belos::LinearProblem::getRHS
Teuchos::RCP< const MV > getRHS() const
A pointer to the right-hand side B.
Definition: BelosLinearProblem.hpp:326
Belos::StatusTestGenResNorm
An implementation of StatusTestResNorm using a family of residual norms.
Definition: BelosStatusTestGenResNorm.hpp:79
Belos::StatusTestGenResNorm::getLOADetected
bool getLOADetected() const
Returns a boolean indicating a loss of accuracy has been detected in computing the residual.
Definition: BelosStatusTestGenResNorm.hpp:239
Belos::Iteration::getNativeResiduals
virtual Teuchos::RCP< const MV > getNativeResiduals(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > *norms) const =0
Get the residuals native to the solver.
Belos::StatusTestGenResNorm::defineResForm
int defineResForm(ResType TypeOfResidual, NormType TypeOfNorm)
Define form of the residual, its norm and optional weighting std::vector.
Definition: BelosStatusTestGenResNorm.hpp:425
Belos::StatusTestResNorm
An abstract class of StatusTest for stopping criteria using residual norms.
Definition: BelosStatusTestResNorm.hpp:62
Belos::StatusTestGenResNorm::ResType
ResType
Select how the residual std::vector is produced.
Definition: BelosStatusTestGenResNorm.hpp:94
Belos::LinearProblem::updateSolution
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.
Definition: BelosLinearProblem.hpp:745
Belos::StatusTestGenResNorm::getScaledNormValue
const std::vector< MagnitudeType > * getScaledNormValue() const
Returns the scaled norm value, .
Definition: BelosStatusTestGenResNorm.hpp:235
Belos::NormType
NormType
The type of vector norm to compute.
Definition: BelosTypes.hpp:97
Belos::StatusTestGenResNorm::reset
void reset()
Resets the internal configuration to the initial state.
Definition: BelosStatusTestGenResNorm.hpp:412
Belos::StatusTestGenResNorm::getTolerance
MagnitudeType getTolerance() const
Returns the value of the tolerance, , set in the constructor.
Definition: BelosStatusTestGenResNorm.hpp:226
Belos::LinearProblem::getInitResVec
Teuchos::RCP< const MV > getInitResVec() const
A pointer to the initial unpreconditioned residual vector.
Definition: BelosLinearProblem.hpp:949
BelosMultiVecTraits.hpp
Declaration of basic traits for the multivector type.
Belos::StatusTestGenResNorm::StatusTestGenResNorm
StatusTestGenResNorm(MagnitudeType Tolerance, int quorum=-1, bool showMaxResNormOnly=false)
Constructor.
Definition: BelosStatusTestGenResNorm.hpp:385
Belos::StatusTestGenResNorm::convIndices
std::vector< int > convIndices()
Returns the std::vector containing the indices of the residuals that passed the test.
Definition: BelosStatusTestGenResNorm.hpp:223
Belos::Iteration::getProblem
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const =0
Get a constant reference to the linear problem.
Belos::StatusTestGenResNorm::checkStatus
StatusType checkStatus(Iteration< ScalarType, MV, OP > *iSolver)
Check convergence status: Passed, Failed, or Undefined.
Definition: BelosStatusTestGenResNorm.hpp:453
Belos::NormOfInitRes
Definition: BelosTypes.hpp:120
Belos::NormOfFullInitRes
Definition: BelosTypes.hpp:124
Belos::StatusTestGenResNorm::Explicit
Definition: BelosStatusTestGenResNorm.hpp:95
Belos::NormOfPrecInitRes
Definition: BelosTypes.hpp:121
Belos::Failed
Definition: BelosTypes.hpp:189
Belos::StatusTestGenResNorm::Implicit
Definition: BelosStatusTestGenResNorm.hpp:94
Belos::MultiVecTraits
Traits class which defines basic operations on multivectors.
Definition: BelosMultiVecTraits.hpp:129
Belos::StatusTestGenResNorm::getTestValue
const std::vector< MagnitudeType > * getTestValue() const
Returns the test value, , computed in most recent call to CheckStatus.
Definition: BelosStatusTestGenResNorm.hpp:229
Belos::StatusTestGenResNorm::getQuorum
int getQuorum() const
Returns the number of residuals that must pass the convergence test before Passed is returned.
Definition: BelosStatusTestGenResNorm.hpp:217
Belos::StatusTestGenResNorm::getStatus
StatusType getStatus() const
Return the result of the most recent CheckStatus call.
Definition: BelosStatusTestGenResNorm.hpp:186
Belos::StatusTestGenResNorm::setShowMaxResNormOnly
int setShowMaxResNormOnly(bool showMaxResNormOnly)
Set whether the only maximum residual norm is displayed when the print() method is called.
Definition: BelosStatusTestGenResNorm.hpp:170

Generated on Thu Feb 27 2020 16:06:46 for Belos by doxygen 1.8.16