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" 
   20#include "Teuchos_DebugDefaultAsserts.hpp" 
   21#include "Teuchos_VerbosityLevel.hpp" 
   22#include "Teuchos_as.hpp" 
   33  :convergenceTestFrequency_(-1),
 
 
   45  const SolveCriteria<Scalar> &solveCriteria,
 
   46  const int convergenceTestFrequency
 
   51  typedef ScalarTraits<ScalarMag> SMT;
 
   55  TEUCHOS_ASSERT_INEQUALITY(solveCriteria.requestedTol, >=, SMT::zero());
 
   56  TEUCHOS_ASSERT_INEQUALITY(solveCriteria.solveMeasureType.numerator, !=, SOLVE_MEASURE_ONE);
 
   57  TEUCHOS_ASSERT(nonnull(solveCriteria.numeratorReductionFunc) ||
 
   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));
 
 
   76ArrayView<const typename ScalarTraits<Scalar>::magnitudeType>
 
   79  return lastAchievedTol_;
 
 
   86template <
class Scalar>
 
   96  TEUCHOS_ASSERT(iSolver);
 
  101  if (currIter == 0 || currIter % convergenceTestFrequency_ != 0) {
 
  106  const RCP<FancyOStream> out = this->getOStream();
 
  107  const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
 
  110  const int numRhs = lp.
getRHS()->domain()->dim();
 
  113  if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_RHS)) {
 
  114    TEUCHOS_TEST_FOR_EXCEPT_MSG(
true, 
"ToDo: Handle ||b||");
 
  118  if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
 
  119    TEUCHOS_TEST_FOR_EXCEPT_MSG(
true, 
"ToDo: Handle ||r0||");
 
  123  RCP<MultiVectorBase<Scalar> > X;
 
  130  RCP<MultiVectorBase<Scalar> > R;
 
  132    R = createMembers(lp.
getOperator()->range(), X->domain());
 
  138  lastNumerator_.resize(numRhs);
 
  139  lastDenominator_.resize(numRhs);
 
  141  for (
int j = 0; j < numRhs; ++j) {
 
  142    const RCP<const VectorBase<Scalar> > x_j = (nonnull(X) ? X->col(j) : null);
 
  143    const RCP<const VectorBase<Scalar> > r_j = (nonnull(R) ? R->col(j) : null);
 
  144    lastNumerator_[j] = computeReductionFunctional(
 
  145      solveCriteria_.solveMeasureType.numerator,
 
  146      solveCriteria_.numeratorReductionFunc.ptr(),
 
  147      x_j.ptr(), r_j.ptr() );
 
  148    lastDenominator_[j] = computeReductionFunctional(
 
  149      solveCriteria_.solveMeasureType.denominator,
 
  150      solveCriteria_.denominatorReductionFunc.ptr(),
 
  151      x_j.ptr(), r_j.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);
 
  163    if (includesVerbLevel(verbLevel, Teuchos::VERB_MEDIUM)) {
 
  164      printRhsStatus(currIter, j, *out);
 
  166    if (!sys_converged_j) {
 
  167      systemsAreConverged = 
false;
 
  172  lastCurrIter_ = currIter;
 
  174  return lastRtnStatus_;
 
 
  178template <
class Scalar>
 
  182  return lastRtnStatus_;
 
 
  186template <
class Scalar>
 
  191  lastNumerator_.clear();
 
  192  lastDenominator_.clear();
 
  193  lastAchievedTol_.clear();
 
 
  199template <
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);
 
 
  214template <
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;
 
  224  ScalarMag rtn = -SMT::one();
 
  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:
 
  237      TEUCHOS_TEST_FOR_EXCEPT_MSG(
true, 
"ToDo: Handle ||r0||!)")
 
  238    case SOLVE_MEASURE_NORM_RHS:
 
  239      TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||b||!)");
 
  240    TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT();
 
  242  if (rtn >= SMT::zero()) {
 
  245  else if (nonnull(v) && rtn < SMT::zero()) {
 
  246    if (nonnull(reductFunc)) {
 
  247      rtn = reductFunc->reduce(*v);
 
  253  TEUCHOS_IF_ELSE_DEBUG_ASSERT();
 
  258template <
class Scalar>
 
  260GeneralSolveCriteriaBelosStatusTest<Scalar>::printRhsStatus(
 
  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")
 
Thyra specializations of MultiVecTraits and OperatorTraits.
virtual int getNumIters() const=0
virtual Teuchos::RCP< MV > getCurrentUpdate() const=0
virtual const LinearProblem< ScalarType, MV, OP > & getProblem() const=0
Teuchos::RCP< const OP > getOperator() const
Teuchos::RCP< const MV > getRHS() const
virtual void computeCurrResVec(MV *R, const MV *X=0, const MV *B=0) const
virtual Teuchos::RCP< MV > updateSolution(const Teuchos::RCP< MV > &update=Teuchos::null, bool updateLP=false, ScalarType scale=Teuchos::ScalarTraits< ScalarType >::one())
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)
virtual Belos::StatusType getStatus() const
virtual Belos::StatusType checkStatus(Belos::Iteration< Scalar, MV, OP > *iSolver)
ScalarTraits< Scalar >::magnitudeType ScalarMag
ArrayView< const ScalarMag > achievedTol() const
virtual void print(std::ostream &os, int indent) const