45 #ifndef THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
46 #define THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
48 #include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp"
49 #include "Thyra_VectorSpaceBase.hpp"
50 #include "Thyra_MultiVectorBase.hpp"
51 #include "Thyra_VectorBase.hpp"
52 #include "Thyra_MultiVectorStdOps.hpp"
53 #include "Thyra_VectorStdOps.hpp"
55 #include "Teuchos_DebugDefaultAsserts.hpp"
56 #include "Teuchos_VerbosityLevel.hpp"
57 #include "Teuchos_as.hpp"
66 template<
class Scalar>
68 :convergenceTestFrequency_(-1),
72 lastRtnStatus_(Belos::Undefined)
78 template<
class Scalar>
81 const int convergenceTestFrequency
86 typedef ScalarTraits<ScalarMag> SMT;
97 solveCriteria_ = solveCriteria;
98 convergenceTestFrequency_ = convergenceTestFrequency;
104 compute_x_ = (compute_r_ ||
110 template<
class Scalar>
111 ArrayView<const typename ScalarTraits<Scalar>::magnitudeType>
114 return lastAchievedTol_;
121 template <
class Scalar>
124 Belos::Iteration<Scalar,MV,OP> *iSolver
134 const int currIter = iSolver->getNumIters();
136 if (currIter == 0 || currIter % convergenceTestFrequency_ != 0) {
138 return Belos::Undefined;
144 const Belos::LinearProblem<Scalar,MV,OP>& lp = iSolver->getProblem();
145 const int numRhs = lp.getRHS()->domain()->dim();
160 RCP<MV> X_update = iSolver->getCurrentUpdate();
161 X = lp.updateSolution(X_update);
167 R = createMembers(lp.getOperator()->range(), X->domain());
168 lp.computeCurrResVec(&*R, &*X);
173 lastNumerator_.resize(numRhs);
174 lastDenominator_.resize(numRhs);
176 for (
int j = 0; j < numRhs; ++j) {
179 lastNumerator_[j] = computeReductionFunctional(
180 solveCriteria_.solveMeasureType.numerator,
181 solveCriteria_.numeratorReductionFunc.ptr(),
183 lastDenominator_[j] = computeReductionFunctional(
184 solveCriteria_.solveMeasureType.denominator,
185 solveCriteria_.denominatorReductionFunc.ptr(),
191 bool systemsAreConverged =
true;
192 lastAchievedTol_.resize(numRhs);
194 for (
int j = 0; j < numRhs; ++j) {
195 const ScalarMag convRatio = lastNumerator_[j] / lastDenominator_[j];
196 lastAchievedTol_[j] = convRatio;
197 const bool sys_converged_j = (convRatio <= solveCriteria_.requestedTol);
199 printRhsStatus(currIter, j, *out);
201 if (!sys_converged_j) {
202 systemsAreConverged =
false;
206 lastRtnStatus_ = (systemsAreConverged ? Belos::Passed : Belos::Failed);
207 lastCurrIter_ = currIter;
209 return lastRtnStatus_;
213 template <
class Scalar>
217 return lastRtnStatus_;
221 template <
class Scalar>
226 lastNumerator_.clear();
227 lastDenominator_.clear();
228 lastAchievedTol_.clear();
230 lastRtnStatus_ = Belos::Undefined;
234 template <
class Scalar>
236 std::ostream& os,
int indent
239 const int numRhs = lastNumerator_.size();
240 for (
int j = 0; j < numRhs; ++j) {
241 printRhsStatus(lastCurrIter_, j, os, indent);
249 template <
class Scalar>
258 typedef ScalarTraits<ScalarMag> SMT;
259 ScalarMag rtn = -SMT::one();
260 Ptr<const VectorBase<Scalar> > v;
261 switch(measureType) {
274 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||b||!)");
277 if (rtn >= SMT::zero()) {
280 else if (
nonnull(v) && rtn < SMT::zero()) {
282 rtn = reductFunc->reduce(*v);
293 template <
class Scalar>
295 GeneralSolveCriteriaBelosStatusTest<Scalar>::printRhsStatus(
296 const int currIter,
const int j, std::ostream &out,
300 const ScalarMag convRatio = lastNumerator_[j] / lastDenominator_[j];
301 const bool sys_converged_j = (convRatio <= solveCriteria_.requestedTol);
302 for (
int i = 0; i < indent; ++i) { out <<
" "; }
304 <<
"["<<currIter<<
"] "
305 <<
"gN(vN("<<j<<
"))/gD(vD("<<j<<
")) = "
306 << lastNumerator_[j] <<
"/" << lastDenominator_[j] <<
" = "
307 << convRatio <<
" <= " << solveCriteria_.requestedTol <<
" : "
308 << (sys_converged_j ?
" true" :
"false")
316 #endif // THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP