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

Generated on Thu Feb 27 2020 16:06:46 for Stratimikos by doxygen 1.8.16