Stratimikos Version of the Day
Loading...
Searching...
No Matches
Thyra_BelosLinearOpWithSolve_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_BELOS_LINEAR_OP_WITH_SOLVE_HPP
11#define THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
12
13#include "Thyra_BelosLinearOpWithSolve_decl.hpp"
14#include "Thyra_GeneralSolveCriteriaBelosStatusTest.hpp"
15#include "Thyra_LinearOpWithSolveHelpers.hpp"
16#include "Teuchos_DebugDefaultAsserts.hpp"
17#include "Teuchos_Assert.hpp"
18#include "Teuchos_TimeMonitor.hpp"
19#include "Stratimikos_Config.h"
20#ifdef HAVE_STRATIMIKOS_THYRATPETRAADAPTERS
21# include "Thyra_TpetraThyraWrappers.hpp"
22# include <MatrixMarket_Tpetra.hpp>
23# include <TpetraExt_MatrixMatrix.hpp>
24#endif
25
26namespace {
27 // Set the Belos solver's parameter list to scale its residual norms
28 // in the specified way.
29 //
30 // We break this out in a separate function because the parameters
31 // to set depend on which parameters the Belos solver supports. Not
32 // all Belos solvers support both the "Implicit Residual Scaling"
33 // and "Explicit Residual Scaling" parameters, so we have to check
34 // the solver's list of valid parameters for the existence of these.
35 //
36 // Scaling options: Belos lets you decide whether the solver will
37 // scale residual norms by the (left-)preconditioned initial
38 // residual norms (residualScalingType = "Norm of Initial
39 // Residual"), or by the unpreconditioned initial residual norms
40 // (residualScalingType = "Norm of RHS"). Usually you want to scale
41 // by the unpreconditioned initial residual norms. This is because
42 // preconditioning is just an optimization, and you really want to
43 // make ||B - AX|| small, rather than ||M B - M (A X)||. If you're
44 // measuring ||B - AX|| and scaling by the initial residual, you
45 // should use the unpreconditioned initial residual to match it.
46 //
47 // Note, however, that the implicit residual test computes
48 // left-preconditioned residuals, if a left preconditioner was
49 // provided. That's OK because when Belos solvers (at least the
50 // GMRES variants) are given a left preconditioner, they first check
51 // the implicit residuals. If those converge, they then check the
52 // explicit residuals. The explicit residual test does _not_ apply
53 // the left preconditioner when computing the residual. The
54 // implicit residual test is just an optimization so that Belos
55 // doesn't have to compute explicit residuals B - A*X at every
56 // iteration. This is why we use the same scaling factor for both
57 // the implicit and explicit residuals.
58 //
59 // Arguments:
60 //
61 // solverParams [in/out] Parameters for the current solve.
62 //
63 // solverValidParams [in] Valid parameter list for the Belos solver.
64 // Result of calling the solver's getValidParameters() method.
65 //
66 // residualScalingType [in] String describing how the solver should
67 // scale residuals. Valid values include "Norm of RHS" and "Norm
68 // of Initial Residual" (these are the only two options this file
69 // currently uses, though Belos offers other options).
70 void
71 setResidualScalingType (const Teuchos::RCP<Teuchos::ParameterList>& solverParams,
72 const Teuchos::RCP<const Teuchos::ParameterList>& solverValidParams,
73 const std::string& residualScalingType)
74 {
75 // Many Belos solvers (especially the GMRES variants) define both
76 // "Implicit Residual Scaling" and "Explicit Residual Scaling"
77 // options.
78 //
79 // "Implicit" means "the left-preconditioned approximate
80 // a.k.a. 'recursive' residual as computed by the Krylov method."
81 //
82 // "Explicit" means ||B - A*X||, the unpreconditioned, "exact"
83 // residual.
84 //
85 // Belos' GMRES implementations chain these two tests in sequence.
86 // Implicit comes first, and explicit is not evaluated unless
87 // implicit passes. In some cases (e.g., no left preconditioner),
88 // GMRES _only_ uses the implicit tests. This means that only
89 // setting "Explicit Residual Scaling" won't change the solver's
90 // behavior. Stratimikos tends to prefer using a right
91 // preconditioner, in which case setting only the "Explicit
92 // Residual Scaling" argument has no effect. Furthermore, if
93 // "Explicit Residual Scaling" is set to something other than the
94 // default (initial residual norm), without "Implicit Residual
95 // Scaling" getting the same setting, then the implicit residual
96 // test will be using a radically different scaling factor than
97 // the user wanted.
98 //
99 // Not all Belos solvers support both options. We check the
100 // solver's valid parameter list first before attempting to set
101 // the option.
102 if (solverValidParams->isParameter ("Implicit Residual Scaling")) {
103 solverParams->set ("Implicit Residual Scaling", residualScalingType);
104 }
105 if (solverValidParams->isParameter ("Explicit Residual Scaling")) {
106 solverParams->set ("Explicit Residual Scaling", residualScalingType);
107 }
108 }
109
110} // namespace (anonymous)
111
112
113namespace Thyra {
114
115
116// Constructors/initializers/accessors
117
118
119template<class Scalar>
121 :convergenceTestFrequency_(-1),
122 isExternalPrec_(false),
123 supportSolveUse_(SUPPORT_SOLVE_UNSPECIFIED),
124 defaultTol_ (-1.0),
125 label_(""),
126 filenameLHS_(""),
127 filenameRHS_(""),
128 init_(false),
129 counter_(0)
130{}
131
132
133template<class Scalar>
136 const RCP<Teuchos::ParameterList> &solverPL,
137 const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver,
138 const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
139 const RCP<const PreconditionerBase<Scalar> > &prec,
140 const bool isExternalPrec_in,
141 const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
142 const ESupportSolveUse &supportSolveUse_in,
143 const int convergenceTestFrequency
144 )
145{
146 using Teuchos::as;
147 using Teuchos::TypeNameTraits;
148 using Teuchos::Exceptions::InvalidParameterType;
149 typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type;
150
151 this->setLinePrefix("BELOS/T");
152 // ToDo: Validate input
153 lp_ = lp;
154 solverPL_ = solverPL;
155 iterativeSolver_ = iterativeSolver;
156 fwdOpSrc_ = fwdOpSrc;
157 prec_ = prec;
158 isExternalPrec_ = isExternalPrec_in;
159 approxFwdOpSrc_ = approxFwdOpSrc;
160 supportSolveUse_ = supportSolveUse_in;
161 convergenceTestFrequency_ = convergenceTestFrequency;
162 init_ = false;
163 // Check if "Convergence Tolerance" is in the solver parameter list. If
164 // not, use the default from the solver.
165 if ( !is_null(solverPL_) ) {
166 if (solverPL_->isParameter("Convergence Tolerance")) {
167
168 // Stratimikos prefers tolerances as double, no matter the
169 // Scalar type. However, we also want it to accept the
170 // tolerance as magnitude_type, for example float if the Scalar
171 // type is float or std::complex<float>.
172 if (solverPL_->isType<double> ("Convergence Tolerance")) {
173 defaultTol_ =
174 as<magnitude_type> (solverPL_->get<double> ("Convergence Tolerance"));
175 }
176 else if (std::is_same_v<double, magnitude_type>) {
177 // magnitude_type == double in this case, and we've already
178 // checked double above.
179 TEUCHOS_TEST_FOR_EXCEPTION(
180 true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: "
181 "The \"Convergence Tolerance\" parameter, which you provided, must "
182 "have type double (the type of the magnitude of Scalar = double).");
183 }
184 else if (solverPL_->isType<magnitude_type> ("Convergence Tolerance")) {
185 defaultTol_ = solverPL_->get<magnitude_type> ("Convergence Tolerance");
186 }
187 else {
188 // Throwing InvalidParameterType ensures that the exception's
189 // type is consistent both with what this method would have
190 // thrown before for an unrecognized type, and with what the
191 // user expects in general when the parameter doesn't have the
192 // right type.
193 TEUCHOS_TEST_FOR_EXCEPTION(
194 true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: "
195 "The \"Convergence Tolerance\" parameter, which you provided, must "
196 "have type double (preferred) or the type of the magnitude of Scalar "
197 "= " << TypeNameTraits<Scalar>::name () << ", which is " <<
198 TypeNameTraits<magnitude_type>::name () << " in this case. You can "
199 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
200 }
201 }
202
203 if (solverPL_->isParameter("Timer Label") && solverPL_->isType<std::string>("Timer Label")) {
204 label_ = solverPL_->get<std::string>("Timer Label");
205 lp_->setLabel(label_);
206 }
207 if (solverPL_->isParameter("Filename LHS") && solverPL_->isType<std::string>("Filename LHS")) {
208 filenameLHS_ = solverPL_->get<std::string>("Filename LHS");
209 }
210 if (solverPL_->isParameter("Filename RHS") && solverPL_->isType<std::string>("Filename RHS")) {
211 filenameRHS_ = solverPL_->get<std::string>("Filename RHS");
212 }
213 }
214 else {
215 RCP<const Teuchos::ParameterList> defaultPL =
216 iterativeSolver->getValidParameters();
217
218 // Stratimikos prefers tolerances as double, no matter the
219 // Scalar type. However, we also want it to accept the
220 // tolerance as magnitude_type, for example float if the Scalar
221 // type is float or std::complex<float>.
222 if (defaultPL->isType<double> ("Convergence Tolerance")) {
223 defaultTol_ =
224 as<magnitude_type> (defaultPL->get<double> ("Convergence Tolerance"));
225 }
226 else if (std::is_same_v<double, magnitude_type>) {
227 // magnitude_type == double in this case, and we've already
228 // checked double above.
229 TEUCHOS_TEST_FOR_EXCEPTION(
230 true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: "
231 "The \"Convergence Tolerance\" parameter, which you provided, must "
232 "have type double (the type of the magnitude of Scalar = double).");
233 }
234 else if (defaultPL->isType<magnitude_type> ("Convergence Tolerance")) {
235 defaultTol_ = defaultPL->get<magnitude_type> ("Convergence Tolerance");
236 }
237 else {
238 // Throwing InvalidParameterType ensures that the exception's
239 // type is consistent both with what this method would have
240 // thrown before for an unrecognized type, and with what the
241 // user expects in general when the parameter doesn't have the
242 // right type.
243 TEUCHOS_TEST_FOR_EXCEPTION(
244 true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: "
245 "The \"Convergence Tolerance\" parameter, which you provided, must "
246 "have type double (preferred) or the type of the magnitude of Scalar "
247 "= " << TypeNameTraits<Scalar>::name () << ", which is " <<
248 TypeNameTraits<magnitude_type>::name () << " in this case. You can "
249 "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
250 }
251 }
252}
253
254
255template<class Scalar>
256RCP<const LinearOpSourceBase<Scalar> >
258{
259 RCP<const LinearOpSourceBase<Scalar> >
260 _fwdOpSrc = fwdOpSrc_;
261 fwdOpSrc_ = Teuchos::null;
262 return _fwdOpSrc;
263}
264
265
266template<class Scalar>
267RCP<const PreconditionerBase<Scalar> >
269{
270 RCP<const PreconditionerBase<Scalar> >
271 _prec = prec_;
272 prec_ = Teuchos::null;
273 return _prec;
274}
275
276
277template<class Scalar>
279{
280 return isExternalPrec_;
281}
282
283
284template<class Scalar>
285RCP<const LinearOpSourceBase<Scalar> >
287{
288 RCP<const LinearOpSourceBase<Scalar> >
289 _approxFwdOpSrc = approxFwdOpSrc_;
290 approxFwdOpSrc_ = Teuchos::null;
291 return _approxFwdOpSrc;
292}
293
294
295template<class Scalar>
297{
298 return supportSolveUse_;
299}
300
301
302template<class Scalar>
305 RCP<Teuchos::ParameterList> *solverPL,
306 RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver,
307 RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc,
308 RCP<const PreconditionerBase<Scalar> > *prec,
309 bool *isExternalPrec_in,
310 RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
311 ESupportSolveUse *supportSolveUse_in
312 )
313{
314 if (lp) *lp = lp_;
315 if (solverPL) *solverPL = solverPL_;
316 if (iterativeSolver) *iterativeSolver = iterativeSolver_;
317 if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
318 if (prec) *prec = prec_;
319 if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_;
320 if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
321 if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_;
322
323 lp_ = Teuchos::null;
324 solverPL_ = Teuchos::null;
325 iterativeSolver_ = Teuchos::null;
326 fwdOpSrc_ = Teuchos::null;
327 prec_ = Teuchos::null;
328 isExternalPrec_ = false;
329 approxFwdOpSrc_ = Teuchos::null;
330 supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED;
331 init_ = false;
332}
333
334
335// Overridden from LinearOpBase
336
337
338template<class Scalar>
339RCP< const VectorSpaceBase<Scalar> >
341{
342 if (!is_null(lp_))
343 return lp_->getOperator()->range();
344 return Teuchos::null;
345}
346
347
348template<class Scalar>
349RCP< const VectorSpaceBase<Scalar> >
351{
352 if (!is_null(lp_))
353 return lp_->getOperator()->domain();
354 return Teuchos::null;
355}
356
357
358template<class Scalar>
359RCP<const LinearOpBase<Scalar> >
361{
362 return Teuchos::null; // Not supported yet but could be
363}
364
365
366// Overridden from Teuchos::Describable
367
368
369template<class Scalar>
371{
372 std::ostringstream oss;
373 oss << Teuchos::Describable::description();
374 if ( !is_null(lp_) && !is_null(lp_->getOperator()) ) {
375 oss << "{";
376 oss << "iterativeSolver=\'"<<iterativeSolver_->description()<<"\'";
377 oss << ",fwdOp=\'"<<lp_->getOperator()->description()<<"\'";
378 if (lp_->getLeftPrec().get())
379 oss << ",leftPrecOp=\'"<<lp_->getLeftPrec()->description()<<"\'";
380 if (lp_->getRightPrec().get())
381 oss << ",rightPrecOp=\'"<<lp_->getRightPrec()->description()<<"\'";
382 oss << "}";
383 }
384 // ToDo: Make Belos::SolverManager derive from Teuchos::Describable so
385 // that we can get better information.
386 return oss.str();
387}
388
389
390template<class Scalar>
392 Teuchos::FancyOStream &out_arg,
393 const Teuchos::EVerbosityLevel verbLevel
394 ) const
395{
396 using Teuchos::FancyOStream;
397 using Teuchos::OSTab;
398 using Teuchos::describe;
399 RCP<FancyOStream> out = rcp(&out_arg,false);
400 OSTab tab(out);
401 switch (verbLevel) {
402 case Teuchos::VERB_DEFAULT: // fall-through
403 case Teuchos::VERB_LOW:
404 case Teuchos::VERB_MEDIUM: // fall-through
405 *out << this->description() << std::endl;
406 break;
407 case Teuchos::VERB_HIGH: // fall-through
408 case Teuchos::VERB_EXTREME:
409 {
410 *out
411 << Teuchos::Describable::description()<< "{"
412 << "rangeDim=" << this->range()->dim()
413 << ",domainDim=" << this->domain()->dim() << "}\n";
414 if (lp_->getOperator().get()) {
415 OSTab tab1(out);
416 *out
417 << "iterativeSolver = "<<describe(*iterativeSolver_,verbLevel)
418 << "fwdOp = " << describe(*lp_->getOperator(),verbLevel);
419 if (lp_->getLeftPrec().get())
420 *out << "leftPrecOp = "<<describe(*lp_->getLeftPrec(),verbLevel);
421 if (lp_->getRightPrec().get())
422 *out << "rightPrecOp = "<<describe(*lp_->getRightPrec(),verbLevel);
423 }
424 break;
425 }
426 default:
427 TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
428 }
429}
430
431
432// protected
433
434
435// Overridden from LinearOpBase
436
437
438template<class Scalar>
440{
441 return ::Thyra::opSupported(*lp_->getOperator(),M_trans);
442}
443
444
445template<class Scalar>
447 const EOpTransp M_trans,
448 const MultiVectorBase<Scalar> &X,
449 const Ptr<MultiVectorBase<Scalar> > &Y,
450 const Scalar alpha,
451 const Scalar beta
452 ) const
453{
454 ::Thyra::apply<Scalar>(*lp_->getOperator(), M_trans, X, Y, alpha, beta);
455}
456
457
458// Overridden from LinearOpWithSolveBase
459
460
461template<class Scalar>
462bool
464{
465 return solveSupportsNewImpl(M_trans, Teuchos::null);
466}
467
468
469template<class Scalar>
470bool
472 const Ptr<const SolveCriteria<Scalar> > /* solveCriteria */) const
473{
474 // Only support forward solve right now!
475 if (real_trans(transp)==NOTRANS) return true;
476 return false; // ToDo: Support adjoint solves!
477 // Otherwise, Thyra/Belos now supports every solve criteria type that exists
478 // because of the class Thyra::GeneralSolveCriteriaBelosStatusTest!
479/*
480 if (real_trans(M_trans)==NOTRANS) {
481 return (
482 solveMeasureType.useDefault()
483 ||
484 solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_RHS)
485 ||
486 solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_INIT_RESIDUAL)
487 );
488 }
489*/
490}
491
492
493template<class Scalar>
494bool
496 EOpTransp M_trans, const SolveMeasureType& solveMeasureType) const
497{
498 SolveCriteria<Scalar> solveCriteria(solveMeasureType, SolveCriteria<Scalar>::unspecifiedTolerance());
499 return solveSupportsNewImpl(M_trans, Teuchos::constOptInArg(solveCriteria));
500}
501
502
503template<class Scalar>
504SolveStatus<Scalar>
506 const EOpTransp M_trans,
507 const MultiVectorBase<Scalar> &B,
508 const Ptr<MultiVectorBase<Scalar> > &X,
509 const Ptr<const SolveCriteria<Scalar> > solveCriteria
510 ) const
511{
512
513 THYRA_FUNC_TIME_MONITOR("Stratimikos: BelosLOWS");
514
515 using Teuchos::rcp;
516 using Teuchos::rcpFromRef;
517 using Teuchos::rcpFromPtr;
518 using Teuchos::FancyOStream;
519 using Teuchos::OSTab;
520 using Teuchos::ParameterList;
521 using Teuchos::parameterList;
522 using Teuchos::describe;
523 typedef Teuchos::ScalarTraits<Scalar> ST;
524 typedef typename ST::magnitudeType ScalarMag;
525 Teuchos::Time totalTimer("Stratimikos::BelosLinearOpWithSolve::totalTime");
526 totalTimer.start(true);
527
528 assertSolveSupports(*this, M_trans, solveCriteria);
529 // 2010/08/22: rabartl: Bug 4915 ToDo: Move the above into the NIV function
530 // solve(...).
531
532 const RCP<FancyOStream> out = this->getOStream();
533 const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
534 OSTab tab = this->getOSTab();
535 if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW)) {
536 *out << "\nStarting iterations with Belos:\n";
537 OSTab tab2(out);
538 *out << "Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel);
539 *out << "Using iterative solver = " << describe(*iterativeSolver_,verbLevel);
540 *out << "With #Eqns="<<B.range()->dim()<<", #RHSs="<<B.domain()->dim()<<" ...\n";
541 }
542
543#ifdef HAVE_STRATIMIKOS_THYRATPETRAADAPTERS
544 //
545 // Write RHS and LHS to file if desired
546 //
547 if (filenameLHS_ != "") {
548 try {
549 auto tmv = Thyra::TpetraOperatorVectorExtraction<Scalar>::getTpetraMultiVector(Teuchos::rcpFromPtr(X));
550 Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameLHS_+ "." + label_ + "." + std::to_string(counter_), *tmv, "", "");
551 } catch (const std::logic_error&) {
552 *out << "Warning: Cannot write LHS multivector to file.\n";
553 }
554 }
555 if (filenameRHS_ != "") {
556 try {
557 auto tmv = Thyra::TpetraOperatorVectorExtraction<Scalar>::getConstTpetraMultiVector(Teuchos::rcpFromRef(B));
558 Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameRHS_+ "." + label_ + "." + std::to_string(counter_), *tmv, "", "");
559 } catch (const std::logic_error&) {
560 *out << "Warning: Cannot write RHS multivector to file.\n";
561 }
562 }
563 ++counter_;
564#endif
565
566 //
567 // Set RHS and LHS
568 //
569
570 bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(B) );
571 TEUCHOS_TEST_FOR_EXCEPTION(
572 ret == false, CatastrophicSolveFailure
573 ,"Error, the Belos::LinearProblem could not be set for the current solve!"
574 );
575
576 //
577 // Set the solution criteria
578 //
579
580 // Parameter list for the current solve.
581 const RCP<ParameterList> tmpPL = Teuchos::parameterList();
582
583 // The solver's valid parameter list.
584 RCP<const ParameterList> validPL = iterativeSolver_->getValidParameters();
585
586 SolveMeasureType solveMeasureType;
587 RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> > generalSolveCriteriaBelosStatusTest;
588 if (nonnull(solveCriteria)) {
589 solveMeasureType = solveCriteria->solveMeasureType;
590 const ScalarMag requestedTol = solveCriteria->requestedTol;
591 if (solveMeasureType.useDefault()) {
592 tmpPL->set("Convergence Tolerance", defaultTol_);
593 }
594 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) {
595 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
596 tmpPL->set("Convergence Tolerance", requestedTol);
597 }
598 else {
599 tmpPL->set("Convergence Tolerance", defaultTol_);
600 }
601 setResidualScalingType (tmpPL, validPL, "Norm of RHS");
602 }
603 else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_INIT_RESIDUAL)) {
604 if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
605 tmpPL->set("Convergence Tolerance", requestedTol);
606 }
607 else {
608 tmpPL->set("Convergence Tolerance", defaultTol_);
609 }
610 setResidualScalingType (tmpPL, validPL, "Norm of Initial Residual");
611 }
612 else {
613 // Set the most generic (and inefficient) solve criteria
614 generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest(
615 *solveCriteria, convergenceTestFrequency_);
616 // Set the verbosity level (one level down)
617 generalSolveCriteriaBelosStatusTest->setOStream(out);
618 generalSolveCriteriaBelosStatusTest->setVerbLevel(incrVerbLevel(verbLevel, -1));
619 // Set the default convergence tolerance to always converged to allow
620 // the above status test to control things.
621 tmpPL->set("Convergence Tolerance", 1.0);
622 }
623 // maximum iterations
624 if (nonnull(solveCriteria->extraParameters)) {
625 if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,"Maximum Iterations")) {
626 tmpPL->set("Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,"Maximum Iterations"));
627 }
628 }
629 // If a preconditioner is on the left, then the implicit residual test
630 // scaling should be the preconditioned initial residual.
631 if (Teuchos::nonnull(lp_->getLeftPrec()) &&
632 validPL->isParameter ("Implicit Residual Scaling"))
633 tmpPL->set("Implicit Residual Scaling",
634 "Norm of Preconditioned Initial Residual");
635 }
636 else {
637 // No solveCriteria was even passed in!
638 tmpPL->set("Convergence Tolerance", defaultTol_);
639 }
640
641 //
642 // Solve the linear system
643 //
644
645 Belos::ReturnType belosSolveStatus;
646 {
647 // Write detailed convergence information if requested for levels >= VERB_LOW
648 RCP<std::ostream>
649 outUsed =
650 ( static_cast<int>(verbLevel) >= static_cast<int>(Teuchos::VERB_LOW)
651 ? out
652 : rcp(new FancyOStream(rcp(new Teuchos::oblackholestream())))
653 );
654 Teuchos::OSTab tab1(outUsed,1,"BELOS");
655 tmpPL->set("Output Stream", outUsed);
656 if (!init_) {
657 iterativeSolver_->setParameters(tmpPL);
658 init_ = !nonnull(solveCriteria);
659 }
660 if (nonnull(generalSolveCriteriaBelosStatusTest)) {
661 iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest);
662 }
663 try {
664 belosSolveStatus = iterativeSolver_->solve();
665 }
666 catch (Belos::BelosError& ex) {
667 TEUCHOS_TEST_FOR_EXCEPTION( true,
668 CatastrophicSolveFailure,
669 ex.what() );
670 }
671 }
672
673 //
674 // Report the solve status
675 //
676
677 totalTimer.stop();
678
679 SolveStatus<Scalar> solveStatus;
680
681 switch (belosSolveStatus) {
682 case Belos::Unconverged: {
683 solveStatus.solveStatus = SOLVE_STATUS_UNCONVERGED;
684 // Set achievedTol even if the solver did not converge. This is
685 // helpful for things like nonlinear solvers, which might be
686 // able to use a partially converged result, and which would
687 // like to know the achieved convergence tolerance for use in
688 // computing bounds. It's also helpful for estimating whether a
689 // small increase in the maximum iteration count might be
690 // helpful next time.
691 try {
692 // Some solvers might not have implemented achievedTol().
693 // The default implementation throws std::runtime_error.
694 solveStatus.achievedTol = iterativeSolver_->achievedTol();
695 } catch (std::runtime_error&) {
696 // Do nothing; use the default value of achievedTol.
697 }
698 break;
699 }
700 case Belos::Converged: {
701 solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
702 if (nonnull(generalSolveCriteriaBelosStatusTest)) {
703 // The user set a custom status test. This means that we
704 // should ask the custom status test itself, rather than the
705 // Belos solver, what the final achieved convergence tolerance
706 // was.
707 const ArrayView<const ScalarMag> achievedTol =
708 generalSolveCriteriaBelosStatusTest->achievedTol();
709 solveStatus.achievedTol = Teuchos::ScalarTraits<ScalarMag>::zero();
710 for (Ordinal i = 0; i < achievedTol.size(); ++i) {
711 solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]);
712 }
713 }
714 else {
715 try {
716 // Some solvers might not have implemented achievedTol().
717 // The default implementation throws std::runtime_error.
718 solveStatus.achievedTol = iterativeSolver_->achievedTol();
719 } catch (std::runtime_error&) {
720 // Use the default convergence tolerance. This is a correct
721 // upper bound, since we did actually converge.
722 solveStatus.achievedTol = tmpPL->get("Convergence Tolerance", defaultTol_);
723 }
724 }
725 break;
726 }
727 TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT();
728 }
729
730 std::ostringstream ossmessage;
731 ossmessage
732 << "The Belos solver " << (label_ != "" ? ("\"" + label_ + "\" ") : "")
733 << "of type \""<<iterativeSolver_->description()
734 <<"\" returned a solve status of \""<< toString(solveStatus.solveStatus) << "\""
735 << " in " << iterativeSolver_->getNumIters() << " iterations"
736 << " with total CPU time of " << totalTimer.totalElapsedTime() << " sec" ;
737 if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
738 *out << "\n" << ossmessage.str() << "\n";
739
740 solveStatus.message = ossmessage.str();
741
742 // Dump the getNumIters() and the achieved convergence tolerance
743 // into solveStatus.extraParameters, as the "Belos/Iteration Count"
744 // resp. "Belos/Achieved Tolerance" parameters.
745 if (solveStatus.extraParameters.is_null()) {
746 solveStatus.extraParameters = parameterList ();
747 }
748 solveStatus.extraParameters->set ("Belos/Iteration Count",
749 iterativeSolver_->getNumIters());\
750 // package independent version of the same
751 solveStatus.extraParameters->set ("Iteration Count",
752 iterativeSolver_->getNumIters());\
753 // NOTE (mfh 13 Dec 2011) Though the most commonly used Belos
754 // solvers do implement achievedTol(), some Belos solvers currently
755 // do not. In the latter case, if the solver did not converge, the
756 // reported achievedTol() value may just be the default "invalid"
757 // value -1, and if the solver did converge, the reported value will
758 // just be the convergence tolerance (a correct upper bound).
759 solveStatus.extraParameters->set ("Belos/Achieved Tolerance",
760 solveStatus.achievedTol);
761
762// This information is in the previous line, which is printed anytime the verbosity
763// is not set to Teuchos::VERB_NONE, so I'm commenting this out for now.
764// if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_NONE))
765// *out << "\nTotal solve time in Belos = "<<totalTimer.totalElapsedTime()<<" sec\n";
766
767 return solveStatus;
768
769}
770
771
772} // end namespace Thyra
773
774
775#endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
void uninitialize(RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > *lp=NULL, RCP< Teuchos::ParameterList > *solverPL=NULL, RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > *iterativeSolver=NULL, RCP< const LinearOpSourceBase< Scalar > > *fwdOpSrc=NULL, RCP< const PreconditionerBase< Scalar > > *prec=NULL, bool *isExternalPrec=NULL, RCP< const LinearOpSourceBase< Scalar > > *approxFwdOpSrc=NULL, ESupportSolveUse *supportSolveUse=NULL)
Uninitializes and returns stored quantities.
RCP< const LinearOpBase< Scalar > > clone() const
virtual bool opSupportedImpl(EOpTransp M_trans) const
RCP< const VectorSpaceBase< Scalar > > domain() const
virtual bool solveSupportsNewImpl(EOpTransp transp, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
RCP< const VectorSpaceBase< Scalar > > range() const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
virtual SolveStatus< Scalar > solveImpl(const EOpTransp transp, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
RCP< const PreconditionerBase< Scalar > > extract_prec()
virtual bool solveSupportsImpl(EOpTransp M_trans) const
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
void initialize(const RCP< Belos::LinearProblem< Scalar, MV_t, LO_t > > &lp, const RCP< Teuchos::ParameterList > &solverPL, const RCP< Belos::SolverManager< Scalar, MV_t, LO_t > > &iterativeSolver, const RCP< const LinearOpSourceBase< Scalar > > &fwdOpSrc, const RCP< const PreconditionerBase< Scalar > > &prec, const bool isExternalPrec, const RCP< const LinearOpSourceBase< Scalar > > &approxFwdOpSrc, const ESupportSolveUse &supportSolveUse, const int convergenceTestFrequency)
Initializes given precreated solver objects.
NOTRANS

Generated for Stratimikos by doxygen 1.9.8