10 #ifndef THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
11 #define THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
13 #include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp"
14 #include "Thyra_VectorSpaceBase.hpp"
15 #include "Thyra_MultiVectorBase.hpp"
16 #include "Thyra_VectorBase.hpp"
17 #include "Thyra_MultiVectorStdOps.hpp"
18 #include "Thyra_VectorStdOps.hpp"
31 template<
class Scalar>
33 :convergenceTestFrequency_(-1),
37 lastRtnStatus_(Belos::Undefined)
43 template<
class Scalar>
45 const SolveCriteria<Scalar> &solveCriteria,
46 const int convergenceTestFrequency
51 typedef ScalarTraits<ScalarMag> SMT;
58 nonnull(solveCriteria.denominatorReductionFunc) );
62 solveCriteria_ = solveCriteria;
63 convergenceTestFrequency_ = convergenceTestFrequency;
67 compute_r_ = solveCriteria.solveMeasureType.contains(SOLVE_MEASURE_NORM_RESIDUAL);
69 compute_x_ = (compute_r_ ||
70 solveCriteria.solveMeasureType.contains(SOLVE_MEASURE_NORM_SOLUTION));
75 template<
class Scalar>
76 ArrayView<const typename ScalarTraits<Scalar>::magnitudeType>
79 return lastAchievedTol_;
86 template <
class Scalar>
89 Belos::Iteration<Scalar,MV,OP> *iSolver
99 const int currIter = iSolver->getNumIters();
101 if (currIter == 0 || currIter % convergenceTestFrequency_ != 0) {
103 return Belos::Undefined;
109 const Belos::LinearProblem<Scalar,MV,OP>& lp = iSolver->getProblem();
110 const int numRhs = lp.getRHS()->domain()->dim();
113 if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_RHS)) {
118 if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
125 RCP<MV> X_update = iSolver->getCurrentUpdate();
126 X = lp.updateSolution(X_update);
132 R = createMembers(lp.getOperator()->range(), X->domain());
133 lp.computeCurrResVec(&*R, &*X);
138 lastNumerator_.resize(numRhs);
139 lastDenominator_.resize(numRhs);
141 for (
int j = 0; j < numRhs; ++j) {
144 lastNumerator_[j] = computeReductionFunctional(
145 solveCriteria_.solveMeasureType.numerator,
146 solveCriteria_.numeratorReductionFunc.ptr(),
148 lastDenominator_[j] = computeReductionFunctional(
149 solveCriteria_.solveMeasureType.denominator,
150 solveCriteria_.denominatorReductionFunc.ptr(),
156 bool systemsAreConverged =
true;
157 lastAchievedTol_.resize(numRhs);
159 for (
int j = 0; j < numRhs; ++j) {
160 const ScalarMag convRatio = lastNumerator_[j] / lastDenominator_[j];
161 lastAchievedTol_[j] = convRatio;
162 const bool sys_converged_j = (convRatio <= solveCriteria_.requestedTol);
164 printRhsStatus(currIter, j, *out);
166 if (!sys_converged_j) {
167 systemsAreConverged =
false;
171 lastRtnStatus_ = (systemsAreConverged ? Belos::Passed : Belos::Failed);
172 lastCurrIter_ = currIter;
174 return lastRtnStatus_;
178 template <
class Scalar>
182 return lastRtnStatus_;
186 template <
class Scalar>
191 lastNumerator_.clear();
192 lastDenominator_.clear();
193 lastAchievedTol_.clear();
195 lastRtnStatus_ = Belos::Undefined;
199 template <
class Scalar>
201 std::ostream& os,
int indent
204 const int numRhs = lastNumerator_.size();
205 for (
int j = 0; j < numRhs; ++j) {
206 printRhsStatus(lastCurrIter_, j, os, indent);
214 template <
class Scalar>
217 ESolveMeasureNormType measureType,
218 const Ptr<
const ReductionFunctional<Scalar> > &reductFunc,
219 const Ptr<
const VectorBase<Scalar> > &x,
220 const Ptr<
const VectorBase<Scalar> > &r
223 typedef ScalarTraits<ScalarMag> SMT;
225 Ptr<const VectorBase<Scalar> > v;
226 switch(measureType) {
227 case SOLVE_MEASURE_ONE:
230 case SOLVE_MEASURE_NORM_RESIDUAL:
233 case SOLVE_MEASURE_NORM_SOLUTION:
236 case SOLVE_MEASURE_NORM_INIT_RESIDUAL:
238 case SOLVE_MEASURE_NORM_RHS:
242 if (rtn >= SMT::zero()) {
245 else if (
nonnull(v) && rtn < SMT::zero()) {
247 rtn = reductFunc->reduce(*v);
258 template <
class Scalar>
261 const int currIter,
const int j, std::ostream &out,
265 const ScalarMag convRatio = lastNumerator_[j] / lastDenominator_[j];
266 const bool sys_converged_j = (convRatio <= solveCriteria_.requestedTol);
267 for (
int i = 0; i < indent; ++i) { out <<
" "; }
269 <<
"["<<currIter<<
"] "
270 <<
"gN(vN("<<j<<
"))/gD(vD("<<j<<
")) = "
271 << lastNumerator_[j] <<
"/" << lastDenominator_[j] <<
" = "
272 << convRatio <<
" <= " << solveCriteria_.requestedTol <<
" : "
273 << (sys_converged_j ?
" true" :
"false")
281 #endif // THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
GeneralSolveCriteriaBelosStatusTest()
ArrayView< const ScalarMag > achievedTol() const
virtual Belos::StatusType getStatus() const
#define TEUCHOS_ASSERT_INEQUALITY(val1, comp, val2)
void printRhsStatus(const int currIter, const int j, std::ostream &out, int indent=0) const
Thyra specializations of MultiVecTraits and OperatorTraits.
virtual Belos::StatusType checkStatus(Belos::Iteration< Scalar, MV, OP > *iSolver)
#define TEUCHOS_IF_ELSE_DEBUG_ASSERT()
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
#define TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT()
ScalarTraits< Scalar >::magnitudeType ScalarMag
ScalarMag computeReductionFunctional(ESolveMeasureNormType measureType, const Ptr< const ReductionFunctional< Scalar > > &reductFunc, const Ptr< const VectorBase< Scalar > > &x, const Ptr< const VectorBase< Scalar > > &r) const
TEUCHOSCORE_LIB_DLL_EXPORT bool includesVerbLevel(const EVerbosityLevel verbLevel, const EVerbosityLevel requestedVerbLevel, const bool isDefaultLevel=false)
TypeTo as(const TypeFrom &t)
bool nonnull(const boost::shared_ptr< T > &p)
void setSolveCriteria(const SolveCriteria< Scalar > &solveCriteria, const int convergenceTestFrequency)
#define TEUCHOS_ASSERT(assertion_test)
virtual void print(std::ostream &os, int indent) const