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),
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>
134 const int currIter = iSolver->getNumIters();
136 if (currIter == 0 || currIter % convergenceTestFrequency_ != 0) {
145 const int numRhs = lp.
getRHS()->domain()->dim();
160 RCP<MV> X_update = iSolver->getCurrentUpdate();
167 R = createMembers(lp.
getOperator()->range(), X->domain());
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;
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();
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
GeneralSolveCriteriaBelosStatusTest()
Subclass of Belos::StatusTest that implements every possible form of SolveCriteria that exists by for...
ArrayView< const ScalarMag > achievedTol() const
#define TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT()
bool contains(ESolveMeasureNormType measure) const
virtual Teuchos::RCP< MV > updateSolution(const Teuchos::RCP< MV > &update=Teuchos::null, bool updateLP=false, ScalarType scale=Teuchos::ScalarTraits< ScalarType >::one())
virtual Belos::StatusType getStatus() const
#define TEUCHOS_ASSERT_INEQUALITY(val1, comp, val2)
Teuchos::RCP< const MV > getRHS() const
Thyra specializations of MultiVecTraits and OperatorTraits.
virtual Belos::StatusType checkStatus(Belos::Iteration< Scalar, MV, OP > *iSolver)
ESolveMeasureNormType numerator
Teuchos::RCP< const OP > getOperator() const
#define TEUCHOS_TEST_FOR_EXCEPT_MSG(throw_exception_test, msg)
#define TEUCHOS_IF_ELSE_DEBUG_ASSERT()
ScalarTraits< Scalar >::magnitudeType ScalarMag
TEUCHOSCORE_LIB_DLL_EXPORT bool includesVerbLevel(const EVerbosityLevel verbLevel, const EVerbosityLevel requestedVerbLevel, const bool isDefaultLevel=false)
bool nonnull(const boost::shared_ptr< T > &p)
virtual void computeCurrResVec(MV *R, const MV *X=0, const MV *B=0) const
RCP< const ReductionFunctional< Scalar > > numeratorReductionFunc
void setSolveCriteria(const SolveCriteria< Scalar > &solveCriteria, const int convergenceTestFrequency)
TypeTo as(const TypeFrom &t)
RCP< const ReductionFunctional< Scalar > > denominatorReductionFunc
SolveMeasureType solveMeasureType
#define TEUCHOS_ASSERT(assertion_test)
virtual void print(std::ostream &os, int indent) const