Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosStatusTestImpResNorm.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_IMPRESNORM_H
43 #define BELOS_STATUS_TEST_IMPRESNORM_H
44 
51 #include "BelosLinearProblem.hpp"
52 #include "BelosMultiVecTraits.hpp"
53 #include "Teuchos_as.hpp"
54 
101 namespace Belos {
102 
103 template <class ScalarType, class MV, class OP>
104 class StatusTestImpResNorm: public StatusTestResNorm<ScalarType,MV,OP> {
105 public:
108 
109 private:
111 
116 
117 public:
119 
120 
133  int quorum = -1,
134  bool showMaxResNormOnly = false);
135 
137  virtual ~StatusTestImpResNorm();
138 
140 
142 
144 
150  int defineResForm( NormType TypeOfNorm);
151 
153 
173  int defineScaleForm( ScaleType TypeOfScaling, NormType TypeOfNorm, MagnitudeType ScaleValue = Teuchos::ScalarTraits<MagnitudeType>::one());
174 
176 
179  int setTolerance (MagnitudeType tolerance) {
180  tolerance_ = tolerance;
181  return 0;
182  }
183 
186  int setQuorum (int quorum) {
187  quorum_ = quorum;
188  return 0;
189  }
190 
192  int setShowMaxResNormOnly (bool showMaxResNormOnly) {
193  showMaxResNormOnly_ = showMaxResNormOnly;
194  return 0;
195  }
196 
198 
200 
201 
209 
211  StatusType getStatus() const {return(status_);};
213 
215 
216 
218  void reset();
219 
221 
223 
224 
226  void print(std::ostream& os, int indent = 0) const;
227 
229  void printStatus(std::ostream& os, StatusType type) const;
231 
233 
234 
236  Teuchos::RCP<MV> getSolution() { return curSoln_; }
237 
240  int getQuorum() const { return quorum_; }
241 
243  bool getShowMaxResNormOnly() { return showMaxResNormOnly_; }
244 
246  std::vector<int> convIndices() { return ind_; }
247 
256  return tolerance_;
257  }
258 
267  return currTolerance_;
268  }
269 
271  const std::vector<MagnitudeType>* getTestValue() const {return(&testvector_);};
272 
274  const std::vector<MagnitudeType>* getResNormValue() const {return(&resvector_);};
275 
277  const std::vector<MagnitudeType>* getScaledNormValue() const {return(&scalevector_);};
278 
280  bool getLOADetected() const { return lossDetected_; }
282 
283 
286 
294 
297 
299  std::string description() const
300  {
301  std::ostringstream oss;
302  oss << "Belos::StatusTestImpResNorm<>: " << resFormStr();
303  oss << ", tol = " << tolerance_;
304  return oss.str();
305  }
307 
308  protected:
309 
310  private:
311 
313 
314 
315  std::string resFormStr() const
316  {
317  std::ostringstream oss;
318  oss << "(";
319  oss << ((resnormtype_==OneNorm) ? "1-Norm" : (resnormtype_==TwoNorm) ? "2-Norm" : "Inf-Norm");
320  oss << " Res Vec) ";
321 
322  // If there is no residual scaling, return current string.
323  if (scaletype_!=None)
324  {
325  // Insert division sign.
326  oss << "/ ";
327 
328  // Determine output string for scaling, if there is any.
329  if (scaletype_==UserProvided)
330  oss << " (User Scale)";
331  else {
332  oss << "(";
333  oss << ((scalenormtype_==OneNorm) ? "1-Norm" : (resnormtype_==TwoNorm) ? "2-Norm" : "Inf-Norm");
334  if (scaletype_==NormOfInitRes)
335  oss << " Res0";
336  else if (scaletype_==NormOfPrecInitRes)
337  oss << " Prec Res0";
338  else
339  oss << " RHS ";
340  oss << ")";
341  }
342  }
343 
344  return oss.str();
345  }
346 
348 
349 
351  MagnitudeType tolerance_, currTolerance_;
352 
354  int quorum_;
355 
357  bool showMaxResNormOnly_;
358 
360  NormType resnormtype_;
361 
363  ScaleType scaletype_;
364 
366  NormType scalenormtype_;
367 
369  MagnitudeType scalevalue_;
370 
372  std::vector<MagnitudeType> scalevector_;
373 
375  std::vector<MagnitudeType> resvector_;
376 
378  std::vector<MagnitudeType> testvector_;
379 
381  Teuchos::RCP<MV> curSoln_;
382 
384  std::vector<int> ind_;
385 
387  StatusType status_;
388 
390  int curBlksz_;
391 
393  int curNumRHS_;
394 
396  std::vector<int> curLSIdx_;
397 
399  int curLSNum_;
400 
402  int numrhs_;
403 
405  bool firstcallCheckStatus_;
406 
408  bool firstcallDefineResForm_;
409 
411  bool firstcallDefineScaleForm_;
412 
414  bool lossDetected_;
416 
417 };
418 
419 template <class ScalarType, class MV, class OP>
421 StatusTestImpResNorm (MagnitudeType Tolerance, int quorum, bool showMaxResNormOnly)
422  : tolerance_(Tolerance),
423  currTolerance_(Tolerance),
424  quorum_(quorum),
425  showMaxResNormOnly_(showMaxResNormOnly),
426  resnormtype_(TwoNorm),
427  scaletype_(NormOfInitRes),
428  scalenormtype_(TwoNorm),
429  scalevalue_(Teuchos::ScalarTraits<MagnitudeType>::one ()),
430  status_(Undefined),
431  curBlksz_(0),
432  curNumRHS_(0),
433  curLSNum_(0),
434  numrhs_(0),
435  firstcallCheckStatus_(true),
436  firstcallDefineResForm_(true),
437  firstcallDefineScaleForm_(true),
438  lossDetected_(false)
439 {
440  // This constructor will compute the residual ||r_i||/||r0_i|| <= tolerance using the 2-norm of
441  // the implicit residual vector.
442 }
443 
444 template <class ScalarType, class MV, class OP>
446 {}
447 
448 template <class ScalarType, class MV, class OP>
450 {
451  status_ = Undefined;
452  curBlksz_ = 0;
453  curLSNum_ = 0;
454  curLSIdx_.resize(0);
455  numrhs_ = 0;
456  ind_.resize(0);
457  currTolerance_ = tolerance_;
458  firstcallCheckStatus_ = true;
459  lossDetected_ = false;
460  curSoln_ = Teuchos::null;
461 }
462 
463 template <class ScalarType, class MV, class OP>
465 {
466  TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineResForm_==false,StatusTestError,
467  "StatusTestResNorm::defineResForm(): The residual form has already been defined.");
468  firstcallDefineResForm_ = false;
469 
470  resnormtype_ = TypeOfNorm;
471 
472  return(0);
473 }
474 
475 template <class ScalarType, class MV, class OP>
477  MagnitudeType ScaleValue )
478 {
479  TEUCHOS_TEST_FOR_EXCEPTION(firstcallDefineScaleForm_==false,StatusTestError,
480  "StatusTestResNorm::defineScaleForm(): The scaling type has already been defined.");
481  firstcallDefineScaleForm_ = false;
482 
483  scaletype_ = TypeOfScaling;
484  scalenormtype_ = TypeOfNorm;
485  scalevalue_ = ScaleValue;
486 
487  return(0);
488 }
489 
490 template <class ScalarType, class MV, class OP>
493 {
494  using Teuchos::as;
495  using Teuchos::RCP;
496 
497  const MagnitudeType zero = STM::zero ();
498  const LinearProblem<ScalarType,MV,OP>& lp = iSolver->getProblem ();
499 
500  // Compute scaling term (done once for each block that's being solved)
501  if (firstcallCheckStatus_) {
502  StatusType status = firstCallCheckStatusSetup (iSolver);
503  if (status == Failed) {
504  status_ = Failed;
505  return status_;
506  }
507  }
508 
509  // mfh 23 Apr 2012: I don't know exactly what this code does. It
510  // has something to do with picking the block of right-hand sides
511  // which we're currently checking.
512  if (curLSNum_ != lp.getLSNumber ()) {
513  //
514  // We have moved on to the next rhs block
515  //
516  curLSNum_ = lp.getLSNumber();
517  curLSIdx_ = lp.getLSIndex();
518  curBlksz_ = (int)curLSIdx_.size();
519  int validLS = 0;
520  for (int i=0; i<curBlksz_; ++i) {
521  if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
522  validLS++;
523  }
524  curNumRHS_ = validLS;
525  curSoln_ = Teuchos::null;
526  } else {
527  //
528  // We are in the same rhs block, return if we are converged
529  //
530  if (status_ == Passed) {
531  return status_;
532  }
533  }
534 
535  //
536  // Get the "native" residual norms from the solver for this block of
537  // right-hand sides. If the solver's getNativeResiduals() method
538  // actually returns a multivector, compute the norms of the columns
539  // of the multivector explicitly. Otherwise, we assume that
540  // resvector_ contains the norms.
541  //
542  // Note that "compute the norms explicitly" doesn't necessarily mean
543  // the "explicit" residual norms (in the sense discussed in this
544  // class' documentation). These are just some vectors returned by
545  // the solver. Some Krylov methods, like CG, compute a residual
546  // vector "recursively." This is an "implicit residual" in the
547  // sense of this class' documentation. It equals the explicit
548  // residual in exact arithmetic, but due to rounding error, it is
549  // usually different than the explicit residual.
550  //
551  // FIXME (mfh 23 Apr 2012) This method does _not_ respect the
552  // OrthoManager used by the solver.
553  //
554  std::vector<MagnitudeType> tmp_resvector( curBlksz_ );
555  RCP<const MV> residMV = iSolver->getNativeResiduals (&tmp_resvector);
556  if (! residMV.is_null ()) {
557  // We got a multivector back. Compute the norms explicitly.
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) {
562  // Check if this index is valid
563  if (*p != -1) {
564  resvector_[*p] = tmp_resvector[i];
565  }
566  }
567  } else {
568  typename std::vector<int>::iterator p = curLSIdx_.begin();
569  for (int i=0; p<curLSIdx_.end(); ++p, ++i) {
570  // Check if this index is valid
571  if (*p != -1) {
572  resvector_[*p] = tmp_resvector[i];
573  }
574  }
575  }
576  //
577  // Scale the unscaled residual norms we computed or obtained above.
578  //
579  if (scalevector_.size () > 0) {
580  // There are per-vector scaling factors to apply.
581  typename std::vector<int>::iterator p = curLSIdx_.begin();
582  for (; p<curLSIdx_.end(); ++p) {
583  // Check if this index is valid
584  if (*p != -1) {
585  // Scale the vector accordingly
586  if ( scalevector_[ *p ] != zero ) {
587  // Don't intentionally divide by zero.
588  testvector_[ *p ] = resvector_[ *p ] / scalevector_[ *p ] / scalevalue_;
589  } else {
590  testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
591  }
592  }
593  }
594  }
595  else { // There are no per-vector scaling factors.
596  typename std::vector<int>::iterator p = curLSIdx_.begin();
597  for (; p<curLSIdx_.end(); ++p) {
598  // Check if this index is valid
599  if (*p != -1) {
600  testvector_[ *p ] = resvector_[ *p ] / scalevalue_;
601  }
602  }
603  }
604 
605  // Count how many scaled residual norms (in testvector_) pass, using
606  // the current tolerance (currTolerance_) rather than the original
607  // tolerance (tolerance_). If at least quorum_ of them pass, we
608  // have a quorum for the whole test to pass.
609  //
610  // We also check here whether any of the scaled residual norms is
611  // NaN, and throw an exception in that case.
612  int have = 0;
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) {
617  // Check if this index is valid
618  if (*p != -1) {
619  if (testvector_[ *p ] > currTolerance_) {
620  // The current residual norm doesn't pass. Do nothing.
621  } else if (testvector_[ *p ] <= currTolerance_) {
622  ind_[have] = *p;
623  lclInd[have] = i;
624  have++; // Yay, the current residual norm passes!
625  } else {
626  // Throw an std::exception if the current residual norm is
627  // NaN. We know that it's NaN because it is not less than,
628  // equal to, or greater than the current tolerance. This is
629  // only possible if either the residual norm or the current
630  // tolerance is NaN; we assume the former. We also mark the
631  // test as failed, in case you want to catch the exception.
632  status_ = Failed;
634  "StatusTestImpResNorm::checkStatus(): One or more of the current "
635  "implicit residual norms is NaN.");
636  }
637  }
638  }
639  // "have" is the number of residual norms that passed.
640  ind_.resize(have);
641  lclInd.resize(have);
642 
643  // Now check the exact residuals
644  if (have) { // At least one residual norm has converged.
645  //
646  // Compute the explicit residual norm(s) from the current solution update.
647  //
648  RCP<MV> cur_update = iSolver->getCurrentUpdate ();
649  curSoln_ = lp.updateSolution (cur_update);
650  RCP<MV> cur_res = MVT::Clone (*curSoln_, MVT::GetNumberVecs (*curSoln_));
651  lp.computeCurrResVec (&*cur_res, &*curSoln_);
652  tmp_resvector.resize (MVT::GetNumberVecs (*cur_res));
653  std::vector<MagnitudeType> tmp_testvector (have);
654  MVT::MvNorm (*cur_res, tmp_resvector, resnormtype_);
655 
656  // Scale the explicit residual norm(s), just like the implicit norm(s).
657  if ( scalevector_.size() > 0 ) {
658  for (int i=0; i<have; ++i) {
659  // Scale the vector accordingly
660  if ( scalevector_[ ind_[i] ] != zero ) {
661  // Don't intentionally divide by zero.
662  tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevector_[ ind_[i] ] / scalevalue_;
663  } else {
664  tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevalue_;
665  }
666  }
667  }
668  else {
669  for (int i=0; i<have; ++i) {
670  tmp_testvector[ i ] = tmp_resvector[ lclInd[i] ] / scalevalue_;
671  }
672  }
673 
674  //
675  // Check whether the explicit residual norms also pass the
676  // convergence test. If not, check whether we want to try
677  // iterating a little more to force both implicit and explicit
678  // residual norms to pass.
679  //
680  int have2 = 0;
681  for (int i = 0; i < have; ++i) {
682  // testvector_ contains the implicit (i.e., recursive, computed
683  // by the algorithm) (possibly scaled) residuals. All of these
684  // pass the convergence test.
685  //
686  // tmp_testvector contains the explicit (i.e., ||B-AX||)
687  // (possibly scaled) residuals. We're checking whether these
688  // pass as well. The explicit residual norms only have to meet
689  // the _original_ tolerance (tolerance_), not the current
690  // tolerance (currTolerance_).
691  if (tmp_testvector[i] <= tolerance_) {
692  ind_[have2] = ind_[i];
693  have2++; // This right-hand side has converged.
694  }
695  else {
696  // Absolute difference between the current explicit and
697  // implicit residual norm.
698  const MagnitudeType diff = STM::magnitude (testvector_[ind_[i]] - tmp_testvector[i]);
699  if (diff > currTolerance_) {
700  // If the above difference is bigger than the current
701  // convergence tolerance, report a loss of accuracy, but
702  // mark this right-hand side converged. (The latter tells
703  // users not to iterate further on this right-hand side,
704  // since it probably won't do much good.) Note that the
705  // current tolerance may have been changed by the branch
706  // below in a previous call to this method.
707  lossDetected_ = true;
708  ind_[have2] = ind_[i];
709  have2++; // Count this right-hand side as converged.
710  }
711  else {
712  // Otherwise, the explicit and implicit residuals are pretty
713  // close together, and the implicit residual has converged,
714  // but the explicit residual hasn't converged. Reduce the
715  // convergence tolerance by some formula related to the
716  // difference, and keep iterating.
717  //
718  // mfh 23 Apr 2012: I have no idea why the currTolerance_
719  // update formula is done the way it's done below. It
720  // doesn't make sense to me. It definitely makes the
721  // current tolerance smaller, though, which is what we want.
722 
723  // We define these constants in this way, rather than simply
724  // writing 1.5 resp. 0.1, to ensure no rounding error from
725  // translating from float to MagnitudeType. Remember that
726  // 0.1 doesn't have an exact representation in binary finite
727  // floating-point arithmetic.
728  const MagnitudeType onePointFive = as<MagnitudeType>(3) / as<MagnitudeType> (2);
729  const MagnitudeType oneTenth = STM::one () / as<MagnitudeType> (10);
730 
731  currTolerance_ = currTolerance_ - onePointFive * diff;
732  while (currTolerance_ < STM::zero ()) {
733  currTolerance_ += oneTenth * diff;
734  }
735  }
736  }
737  }
738  have = have2;
739  ind_.resize(have);
740  }
741 
742  // Check whether we've met the quorum of vectors necessary for the
743  // whole test to pass.
744  int need = (quorum_ == -1) ? curNumRHS_: quorum_;
745  status_ = (have >= need) ? Passed : Failed;
746 
747  // Return the current status
748  return status_;
749 }
750 
751 template <class ScalarType, class MV, class OP>
752 void StatusTestImpResNorm<ScalarType,MV,OP>::print(std::ostream& os, int indent) const
753 {
754  for (int j = 0; j < indent; j ++)
755  os << ' ';
756  printStatus(os, status_);
757  os << resFormStr();
758  if (status_==Undefined)
759  os << ", tol = " << tolerance_ << std::endl;
760  else {
761  os << std::endl;
762  if(showMaxResNormOnly_ && curBlksz_ > 1) {
763  const MagnitudeType maxRelRes = *std::max_element(
764  testvector_.begin()+curLSIdx_[0],testvector_.begin()+curLSIdx_[curBlksz_-1]
765  );
766  for (int j = 0; j < indent + 13; j ++)
767  os << ' ';
768  os << "max{residual["<<curLSIdx_[0]<<"..."<<curLSIdx_[curBlksz_-1]<<"]} = " << maxRelRes
769  << ( maxRelRes <= tolerance_ ? " <= " : " > " ) << tolerance_ << std::endl;
770  }
771  else {
772  for ( int i=0; i<numrhs_; i++ ) {
773  for (int j = 0; j < indent + 13; j ++)
774  os << ' ';
775  os << "residual [ " << i << " ] = " << testvector_[ i ];
776  os << ((testvector_[i]<tolerance_) ? " < " : (testvector_[i]==tolerance_) ? " == " : (testvector_[i]>tolerance_) ? " > " : " " ) << tolerance_ << std::endl;
777  }
778  }
779  }
780  os << std::endl;
781 }
782 
783 template <class ScalarType, class MV, class OP>
785 {
786  os << std::left << std::setw(13) << std::setfill('.');
787  switch (type) {
788  case Passed:
789  os << "Converged";
790  break;
791  case Failed:
792  if (lossDetected_)
793  os << "Unconverged (LoA)";
794  else
795  os << "Unconverged";
796  break;
797  case Undefined:
798  default:
799  os << "**";
800  break;
801  }
802  os << std::left << std::setfill(' ');
803  return;
804 }
805 
806 template <class ScalarType, class MV, class OP>
809 {
810  int i;
811  const MagnitudeType zero = STM::zero ();
812  const MagnitudeType one = STM::one ();
813  const LinearProblem<ScalarType,MV,OP>& lp = iSolver->getProblem();
814  // Compute scaling term (done once for each block that's being solved)
815  if (firstcallCheckStatus_) {
816  //
817  // Get some current solver information.
818  //
819  firstcallCheckStatus_ = false;
820 
821  if (scaletype_== NormOfRHS) {
822  Teuchos::RCP<const MV> rhs = lp.getRHS();
823  numrhs_ = MVT::GetNumberVecs( *rhs );
824  scalevector_.resize( numrhs_ );
825  MVT::MvNorm( *rhs, scalevector_, scalenormtype_ );
826  }
827  else if (scaletype_==NormOfInitRes) {
828  Teuchos::RCP<const MV> init_res = lp.getInitResVec();
829  numrhs_ = MVT::GetNumberVecs( *init_res );
830  scalevector_.resize( numrhs_ );
831  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
832  }
833  else if (scaletype_==NormOfPrecInitRes) {
835  numrhs_ = MVT::GetNumberVecs( *init_res );
836  scalevector_.resize( numrhs_ );
837  MVT::MvNorm( *init_res, scalevector_, scalenormtype_ );
838  }
839  else {
840  numrhs_ = MVT::GetNumberVecs( *(lp.getRHS()) );
841  }
842 
843  resvector_.resize( numrhs_ );
844  testvector_.resize( numrhs_ );
845 
846  curLSNum_ = lp.getLSNumber();
847  curLSIdx_ = lp.getLSIndex();
848  curBlksz_ = (int)curLSIdx_.size();
849  int validLS = 0;
850  for (i=0; i<curBlksz_; ++i) {
851  if (curLSIdx_[i] > -1 && curLSIdx_[i] < numrhs_)
852  validLS++;
853  }
854  curNumRHS_ = validLS;
855  //
856  // Initialize the testvector.
857  for (i=0; i<numrhs_; i++) { testvector_[i] = one; }
858 
859  // Return an error if the scaling is zero.
860  if (scalevalue_ == zero) {
861  return Failed;
862  }
863  }
864  return Undefined;
865 }
866 
867 } // end namespace Belos
868 
869 #endif /* BELOS_STATUS_TEST_IMPRESNORM_H */
virtual Teuchos::RCP< const MV > getNativeResiduals(std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > *norms) const =0
Get the residuals native to the solver.
Teuchos::RCP< MV > getSolution()
Returns the current solution estimate that was computed for the most recent residual test...
virtual ~StatusTestImpResNorm()
Destructor (virtual for memory safety).
const std::vector< MagnitudeType > * getTestValue() const
Returns the test value, , computed in most recent call to CheckStatus.
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:120
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
Convergence test using the implicit residual norm(s), with an explicit residual norm(s) check for los...
StatusType checkStatus(Iteration< ScalarType, MV, OP > *iSolver)
Check convergence status: Passed, Failed, or Undefined.
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.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual Teuchos::RCP< MV > getCurrentUpdate() const =0
Get the current update to the linear system.
Declaration of basic traits for the multivector type.
Teuchos::RCP< const MV > getInitPrecResVec() const
A pointer to the preconditioned initial residual vector.
An abstract class of StatusTest for stopping criteria using residual norms.
Teuchos::RCP< const MV > getRHS() const
A pointer to the right-hand side B.
StatusType
Whether the StatusTest wants iteration to stop.
Definition: BelosTypes.hpp:189
const std::vector< MagnitudeType > * getResNormValue() const
Returns the residual norm value, , computed in most recent call to CheckStatus.
Traits class which defines basic operations on multivectors.
int getQuorum() const
Returns the number of residuals that must pass the convergence test before Passed is returned...
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
The type of the magnitude (absolute value) of a ScalarType.
const std::vector< int > & getLSIndex() const
(Zero-based) indices of the linear system(s) currently being solved.
int setShowMaxResNormOnly(bool showMaxResNormOnly)
Set whether the only maximum residual norm is displayed when the print() method is called...
const std::vector< MagnitudeType > * getScaledNormValue() const
Returns the scaled norm value, .
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
int getLSNumber() const
The number of linear systems that have been set.
int setQuorum(int quorum)
Sets the number of residuals that must pass the convergence test before Passed is returned...
StatusType firstCallCheckStatusSetup(Iteration< ScalarType, MV, OP > *iSolver)
Call to setup initial scaling vector.
MagnitudeType getCurrTolerance() const
Current convergence tolerance; may be changed to prevent loss of accuracy.
int defineResForm(NormType TypeOfNorm)
Define form of the residual, its norm and optional weighting vector.
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const =0
Get a constant reference to the linear problem.
Teuchos::RCP< const MV > getInitResVec() const
A pointer to the initial unpreconditioned residual vector.
Belos::StatusTest abstract class for specifying a residual norm stopping criteria.
bool getLOADetected() const
Returns a boolean indicating a loss of accuracy has been detected in computing the residual...
NormType
The type of vector norm to compute.
Definition: BelosTypes.hpp:97
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. ...
StatusType getStatus() const
Return the result of the most recent CheckStatus call.
MagnitudeType getTolerance() const
&quot;Original&quot; convergence tolerance as set by user.
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, define an explicit value.
TypeTo as(const TypeFrom &t)
StatusTestImpResNorm(MagnitudeType Tolerance, int quorum=-1, bool showMaxResNormOnly=false)
Constructor.
int setTolerance(MagnitudeType tolerance)
Set the value of the tolerance.
std::string description() const
Method to return description of the maximum iteration status test.
void printStatus(std::ostream &os, StatusType type) const
Print message for each status specific to this stopping test.
void reset()
Resets the internal configuration to the initial state.
void print(std::ostream &os, int indent=0) const
Output formatted description of stopping test to output stream.
std::vector< int > convIndices()
Returns the vector containing the indices of the residuals that passed the test.
bool getShowMaxResNormOnly()
Returns whether the only maximum residual norm is displayed when the print() method is called...

Generated on Thu Nov 21 2024 09:28:19 for Belos by doxygen 1.8.5