Stratimikos Version of the Day
Loading...
Searching...
No Matches
Thyra_GeneralSolveCriteriaBelosStatusTest_decl.hpp
1// @HEADER
2// *****************************************************************************
3// Stratimikos: Thyra-based strategies for linear solvers
4//
5// Copyright 2006 NTESS and the Stratimikos contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DECL_HPP
11#define THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DECL_HPP
12
13#include "Thyra_SolveSupportTypes.hpp"
14#include "BelosStatusTest.hpp"
17#include "Teuchos_VerboseObject.hpp"
18
19
20namespace Thyra {
21
22
27template<class Scalar>
29 : public Belos::StatusTest<Scalar, MultiVectorBase<Scalar>, LinearOpBase<Scalar> >,
30 public Teuchos::VerboseObject<GeneralSolveCriteriaBelosStatusTest<Scalar> >
31{
32public:
33
37 typedef MultiVectorBase<Scalar> MV;
39 typedef LinearOpBase<Scalar> OP;
41 typedef typename ScalarTraits<Scalar>::magnitudeType ScalarMag;
43
46
49
51 void setSolveCriteria(const SolveCriteria<Scalar> &solveCriteria,
52 const int convergenceTestFrequency);
53
55 ArrayView<const ScalarMag> achievedTol() const;
56
58
64 virtual Belos::StatusType getStatus() const;
66 virtual void reset();
68 virtual void print(std::ostream& os, int indent) const;
70
71private:
72
73 SolveCriteria<Scalar> solveCriteria_;
74 int convergenceTestFrequency_;
75
76 bool compute_x_;
77 bool compute_r_;
78
79 Array<ScalarMag> r0_nrm_;
80 Array<ScalarMag> b_nrm_;
81 Array<ScalarMag> lastNumerator_;
82 Array<ScalarMag> lastDenominator_;
83 Array<ScalarMag> lastAchievedTol_;
84 int lastCurrIter_;
85 Belos::StatusType lastRtnStatus_;
86
87 // Private member functions
88
89 ScalarMag computeReductionFunctional(ESolveMeasureNormType measureType,
90 const Ptr<const ReductionFunctional<Scalar> > &reductFunc,
91 const Ptr<const VectorBase<Scalar> > &x,
92 const Ptr<const VectorBase<Scalar> > &r
93 ) const;
94
95 void printRhsStatus(const int currIter, const int j, std::ostream &out,
96 int indent = 0) const;
97
98};
99
100
105template<class Scalar>
106RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> >
108 const SolveCriteria<Scalar> &solveCriteria,
109 const int convergenceTestFrequency
110 )
111{
112 RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> >
113 gscbst = Teuchos::rcp(new GeneralSolveCriteriaBelosStatusTest<Scalar>);
114 gscbst->setSolveCriteria(solveCriteria, convergenceTestFrequency);
115 return gscbst;
116}
117
118
119} // namespace Thyra
120
121
122#endif // THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DECL_HPP
Subclass of Belos::StatusTest that implements every possible form of SolveCriteria that exists by for...
void setSolveCriteria(const SolveCriteria< Scalar > &solveCriteria, const int convergenceTestFrequency)
RCP< GeneralSolveCriteriaBelosStatusTest< Scalar > > createGeneralSolveCriteriaBelosStatusTest(const SolveCriteria< Scalar > &solveCriteria, const int convergenceTestFrequency)
Nonmember constructor.
virtual Belos::StatusType checkStatus(Belos::Iteration< Scalar, MV, OP > *iSolver)

Generated for Stratimikos by doxygen 1.9.8