45#ifndef THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DECL_HPP
46#define THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DECL_HPP
48#include "Thyra_SolveSupportTypes.hpp"
49#include "BelosStatusTest.hpp"
50#include "BelosLinearProblem.hpp"
51#include "BelosMultiVecTraits.hpp"
52#include "Teuchos_VerboseObject.hpp"
64 :
public Belos::StatusTest<Scalar, MultiVectorBase<Scalar>, LinearOpBase<Scalar> >,
65 public Teuchos::VerboseObject<GeneralSolveCriteriaBelosStatusTest<Scalar> >
72 typedef MultiVectorBase<Scalar>
MV;
74 typedef LinearOpBase<Scalar>
OP;
76 typedef typename ScalarTraits<Scalar>::magnitudeType
ScalarMag;
87 const int convergenceTestFrequency);
97 virtual Belos::StatusType
checkStatus(Belos::Iteration<Scalar,MV,OP> *iSolver);
99 virtual Belos::StatusType
getStatus()
const;
101 virtual void reset();
103 virtual void print(std::ostream& os,
int indent)
const;
125 const Ptr<
const ReductionFunctional<Scalar> > &reductFunc,
126 const Ptr<
const VectorBase<Scalar> > &x,
127 const Ptr<
const VectorBase<Scalar> > &r
130 void printRhsStatus(
const int currIter,
const int j, std::ostream &out,
131 int indent = 0)
const;
140template<
class Scalar>
141RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> >
143 const SolveCriteria<Scalar> &solveCriteria,
144 const int convergenceTestFrequency
147 RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> >
149 gscbst->setSolveCriteria(solveCriteria, convergenceTestFrequency);
Subclass of Belos::StatusTest that implements every possible form of SolveCriteria that exists by for...
GeneralSolveCriteriaBelosStatusTest()
void setSolveCriteria(const SolveCriteria< Scalar > &solveCriteria, const int convergenceTestFrequency)
void printRhsStatus(const int currIter, const int j, std::ostream &out, int indent=0) const
Array< ScalarMag > lastNumerator_
RCP< GeneralSolveCriteriaBelosStatusTest< Scalar > > createGeneralSolveCriteriaBelosStatusTest(const SolveCriteria< Scalar > &solveCriteria, const int convergenceTestFrequency)
Nonmember constructor.
MultiVectorBase< Scalar > MV
int convergenceTestFrequency_
Array< ScalarMag > b_nrm_
virtual Belos::StatusType getStatus() const
virtual Belos::StatusType checkStatus(Belos::Iteration< Scalar, MV, OP > *iSolver)
LinearOpBase< Scalar > OP
ScalarTraits< Scalar >::magnitudeType ScalarMag
ArrayView< const ScalarMag > achievedTol() const
Array< ScalarMag > r0_nrm_
SolveCriteria< Scalar > solveCriteria_
virtual void print(std::ostream &os, int indent) const
Belos::StatusType lastRtnStatus_
Array< ScalarMag > lastDenominator_
Array< ScalarMag > lastAchievedTol_
ScalarMag computeReductionFunctional(ESolveMeasureNormType measureType, const Ptr< const ReductionFunctional< Scalar > > &reductFunc, const Ptr< const VectorBase< Scalar > > &x, const Ptr< const VectorBase< Scalar > > &r) const