Stratimikos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
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 
26 namespace {
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 
113 namespace Thyra {
114 
115 
116 // Constructors/initializers/accessors
117 
118 
119 template<class Scalar>
121  :convergenceTestFrequency_(-1),
122  isExternalPrec_(false),
123  supportSolveUse_(SUPPORT_SOLVE_UNSPECIFIED),
124  defaultTol_ (-1.0),
125  label_(""),
126  filenameLHS_(""),
127  filenameRHS_(""),
128  counter_(0)
129 {}
130 
131 
132 template<class Scalar>
135  const RCP<Teuchos::ParameterList> &solverPL,
136  const RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > &iterativeSolver,
137  const RCP<const LinearOpSourceBase<Scalar> > &fwdOpSrc,
138  const RCP<const PreconditionerBase<Scalar> > &prec,
139  const bool isExternalPrec_in,
140  const RCP<const LinearOpSourceBase<Scalar> > &approxFwdOpSrc,
141  const ESupportSolveUse &supportSolveUse_in,
142  const int convergenceTestFrequency
143  )
144 {
145  using Teuchos::as;
148  typedef typename Teuchos::ScalarTraits<Scalar>::magnitudeType magnitude_type;
149 
150  this->setLinePrefix("BELOS/T");
151  // ToDo: Validate input
152  lp_ = lp;
153  solverPL_ = solverPL;
154  iterativeSolver_ = iterativeSolver;
155  fwdOpSrc_ = fwdOpSrc;
156  prec_ = prec;
157  isExternalPrec_ = isExternalPrec_in;
158  approxFwdOpSrc_ = approxFwdOpSrc;
159  supportSolveUse_ = supportSolveUse_in;
160  convergenceTestFrequency_ = convergenceTestFrequency;
161  // Check if "Convergence Tolerance" is in the solver parameter list. If
162  // not, use the default from the solver.
163  if ( !is_null(solverPL_) ) {
164  if (solverPL_->isParameter("Convergence Tolerance")) {
165 
166  // Stratimikos prefers tolerances as double, no matter the
167  // Scalar type. However, we also want it to accept the
168  // tolerance as magnitude_type, for example float if the Scalar
169  // type is float or std::complex<float>.
170  if (solverPL_->isType<double> ("Convergence Tolerance")) {
171  defaultTol_ =
172  as<magnitude_type> (solverPL_->get<double> ("Convergence Tolerance"));
173  }
174  else if (std::is_same_v<double, magnitude_type>) {
175  // magnitude_type == double in this case, and we've already
176  // checked double above.
178  true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: "
179  "The \"Convergence Tolerance\" parameter, which you provided, must "
180  "have type double (the type of the magnitude of Scalar = double).");
181  }
182  else if (solverPL_->isType<magnitude_type> ("Convergence Tolerance")) {
183  defaultTol_ = solverPL_->get<magnitude_type> ("Convergence Tolerance");
184  }
185  else {
186  // Throwing InvalidParameterType ensures that the exception's
187  // type is consistent both with what this method would have
188  // thrown before for an unrecognized type, and with what the
189  // user expects in general when the parameter doesn't have the
190  // right type.
192  true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: "
193  "The \"Convergence Tolerance\" parameter, which you provided, must "
194  "have type double (preferred) or the type of the magnitude of Scalar "
195  "= " << TypeNameTraits<Scalar>::name () << ", which is " <<
196  TypeNameTraits<magnitude_type>::name () << " in this case. You can "
197  "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
198  }
199  }
200 
201  if (solverPL_->isParameter("Timer Label") && solverPL_->isType<std::string>("Timer Label")) {
202  label_ = solverPL_->get<std::string>("Timer Label");
203  lp_->setLabel(label_);
204  }
205  if (solverPL_->isParameter("Filename LHS") && solverPL_->isType<std::string>("Filename LHS")) {
206  filenameLHS_ = solverPL_->get<std::string>("Filename LHS");
207  }
208  if (solverPL_->isParameter("Filename RHS") && solverPL_->isType<std::string>("Filename RHS")) {
209  filenameRHS_ = solverPL_->get<std::string>("Filename RHS");
210  }
211  }
212  else {
214  iterativeSolver->getValidParameters();
215 
216  // Stratimikos prefers tolerances as double, no matter the
217  // Scalar type. However, we also want it to accept the
218  // tolerance as magnitude_type, for example float if the Scalar
219  // type is float or std::complex<float>.
220  if (defaultPL->isType<double> ("Convergence Tolerance")) {
221  defaultTol_ =
222  as<magnitude_type> (defaultPL->get<double> ("Convergence Tolerance"));
223  }
224  else if (std::is_same_v<double, magnitude_type>) {
225  // magnitude_type == double in this case, and we've already
226  // checked double above.
228  true, std::invalid_argument, "BelosLinearOpWithSolve::initialize: "
229  "The \"Convergence Tolerance\" parameter, which you provided, must "
230  "have type double (the type of the magnitude of Scalar = double).");
231  }
232  else if (defaultPL->isType<magnitude_type> ("Convergence Tolerance")) {
233  defaultTol_ = defaultPL->get<magnitude_type> ("Convergence Tolerance");
234  }
235  else {
236  // Throwing InvalidParameterType ensures that the exception's
237  // type is consistent both with what this method would have
238  // thrown before for an unrecognized type, and with what the
239  // user expects in general when the parameter doesn't have the
240  // right type.
242  true, InvalidParameterType, "BelosLinearOpWithSolve::initialize: "
243  "The \"Convergence Tolerance\" parameter, which you provided, must "
244  "have type double (preferred) or the type of the magnitude of Scalar "
245  "= " << TypeNameTraits<Scalar>::name () << ", which is " <<
246  TypeNameTraits<magnitude_type>::name () << " in this case. You can "
247  "find that type using Teuchos::ScalarTraits<Scalar>::magnitudeType.");
248  }
249  }
250 }
251 
252 
253 template<class Scalar>
256 {
258  _fwdOpSrc = fwdOpSrc_;
259  fwdOpSrc_ = Teuchos::null;
260  return _fwdOpSrc;
261 }
262 
263 
264 template<class Scalar>
267 {
269  _prec = prec_;
270  prec_ = Teuchos::null;
271  return _prec;
272 }
273 
274 
275 template<class Scalar>
277 {
278  return isExternalPrec_;
279 }
280 
281 
282 template<class Scalar>
285 {
287  _approxFwdOpSrc = approxFwdOpSrc_;
288  approxFwdOpSrc_ = Teuchos::null;
289  return _approxFwdOpSrc;
290 }
291 
292 
293 template<class Scalar>
295 {
296  return supportSolveUse_;
297 }
298 
299 
300 template<class Scalar>
303  RCP<Teuchos::ParameterList> *solverPL,
304  RCP<Belos::SolverManager<Scalar,MV_t,LO_t> > *iterativeSolver,
305  RCP<const LinearOpSourceBase<Scalar> > *fwdOpSrc,
306  RCP<const PreconditionerBase<Scalar> > *prec,
307  bool *isExternalPrec_in,
308  RCP<const LinearOpSourceBase<Scalar> > *approxFwdOpSrc,
309  ESupportSolveUse *supportSolveUse_in
310  )
311 {
312  if (lp) *lp = lp_;
313  if (solverPL) *solverPL = solverPL_;
314  if (iterativeSolver) *iterativeSolver = iterativeSolver_;
315  if (fwdOpSrc) *fwdOpSrc = fwdOpSrc_;
316  if (prec) *prec = prec_;
317  if (isExternalPrec_in) *isExternalPrec_in = isExternalPrec_;
318  if (approxFwdOpSrc) *approxFwdOpSrc = approxFwdOpSrc_;
319  if (supportSolveUse_in) *supportSolveUse_in = supportSolveUse_;
320 
321  lp_ = Teuchos::null;
322  solverPL_ = Teuchos::null;
323  iterativeSolver_ = Teuchos::null;
324  fwdOpSrc_ = Teuchos::null;
325  prec_ = Teuchos::null;
326  isExternalPrec_ = false;
327  approxFwdOpSrc_ = Teuchos::null;
328  supportSolveUse_ = SUPPORT_SOLVE_UNSPECIFIED;
329 }
330 
331 
332 // Overridden from LinearOpBase
333 
334 
335 template<class Scalar>
338 {
339  if (!is_null(lp_))
340  return lp_->getOperator()->range();
341  return Teuchos::null;
342 }
343 
344 
345 template<class Scalar>
348 {
349  if (!is_null(lp_))
350  return lp_->getOperator()->domain();
351  return Teuchos::null;
352 }
353 
354 
355 template<class Scalar>
358 {
359  return Teuchos::null; // Not supported yet but could be
360 }
361 
362 
363 // Overridden from Teuchos::Describable
364 
365 
366 template<class Scalar>
368 {
369  std::ostringstream oss;
371  if ( !is_null(lp_) && !is_null(lp_->getOperator()) ) {
372  oss << "{";
373  oss << "iterativeSolver=\'"<<iterativeSolver_->description()<<"\'";
374  oss << ",fwdOp=\'"<<lp_->getOperator()->description()<<"\'";
375  if (lp_->getLeftPrec().get())
376  oss << ",leftPrecOp=\'"<<lp_->getLeftPrec()->description()<<"\'";
377  if (lp_->getRightPrec().get())
378  oss << ",rightPrecOp=\'"<<lp_->getRightPrec()->description()<<"\'";
379  oss << "}";
380  }
381  // ToDo: Make Belos::SolverManager derive from Teuchos::Describable so
382  // that we can get better information.
383  return oss.str();
384 }
385 
386 
387 template<class Scalar>
389  Teuchos::FancyOStream &out_arg,
390  const Teuchos::EVerbosityLevel verbLevel
391  ) const
392 {
393  using Teuchos::FancyOStream;
394  using Teuchos::OSTab;
395  using Teuchos::describe;
396  RCP<FancyOStream> out = rcp(&out_arg,false);
397  OSTab tab(out);
398  switch (verbLevel) {
399  case Teuchos::VERB_DEFAULT: // fall-through
400  case Teuchos::VERB_LOW:
401  case Teuchos::VERB_MEDIUM: // fall-through
402  *out << this->description() << std::endl;
403  break;
404  case Teuchos::VERB_HIGH: // fall-through
406  {
407  *out
409  << "rangeDim=" << this->range()->dim()
410  << ",domainDim=" << this->domain()->dim() << "}\n";
411  if (lp_->getOperator().get()) {
412  OSTab tab1(out);
413  *out
414  << "iterativeSolver = "<<describe(*iterativeSolver_,verbLevel)
415  << "fwdOp = " << describe(*lp_->getOperator(),verbLevel);
416  if (lp_->getLeftPrec().get())
417  *out << "leftPrecOp = "<<describe(*lp_->getLeftPrec(),verbLevel);
418  if (lp_->getRightPrec().get())
419  *out << "rightPrecOp = "<<describe(*lp_->getRightPrec(),verbLevel);
420  }
421  break;
422  }
423  default:
424  TEUCHOS_TEST_FOR_EXCEPT(true); // Should never get here!
425  }
426 }
427 
428 
429 // protected
430 
431 
432 // Overridden from LinearOpBase
433 
434 
435 template<class Scalar>
437 {
438  return ::Thyra::opSupported(*lp_->getOperator(),M_trans);
439 }
440 
441 
442 template<class Scalar>
444  const EOpTransp M_trans,
445  const MultiVectorBase<Scalar> &X,
446  const Ptr<MultiVectorBase<Scalar> > &Y,
447  const Scalar alpha,
448  const Scalar beta
449  ) const
450 {
451  ::Thyra::apply<Scalar>(*lp_->getOperator(), M_trans, X, Y, alpha, beta);
452 }
453 
454 
455 // Overridden from LinearOpWithSolveBase
456 
457 
458 template<class Scalar>
459 bool
461 {
462  return solveSupportsNewImpl(M_trans, Teuchos::null);
463 }
464 
465 
466 template<class Scalar>
467 bool
469  const Ptr<const SolveCriteria<Scalar> > /* solveCriteria */) const
470 {
471  // Only support forward solve right now!
472  if (real_trans(transp)==NOTRANS) return true;
473  return false; // ToDo: Support adjoint solves!
474  // Otherwise, Thyra/Belos now supports every solve criteria type that exists
475  // because of the class Thyra::GeneralSolveCriteriaBelosStatusTest!
476 /*
477  if (real_trans(M_trans)==NOTRANS) {
478  return (
479  solveMeasureType.useDefault()
480  ||
481  solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_RHS)
482  ||
483  solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL,SOLVE_MEASURE_NORM_INIT_RESIDUAL)
484  );
485  }
486 */
487 }
488 
489 
490 template<class Scalar>
491 bool
493  EOpTransp M_trans, const SolveMeasureType& solveMeasureType) const
494 {
495  SolveCriteria<Scalar> solveCriteria(solveMeasureType, SolveCriteria<Scalar>::unspecifiedTolerance());
496  return solveSupportsNewImpl(M_trans, Teuchos::constOptInArg(solveCriteria));
497 }
498 
499 
500 template<class Scalar>
503  const EOpTransp M_trans,
504  const MultiVectorBase<Scalar> &B,
505  const Ptr<MultiVectorBase<Scalar> > &X,
506  const Ptr<const SolveCriteria<Scalar> > solveCriteria
507  ) const
508 {
509 
510  THYRA_FUNC_TIME_MONITOR("Stratimikos: BelosLOWS");
511 
512  using Teuchos::rcp;
513  using Teuchos::rcpFromRef;
514  using Teuchos::rcpFromPtr;
515  using Teuchos::FancyOStream;
516  using Teuchos::OSTab;
518  using Teuchos::parameterList;
519  using Teuchos::describe;
521  typedef typename ST::magnitudeType ScalarMag;
522  Teuchos::Time totalTimer(""), timer("");
523  totalTimer.start(true);
524 
525  assertSolveSupports(*this, M_trans, solveCriteria);
526  // 2010/08/22: rabartl: Bug 4915 ToDo: Move the above into the NIV function
527  // solve(...).
528 
529  const RCP<FancyOStream> out = this->getOStream();
530  const Teuchos::EVerbosityLevel verbLevel = this->getVerbLevel();
531  OSTab tab = this->getOSTab();
532  if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW)) {
533  *out << "\nStarting iterations with Belos:\n";
534  OSTab tab2(out);
535  *out << "Using forward operator = " << describe(*fwdOpSrc_->getOp(),verbLevel);
536  *out << "Using iterative solver = " << describe(*iterativeSolver_,verbLevel);
537  *out << "With #Eqns="<<B.range()->dim()<<", #RHSs="<<B.domain()->dim()<<" ...\n";
538  }
539 
540 #ifdef HAVE_STRATIMIKOS_THYRATPETRAADAPTERS
541  //
542  // Write RHS and LHS to file if desired
543  //
544  if (filenameLHS_ != "") {
545  try {
547  Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameLHS_+ "." + label_ + "." + std::to_string(counter_), *tmv, "", "");
548  } catch (const std::logic_error&) {
549  *out << "Warning: Cannot write LHS multivector to file.\n";
550  }
551  }
552  if (filenameRHS_ != "") {
553  try {
555  Tpetra::MatrixMarket::Writer<::Tpetra::CrsMatrix<Scalar> >::writeDenseFile(filenameRHS_+ "." + label_ + "." + std::to_string(counter_), *tmv, "", "");
556  } catch (const std::logic_error&) {
557  *out << "Warning: Cannot write RHS multivector to file.\n";
558  }
559  }
560  ++counter_;
561 #endif
562 
563  //
564  // Set RHS and LHS
565  //
566 
567  bool ret = lp_->setProblem( rcpFromPtr(X), rcpFromRef(B) );
569  ret == false, CatastrophicSolveFailure
570  ,"Error, the Belos::LinearProblem could not be set for the current solve!"
571  );
572 
573  //
574  // Set the solution criteria
575  //
576 
577  // Parameter list for the current solve.
578  const RCP<ParameterList> tmpPL = Teuchos::parameterList();
579 
580  // The solver's valid parameter list.
581  RCP<const ParameterList> validPL = iterativeSolver_->getValidParameters();
582 
583  SolveMeasureType solveMeasureType;
584  RCP<GeneralSolveCriteriaBelosStatusTest<Scalar> > generalSolveCriteriaBelosStatusTest;
585  if (nonnull(solveCriteria)) {
586  solveMeasureType = solveCriteria->solveMeasureType;
587  const ScalarMag requestedTol = solveCriteria->requestedTol;
588  if (solveMeasureType.useDefault()) {
589  tmpPL->set("Convergence Tolerance", defaultTol_);
590  }
591  else if (solveMeasureType(SOLVE_MEASURE_NORM_RESIDUAL, SOLVE_MEASURE_NORM_RHS)) {
592  if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
593  tmpPL->set("Convergence Tolerance", requestedTol);
594  }
595  else {
596  tmpPL->set("Convergence Tolerance", defaultTol_);
597  }
598  setResidualScalingType (tmpPL, validPL, "Norm of RHS");
599  }
601  if (requestedTol != SolveCriteria<Scalar>::unspecifiedTolerance()) {
602  tmpPL->set("Convergence Tolerance", requestedTol);
603  }
604  else {
605  tmpPL->set("Convergence Tolerance", defaultTol_);
606  }
607  setResidualScalingType (tmpPL, validPL, "Norm of Initial Residual");
608  }
609  else {
610  // Set the most generic (and inefficient) solve criteria
611  generalSolveCriteriaBelosStatusTest = createGeneralSolveCriteriaBelosStatusTest(
612  *solveCriteria, convergenceTestFrequency_);
613  // Set the verbosity level (one level down)
614  generalSolveCriteriaBelosStatusTest->setOStream(out);
615  generalSolveCriteriaBelosStatusTest->setVerbLevel(incrVerbLevel(verbLevel, -1));
616  // Set the default convergence tolerance to always converged to allow
617  // the above status test to control things.
618  tmpPL->set("Convergence Tolerance", 1.0);
619  }
620  // maximum iterations
621  if (nonnull(solveCriteria->extraParameters)) {
622  if (Teuchos::isParameterType<int>(*solveCriteria->extraParameters,"Maximum Iterations")) {
623  tmpPL->set("Maximum Iterations", Teuchos::get<int>(*solveCriteria->extraParameters,"Maximum Iterations"));
624  }
625  }
626  // If a preconditioner is on the left, then the implicit residual test
627  // scaling should be the preconditioned initial residual.
628  if (Teuchos::nonnull(lp_->getLeftPrec()) &&
629  validPL->isParameter ("Implicit Residual Scaling"))
630  tmpPL->set("Implicit Residual Scaling",
631  "Norm of Preconditioned Initial Residual");
632  }
633  else {
634  // No solveCriteria was even passed in!
635  tmpPL->set("Convergence Tolerance", defaultTol_);
636  }
637 
638  //
639  // Solve the linear system
640  //
641 
642  Belos::ReturnType belosSolveStatus;
643  {
644  // Write detailed convergence information if requested for levels >= VERB_LOW
646  outUsed =
647  ( static_cast<int>(verbLevel) >= static_cast<int>(Teuchos::VERB_LOW)
648  ? out
650  );
651  Teuchos::OSTab tab1(outUsed,1,"BELOS");
652  tmpPL->set("Output Stream", outUsed);
653  iterativeSolver_->setParameters(tmpPL);
654  if (nonnull(generalSolveCriteriaBelosStatusTest)) {
655  iterativeSolver_->setUserConvStatusTest(generalSolveCriteriaBelosStatusTest);
656  }
657  try {
658  belosSolveStatus = iterativeSolver_->solve();
659  }
660  catch (Belos::BelosError& ex) {
663  ex.what() );
664  }
665  }
666 
667  //
668  // Report the solve status
669  //
670 
671  totalTimer.stop();
672 
673  SolveStatus<Scalar> solveStatus;
674 
675  switch (belosSolveStatus) {
676  case Belos::Unconverged: {
678  // Set achievedTol even if the solver did not converge. This is
679  // helpful for things like nonlinear solvers, which might be
680  // able to use a partially converged result, and which would
681  // like to know the achieved convergence tolerance for use in
682  // computing bounds. It's also helpful for estimating whether a
683  // small increase in the maximum iteration count might be
684  // helpful next time.
685  try {
686  // Some solvers might not have implemented achievedTol().
687  // The default implementation throws std::runtime_error.
688  solveStatus.achievedTol = iterativeSolver_->achievedTol();
689  } catch (std::runtime_error&) {
690  // Do nothing; use the default value of achievedTol.
691  }
692  break;
693  }
694  case Belos::Converged: {
695  solveStatus.solveStatus = SOLVE_STATUS_CONVERGED;
696  if (nonnull(generalSolveCriteriaBelosStatusTest)) {
697  // The user set a custom status test. This means that we
698  // should ask the custom status test itself, rather than the
699  // Belos solver, what the final achieved convergence tolerance
700  // was.
701  const ArrayView<const ScalarMag> achievedTol =
702  generalSolveCriteriaBelosStatusTest->achievedTol();
704  for (Ordinal i = 0; i < achievedTol.size(); ++i) {
705  solveStatus.achievedTol = std::max(solveStatus.achievedTol, achievedTol[i]);
706  }
707  }
708  else {
709  try {
710  // Some solvers might not have implemented achievedTol().
711  // The default implementation throws std::runtime_error.
712  solveStatus.achievedTol = iterativeSolver_->achievedTol();
713  } catch (std::runtime_error&) {
714  // Use the default convergence tolerance. This is a correct
715  // upper bound, since we did actually converge.
716  solveStatus.achievedTol = tmpPL->get("Convergence Tolerance", defaultTol_);
717  }
718  }
719  break;
720  }
722  }
723 
724  std::ostringstream ossmessage;
725  ossmessage
726  << "The Belos solver " << (label_ != "" ? ("\"" + label_ + "\" ") : "")
727  << "of type \""<<iterativeSolver_->description()
728  <<"\" returned a solve status of \""<< toString(solveStatus.solveStatus) << "\""
729  << " in " << iterativeSolver_->getNumIters() << " iterations"
730  << " with total CPU time of " << totalTimer.totalElapsedTime() << " sec" ;
731  if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_LOW))
732  *out << "\n" << ossmessage.str() << "\n";
733 
734  solveStatus.message = ossmessage.str();
735 
736  // Dump the getNumIters() and the achieved convergence tolerance
737  // into solveStatus.extraParameters, as the "Belos/Iteration Count"
738  // resp. "Belos/Achieved Tolerance" parameters.
739  if (solveStatus.extraParameters.is_null()) {
740  solveStatus.extraParameters = parameterList ();
741  }
742  solveStatus.extraParameters->set ("Belos/Iteration Count",
743  iterativeSolver_->getNumIters());\
744  // package independent version of the same
745  solveStatus.extraParameters->set ("Iteration Count",
746  iterativeSolver_->getNumIters());\
747  // NOTE (mfh 13 Dec 2011) Though the most commonly used Belos
748  // solvers do implement achievedTol(), some Belos solvers currently
749  // do not. In the latter case, if the solver did not converge, the
750  // reported achievedTol() value may just be the default "invalid"
751  // value -1, and if the solver did converge, the reported value will
752  // just be the convergence tolerance (a correct upper bound).
753  solveStatus.extraParameters->set ("Belos/Achieved Tolerance",
754  solveStatus.achievedTol);
755 
756 // This information is in the previous line, which is printed anytime the verbosity
757 // is not set to Teuchos::VERB_NONE, so I'm commenting this out for now.
758 // if (out.get() && static_cast<int>(verbLevel) > static_cast<int>(Teuchos::VERB_NONE))
759 // *out << "\nTotal solve time in Belos = "<<totalTimer.totalElapsedTime()<<" sec\n";
760 
761  return solveStatus;
762 
763 }
764 
765 
766 } // end namespace Thyra
767 
768 
769 #endif // THYRA_BELOS_LINEAR_OP_WITH_SOLVE_HPP
virtual bool solveSupportsImpl(EOpTransp M_trans) const
RCP< const VectorSpaceBase< Scalar > > domain() const
bool is_null(const boost::shared_ptr< T > &p)
virtual RCP< const VectorSpaceBase< Scalar > > range() const =0
EOpTransp
bool useDefault() const
#define TEUCHOS_SWITCH_DEFAULT_DEBUG_ASSERT()
basic_OSTab< char > OSTab
std::string message
virtual SolveStatus< Scalar > solveImpl(const EOpTransp transp, const MultiVectorBase< Scalar > &B, const Ptr< MultiVectorBase< Scalar > > &X, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
RCP< const LinearOpSourceBase< Scalar > > extract_approxFwdOpSrc()
T & get(const std::string &name, T def_value)
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.
bool nonnull(const std::shared_ptr< T > &p)
basic_FancyOStream< char > FancyOStream
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
RCP< const PreconditionerBase< Scalar > > extract_prec()
virtual bool solveSupportsNewImpl(EOpTransp transp, const Ptr< const SolveCriteria< Scalar > > solveCriteria) const
EOpTransp real_trans(EOpTransp transp)
T * get() const
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
RCP< const VectorSpaceBase< Scalar > > range() const
virtual void applyImpl(const EOpTransp M_trans, const MultiVectorBase< Scalar > &X, const Ptr< MultiVectorBase< Scalar > > &Y, const Scalar alpha, const Scalar beta) const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
bool isParameter(const std::string &name) const
static RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetraMultiVector(const RCP< MultiVectorBase< Scalar > > &mv)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::Ordinal Ordinal
virtual std::string description() const
TEUCHOSCORE_LIB_DLL_EXPORT std::string toString(const EVerbosityLevel verbLevel)
ParameterList & setParameters(const ParameterList &source)
ESupportSolveUse
ESolveStatus solveStatus
RCP< const LinearOpSourceBase< Scalar > > extract_fwdOpSrc()
TEUCHOSCORE_LIB_DLL_EXPORT EVerbosityLevel incrVerbLevel(const EVerbosityLevel inputVerbLevel, const int numLevels)
virtual RCP< const VectorSpaceBase< Scalar > > domain() const =0
bool nonnull(const boost::shared_ptr< T > &p)
ScalarMag achievedTol
TypeTo as(const TypeFrom &t)
bool isType(const std::string &name) const
static RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getConstTpetraMultiVector(const RCP< const MultiVectorBase< Scalar > > &mv)
BelosLinearOpWithSolve()
Construct to unintialize.
RCP< ParameterList > extraParameters
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
virtual bool solveSupportsSolveMeasureTypeImpl(EOpTransp M_trans, const SolveMeasureType &solveMeasureType) const
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)

Generated on Fri Dec 20 2024 09:21:59 for Stratimikos by doxygen 1.8.5