Stratimikos Version of the Day
Loading...
Searching...
No Matches
Thyra_GeneralSolveCriteriaBelosStatusTest_def.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_DEF_HPP
11#define THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
12
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"
19#include "BelosThyraAdapter.hpp"
20#include "Teuchos_DebugDefaultAsserts.hpp"
21#include "Teuchos_VerbosityLevel.hpp"
22#include "Teuchos_as.hpp"
23
24
25namespace Thyra {
26
27
28// Constructors/initializers/accessors
29
30
31template<class Scalar>
33 :convergenceTestFrequency_(-1),
34 compute_x_(false),
35 compute_r_(false),
36 lastCurrIter_(-1),
37 lastRtnStatus_(Belos::Undefined)
38{
40}
41
42
43template<class Scalar>
45 const SolveCriteria<Scalar> &solveCriteria,
46 const int convergenceTestFrequency
47 )
48{
49
50 using Teuchos::as;
51 typedef ScalarTraits<ScalarMag> SMT;
52
53 // Make sure the solve criteria is valid
54
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) );
59
60 // Remember the solve criteria sett
61
62 solveCriteria_ = solveCriteria;
63 convergenceTestFrequency_ = convergenceTestFrequency;
64
65 // Determine what needs to be computed on each check
66
67 compute_r_ = solveCriteria.solveMeasureType.contains(SOLVE_MEASURE_NORM_RESIDUAL);
68
69 compute_x_ = (compute_r_ ||
70 solveCriteria.solveMeasureType.contains(SOLVE_MEASURE_NORM_SOLUTION));
71
72}
73
74
75template<class Scalar>
76ArrayView<const typename ScalarTraits<Scalar>::magnitudeType>
78{
79 return lastAchievedTol_;
80}
81
82
83// Overridden public functions from Belos::StatusTest
84
85
86template <class Scalar>
90 )
91{
92
93 using Teuchos::null;
94
95#ifdef THYRA_DEBUG
96 TEUCHOS_ASSERT(iSolver);
97#endif
98
99 const int currIter = iSolver->getNumIters();
100
101 if (currIter == 0 || currIter % convergenceTestFrequency_ != 0) {
102 // We will skip this check!
103 return Belos::Undefined;
104 }
105
106 const RCP<FancyOStream> out = this->getOStream();
107 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
108
109 const Belos::LinearProblem<Scalar,MV,OP>& lp = iSolver->getProblem();
110 const int numRhs = lp.getRHS()->domain()->dim();
111
112 // Compute the rhs norm if requested and not already computed
113 if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_RHS)) {
114 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||b||");
115 }
116
117 // Compute the initial residual norm if requested and not already computed
118 if (solveCriteria_.solveMeasureType.contains(SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
119 TEUCHOS_TEST_FOR_EXCEPT_MSG(true, "ToDo: Handle ||r0||");
120 }
121
122 // Compute X if requested
123 RCP<MultiVectorBase<Scalar> > X;
124 if (compute_x_) {
125 RCP<MV> X_update = iSolver->getCurrentUpdate();
126 X = lp.updateSolution(X_update);
127 }
128
129 // Compute R if requested
130 RCP<MultiVectorBase<Scalar> > R;
131 if (compute_r_) {
132 R = createMembers(lp.getOperator()->range(), X->domain());
133 lp.computeCurrResVec(&*R, &*X);
134 }
135
136 // Form numerators and denominaotors gN(vN)/gD(vD) for each system
137
138 lastNumerator_.resize(numRhs);
139 lastDenominator_.resize(numRhs);
140
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() );
152 }
153
154 // Determine if convRatio <= requestedTol
155
156 bool systemsAreConverged = true;
157 lastAchievedTol_.resize(numRhs);
158
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);
165 }
166 if (!sys_converged_j) {
167 systemsAreConverged = false;
168 }
169 }
170
171 lastRtnStatus_ = (systemsAreConverged ? Belos::Passed : Belos::Failed);
172 lastCurrIter_ = currIter;
173
174 return lastRtnStatus_;
175
176}
177
178template <class Scalar>
181{
182 return lastRtnStatus_;
183}
184
185
186template <class Scalar>
188{
189 r0_nrm_.clear();
190 b_nrm_.clear();
191 lastNumerator_.clear();
192 lastDenominator_.clear();
193 lastAchievedTol_.clear();
194 lastCurrIter_ = -1;
195 lastRtnStatus_ = Belos::Undefined;
196}
197
198
199template <class Scalar>
201 std::ostream& os, int indent
202 ) const
203{
204 const int numRhs = lastNumerator_.size();
205 for (int j = 0; j < numRhs; ++j) {
206 printRhsStatus(lastCurrIter_, j, os, indent);
207 }
208}
209
210
211// private
212
213
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
221 ) const
222{
223 typedef ScalarTraits<ScalarMag> SMT;
224 ScalarMag rtn = -SMT::one();
225 Ptr<const VectorBase<Scalar> > v;
226 switch(measureType) {
227 case SOLVE_MEASURE_ONE:
228 rtn = SMT::one();
229 break;
230 case SOLVE_MEASURE_NORM_RESIDUAL:
231 v = r;
232 break;
233 case SOLVE_MEASURE_NORM_SOLUTION:
234 v = x;
235 break;
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();
241 }
242 if (rtn >= SMT::zero()) {
243 // We already have what we need!
244 }
245 else if (nonnull(v) && rtn < SMT::zero()) {
246 if (nonnull(reductFunc)) {
247 rtn = reductFunc->reduce(*v);
248 }
249 else {
250 rtn = norm(*v);
251 }
252 }
253 TEUCHOS_IF_ELSE_DEBUG_ASSERT();
254 return rtn;
255}
256
257
258template <class Scalar>
259void
260GeneralSolveCriteriaBelosStatusTest<Scalar>::printRhsStatus(
261 const int currIter, const int j, std::ostream &out,
262 int indent
263 ) const
264{
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 << " "; }
268 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")
274 << "\n";
275}
276
277
278} // namespace Thyra
279
280
281#endif // THYRA_GENERAL_SOLVE_CRITERIA_BELOS_STATUS_TEST_DEF_HPP
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...
void setSolveCriteria(const SolveCriteria< Scalar > &solveCriteria, const int convergenceTestFrequency)
virtual Belos::StatusType checkStatus(Belos::Iteration< Scalar, MV, OP > *iSolver)
Undefined

Generated for Stratimikos by doxygen 1.9.8