Belos  Version of the Day
BelosGmresPolySolMgr.hpp
Go to the documentation of this file.
1 //@HEADER
2 // ************************************************************************
3 //
4 // Belos: Block Linear Solvers Package
5 // Copyright 2004 Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 //@HEADER
41 
42 #ifndef BELOS_GMRES_POLY_SOLMGR_HPP
43 #define BELOS_GMRES_POLY_SOLMGR_HPP
44 
48 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosLinearProblem.hpp"
53 #include "BelosSolverManager.hpp"
54 
55 #include "BelosGmresPolyOp.hpp"
56 #include "BelosGmresIteration.hpp"
57 #include "BelosBlockGmresIter.hpp"
64 #include "BelosStatusTestCombo.hpp"
66 #include "BelosOutputManager.hpp"
67 #include "Teuchos_BLAS.hpp"
68 #include "Teuchos_as.hpp"
69 #ifdef BELOS_TEUCHOS_TIME_MONITOR
70 #include "Teuchos_TimeMonitor.hpp"
71 #endif
72 
73 
74 namespace Belos {
75 
77 
78 
86  GmresPolySolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
87  {}};
88 
96  GmresPolySolMgrPolynomialFailure(const std::string& what_arg) : BelosError(what_arg)
97  {}};
98 
106  GmresPolySolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
107  {}};
108 
154 template<class ScalarType, class MV, class OP>
155 class GmresPolySolMgr : public SolverManager<ScalarType,MV,OP> {
156 private:
159  typedef Teuchos::ScalarTraits<ScalarType> STS;
160  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
161  typedef Teuchos::ScalarTraits<MagnitudeType> MT;
162 
163 public:
164 
166 
167 
173  GmresPolySolMgr();
174 
188  GmresPolySolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
189  const Teuchos::RCP<Teuchos::ParameterList> &pl );
190 
192  virtual ~GmresPolySolMgr() {};
194 
196 
197 
201  return *problem_;
202  }
203 
206  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
207 
210  Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
211 
217  Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
218  return Teuchos::tuple(timerSolve_, timerPoly_);
219  }
220 
222  int getNumIters() const {
223  return numIters_;
224  }
225 
229  bool isLOADetected() const { return loaDetected_; }
230 
232 
234 
235 
237  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) { problem_ = problem; isSTSet_ = false; }
238 
240  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params );
241 
243 
245 
254  void reset( const ResetType type ) {
255  if ((type & Belos::Problem) && ! problem_.is_null ()) {
256  problem_->setProblem ();
257  isPolyBuilt_ = false; // Rebuild the GMRES polynomial
258  }
259  }
260 
262 
264 
282  ReturnType solve();
283 
285 
288 
290  std::string description() const;
291 
293 
294 private:
295 
296  // Method for checking current status test against defined linear problem.
297  bool checkStatusTest();
298 
299  // Method for generating GMRES polynomial.
300  bool generatePoly();
301 
302  // Linear problem.
303  Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
304 
305  // Output manager.
306  Teuchos::RCP<OutputManager<ScalarType> > printer_;
307  Teuchos::RCP<std::ostream> outputStream_;
308 
309  // Status test.
310  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
311  Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
312  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
313  Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
314  Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
315 
316  // Orthogonalization manager.
317  Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
318 
319  // Current parameter list.
320  Teuchos::RCP<Teuchos::ParameterList> params_;
321 
322  // Default solver values.
323  static const MagnitudeType polytol_default_;
324  static const MagnitudeType convtol_default_;
325  static const MagnitudeType orthoKappa_default_;
326  static const int maxDegree_default_;
327  static const int maxRestarts_default_;
328  static const int maxIters_default_;
329  static const bool strictConvTol_default_;
330  static const bool showMaxResNormOnly_default_;
331  static const int blockSize_default_;
332  static const int numBlocks_default_;
333  static const int verbosity_default_;
334  static const int outputStyle_default_;
335  static const int outputFreq_default_;
336  static const std::string impResScale_default_;
337  static const std::string expResScale_default_;
338  static const std::string label_default_;
339  static const std::string orthoType_default_;
340  static const Teuchos::RCP<std::ostream> outputStream_default_;
341 
342  // Current solver values.
343  MagnitudeType polytol_, convtol_, orthoKappa_;
344  int maxDegree_, maxRestarts_, maxIters_, numIters_;
345  int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_;
346  bool strictConvTol_, showMaxResNormOnly_;
347  std::string orthoType_;
348  std::string impResScale_, expResScale_;
349 
350  // Polynomial storage
351  int poly_dim_;
352  Teuchos::RCP<Teuchos::SerialDenseMatrix<int, ScalarType> > poly_H_, poly_y_;
353  Teuchos::RCP<Teuchos::SerialDenseVector<int, ScalarType> > poly_r0_;
354  Teuchos::RCP<Belos::GmresPolyOp<ScalarType, MV, OP> > poly_Op_;
355 
356  // Timers.
357  std::string label_;
358  Teuchos::RCP<Teuchos::Time> timerSolve_, timerPoly_;
359 
360  // Internal state variables.
361  bool isPolyBuilt_;
362  bool isSet_, isSTSet_, expResTest_;
363  bool loaDetected_;
364 
366  mutable Teuchos::RCP<const Teuchos::ParameterList> validPL_;
367 };
368 
369 
370 // Default solver values.
371 template<class ScalarType, class MV, class OP>
372 const typename GmresPolySolMgr<ScalarType,MV,OP>::MagnitudeType
374 
375 template<class ScalarType, class MV, class OP>
376 const typename GmresPolySolMgr<ScalarType,MV,OP>::MagnitudeType
378 
379 template<class ScalarType, class MV, class OP>
380 const typename GmresPolySolMgr<ScalarType,MV,OP>::MagnitudeType
382  -Teuchos::ScalarTraits<MagnitudeType>::one();
383 
384 template<class ScalarType, class MV, class OP>
386 
387 template<class ScalarType, class MV, class OP>
389 
390 template<class ScalarType, class MV, class OP>
392 
393 template<class ScalarType, class MV, class OP>
395 
396 template<class ScalarType, class MV, class OP>
398 
399 template<class ScalarType, class MV, class OP>
401 
402 template<class ScalarType, class MV, class OP>
404 
405 template<class ScalarType, class MV, class OP>
407 
408 template<class ScalarType, class MV, class OP>
410 
411 template<class ScalarType, class MV, class OP>
413 
414 template<class ScalarType, class MV, class OP>
415 const std::string GmresPolySolMgr<ScalarType,MV,OP>::impResScale_default_ = "Norm of RHS";
416 
417 template<class ScalarType, class MV, class OP>
418 const std::string GmresPolySolMgr<ScalarType,MV,OP>::expResScale_default_ = "Norm of RHS";
419 
420 template<class ScalarType, class MV, class OP>
421 const std::string GmresPolySolMgr<ScalarType,MV,OP>::label_default_ = "Belos";
422 
423 template<class ScalarType, class MV, class OP>
425 
426 template<class ScalarType, class MV, class OP>
427 const Teuchos::RCP<std::ostream>
428 GmresPolySolMgr<ScalarType,MV,OP>::outputStream_default_ = Teuchos::rcpFromRef (std::cout);
429 
430 
431 template<class ScalarType, class MV, class OP>
433  outputStream_ (outputStream_default_),
434  polytol_ (polytol_default_),
435  convtol_ (convtol_default_),
436  orthoKappa_ (orthoKappa_default_),
437  maxDegree_ (maxDegree_default_),
438  maxRestarts_ (maxRestarts_default_),
439  maxIters_ (maxIters_default_),
440  numIters_ (0),
441  blockSize_ (blockSize_default_),
442  numBlocks_ (numBlocks_default_),
443  verbosity_ (verbosity_default_),
444  outputStyle_ (outputStyle_default_),
445  outputFreq_ (outputFreq_default_),
446  strictConvTol_ (strictConvTol_default_),
447  showMaxResNormOnly_ (showMaxResNormOnly_default_),
448  orthoType_ (orthoType_default_),
449  impResScale_ (impResScale_default_),
450  expResScale_ (expResScale_default_),
451  poly_dim_ (0),
452  label_ (label_default_),
453  isPolyBuilt_ (false),
454  isSet_ (false),
455  isSTSet_ (false),
456  expResTest_ (false),
457  loaDetected_ (false)
458 {}
459 
460 
461 template<class ScalarType, class MV, class OP>
463 GmresPolySolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
464  const Teuchos::RCP<Teuchos::ParameterList> &pl) :
465  problem_ (problem),
466  outputStream_ (outputStream_default_),
467  polytol_ (polytol_default_),
468  convtol_ (convtol_default_),
469  orthoKappa_ (orthoKappa_default_),
470  maxDegree_ (maxDegree_default_),
471  maxRestarts_ (maxRestarts_default_),
472  maxIters_ (maxIters_default_),
473  numIters_ (0),
474  blockSize_ (blockSize_default_),
475  numBlocks_ (numBlocks_default_),
476  verbosity_ (verbosity_default_),
477  outputStyle_ (outputStyle_default_),
478  outputFreq_ (outputFreq_default_),
479  strictConvTol_ (strictConvTol_default_),
480  showMaxResNormOnly_ (showMaxResNormOnly_default_),
481  orthoType_ (orthoType_default_),
482  impResScale_ (impResScale_default_),
483  expResScale_ (expResScale_default_),
484  poly_dim_ (0),
485  label_ (label_default_),
486  isPolyBuilt_ (false),
487  isSet_ (false),
488  isSTSet_ (false),
489  expResTest_ (false),
490  loaDetected_ (false)
491 {
492  TEUCHOS_TEST_FOR_EXCEPTION(
493  problem_.is_null (), std::invalid_argument,
494  "Belos::GmresPolySolMgr: The given linear problem is null. "
495  "Please call this constructor with a nonnull LinearProblem argument, "
496  "or call the constructor that does not take a LinearProblem.");
497 
498  // If the input parameter list is null, then the parameters take
499  // default values.
500  if (! pl.is_null ()) {
501  setParameters (pl);
502  }
503 }
504 
505 
506 template<class ScalarType, class MV, class OP>
507 Teuchos::RCP<const Teuchos::ParameterList>
509 {
510  if (validPL_.is_null ()) {
511  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList ();
512  pl->set("Polynomial Tolerance", polytol_default_,
513  "The relative residual tolerance that used to construct the GMRES polynomial.");
514  pl->set("Maximum Degree", maxDegree_default_,
515  "The maximum degree allowed for any GMRES polynomial.");
516  pl->set("Convergence Tolerance", convtol_default_,
517  "The relative residual tolerance that needs to be achieved by the\n"
518  "iterative solver in order for the linear system to be declared converged." );
519  pl->set("Maximum Restarts", maxRestarts_default_,
520  "The maximum number of restarts allowed for each\n"
521  "set of RHS solved.");
522  pl->set("Maximum Iterations", maxIters_default_,
523  "The maximum number of block iterations allowed for each\n"
524  "set of RHS solved.");
525  pl->set("Num Blocks", numBlocks_default_,
526  "The maximum number of blocks allowed in the Krylov subspace\n"
527  "for each set of RHS solved.");
528  pl->set("Block Size", blockSize_default_,
529  "The number of vectors in each block. This number times the\n"
530  "number of blocks is the total Krylov subspace dimension.");
531  pl->set("Verbosity", verbosity_default_,
532  "What type(s) of solver information should be outputted\n"
533  "to the output stream.");
534  pl->set("Output Style", outputStyle_default_,
535  "What style is used for the solver information outputted\n"
536  "to the output stream.");
537  pl->set("Output Frequency", outputFreq_default_,
538  "How often convergence information should be outputted\n"
539  "to the output stream.");
540  pl->set("Output Stream", outputStream_default_,
541  "A reference-counted pointer to the output stream where all\n"
542  "solver output is sent.");
543  pl->set("Strict Convergence", strictConvTol_default_,
544  "After polynomial is applied, whether solver should try to achieve\n"
545  "the relative residual tolerance.");
546  pl->set("Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
547  "When convergence information is printed, only show the maximum\n"
548  "relative residual norm when the block size is greater than one.");
549  pl->set("Implicit Residual Scaling", impResScale_default_,
550  "The type of scaling used in the implicit residual convergence test.");
551  pl->set("Explicit Residual Scaling", expResScale_default_,
552  "The type of scaling used in the explicit residual convergence test.");
553  pl->set("Timer Label", label_default_,
554  "The string to use as a prefix for the timer labels.");
555  pl->set("Orthogonalization", orthoType_default_,
556  "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
557  pl->set("Orthogonalization Constant",orthoKappa_default_,
558  "The constant used by DGKS orthogonalization to determine\n"
559  "whether another step of classical Gram-Schmidt is necessary.");
560  validPL_ = pl;
561  }
562  return validPL_;
563 }
564 
565 
566 template<class ScalarType, class MV, class OP>
568 setParameters (const Teuchos::RCP<Teuchos::ParameterList>& params)
569 {
570  // Create the internal parameter list if ones doesn't already exist.
571  if (params_.is_null ()) {
572  params_ = Teuchos::parameterList (*getValidParameters ());
573  }
574  else {
575  params->validateParameters (*getValidParameters ());
576  }
577 
578  // Check for maximum polynomial degree
579  if (params->isParameter("Maximum Degree")) {
580  maxDegree_ = params->get("Maximum Degree",maxDegree_default_);
581 
582  // Update parameter in our list.
583  params_->set("Maximum Degree", maxDegree_);
584  }
585 
586  // Check for maximum number of restarts
587  if (params->isParameter("Maximum Restarts")) {
588  maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
589 
590  // Update parameter in our list.
591  params_->set("Maximum Restarts", maxRestarts_);
592  }
593 
594  // Check for maximum number of iterations
595  if (params->isParameter("Maximum Iterations")) {
596  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
597 
598  // Update parameter in our list and in status test.
599  params_->set("Maximum Iterations", maxIters_);
600  if (maxIterTest_!=Teuchos::null)
601  maxIterTest_->setMaxIters( maxIters_ );
602  }
603 
604  // Check for blocksize
605  if (params->isParameter("Block Size")) {
606  blockSize_ = params->get("Block Size",blockSize_default_);
607  TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
608  "Belos::GmresPolySolMgr: \"Block Size\" must be strictly positive.");
609 
610  // Update parameter in our list.
611  params_->set("Block Size", blockSize_);
612  }
613 
614  // Check for the maximum number of blocks.
615  if (params->isParameter("Num Blocks")) {
616  numBlocks_ = params->get("Num Blocks",numBlocks_default_);
617  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
618  "Belos::GmresPolySolMgr: \"Num Blocks\" must be strictly positive.");
619 
620  // Update parameter in our list.
621  params_->set("Num Blocks", numBlocks_);
622  }
623 
624  // Check to see if the timer label changed.
625  if (params->isParameter("Timer Label")) {
626  std::string tempLabel = params->get("Timer Label", label_default_);
627 
628  // Update parameter in our list and solver timer
629  if (tempLabel != label_) {
630  label_ = tempLabel;
631  params_->set("Timer Label", label_);
632  std::string solveLabel = label_ + ": GmresPolySolMgr total solve time";
633 #ifdef BELOS_TEUCHOS_TIME_MONITOR
634  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
635 #endif
636  std::string polyLabel = label_ + ": GmresPolySolMgr polynomial creation time";
637 #ifdef BELOS_TEUCHOS_TIME_MONITOR
638  timerPoly_ = Teuchos::TimeMonitor::getNewCounter(polyLabel);
639 #endif
640  if (ortho_ != Teuchos::null) {
641  ortho_->setLabel( label_ );
642  }
643  }
644  }
645 
646  // Check if the orthogonalization changed.
647  if (params->isParameter("Orthogonalization")) {
648  std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
649  TEUCHOS_TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS",
650  std::invalid_argument,
651  "Belos::GmresPolySolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
652  if (tempOrthoType != orthoType_) {
653  params_->set("Orthogonalization", orthoType_);
654  orthoType_ = tempOrthoType;
655  // Create orthogonalization manager
656  if (orthoType_=="DGKS") {
657  if (orthoKappa_ <= 0) {
658  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
659  }
660  else {
661  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
662  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
663  }
664  }
665  else if (orthoType_=="ICGS") {
666  ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
667  }
668  else if (orthoType_=="IMGS") {
669  ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
670  }
671  }
672  }
673 
674  // Check which orthogonalization constant to use.
675  if (params->isParameter("Orthogonalization Constant")) {
676  orthoKappa_ = params->get("Orthogonalization Constant",orthoKappa_default_);
677 
678  // Update parameter in our list.
679  params_->set("Orthogonalization Constant",orthoKappa_);
680  if (orthoType_=="DGKS") {
681  if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
682  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
683  }
684  }
685  }
686 
687  // Check for a change in verbosity level
688  if (params->isParameter("Verbosity")) {
689  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
690  verbosity_ = params->get("Verbosity", verbosity_default_);
691  } else {
692  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
693  }
694 
695  // Update parameter in our list.
696  params_->set("Verbosity", verbosity_);
697  if (printer_ != Teuchos::null)
698  printer_->setVerbosity(verbosity_);
699  }
700 
701  // Check for a change in output style
702  if (params->isParameter("Output Style")) {
703  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
704  outputStyle_ = params->get("Output Style", outputStyle_default_);
705  } else {
706  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
707  }
708 
709  // Reconstruct the convergence test if the explicit residual test is not being used.
710  params_->set("Output Style", outputStyle_);
711  if (outputTest_ != Teuchos::null) {
712  isSTSet_ = false;
713  }
714  }
715 
716  // output stream
717  if (params->isParameter("Output Stream")) {
718  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
719 
720  // Update parameter in our list.
721  params_->set("Output Stream", outputStream_);
722  if (printer_ != Teuchos::null)
723  printer_->setOStream( outputStream_ );
724  }
725 
726  // frequency level
727  if (verbosity_ & Belos::StatusTestDetails) {
728  if (params->isParameter("Output Frequency")) {
729  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
730  }
731 
732  // Update parameter in out list and output status test.
733  params_->set("Output Frequency", outputFreq_);
734  if (outputTest_ != Teuchos::null)
735  outputTest_->setOutputFrequency( outputFreq_ );
736  }
737 
738  // Create output manager if we need to.
739  if (printer_ == Teuchos::null) {
740  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
741  }
742 
743  // Convergence
744  //typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t; // unused
745  //typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t; // unused
746 
747  // Check for polynomial convergence tolerance
748  if (params->isParameter("Polynomial Tolerance")) {
749  polytol_ = params->get("Polynomial Tolerance",polytol_default_);
750 
751  // Update parameter in our list and residual tests.
752  params_->set("Polynomial Tolerance", polytol_);
753  }
754 
755  // Check for convergence tolerance
756  if (params->isParameter("Convergence Tolerance")) {
757  convtol_ = params->get("Convergence Tolerance",convtol_default_);
758 
759  // Update parameter in our list and residual tests.
760  params_->set("Convergence Tolerance", convtol_);
761  if (impConvTest_ != Teuchos::null)
762  impConvTest_->setTolerance( convtol_ );
763  if (expConvTest_ != Teuchos::null)
764  expConvTest_->setTolerance( convtol_ );
765  }
766 
767  // Check if user requires solver to reach convergence tolerance
768  if (params->isParameter("Strict Convergence")) {
769  strictConvTol_ = params->get("Strict Convergence",strictConvTol_default_);
770 
771  // Update parameter in our list and residual tests
772  params_->set("Strict Convergence", strictConvTol_);
773  }
774 
775  // Check for a change in scaling, if so we need to build new residual tests.
776  if (params->isParameter("Implicit Residual Scaling")) {
777  std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
778 
779  // Only update the scaling if it's different.
780  if (impResScale_ != tempImpResScale) {
781  Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
782  impResScale_ = tempImpResScale;
783 
784  // Update parameter in our list and residual tests
785  params_->set("Implicit Residual Scaling", impResScale_);
786  if (impConvTest_ != Teuchos::null) {
787  try {
788  impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
789  }
790  catch (std::exception& e) {
791  // Make sure the convergence test gets constructed again.
792  isSTSet_ = false;
793  }
794  }
795  }
796  }
797 
798  if (params->isParameter("Explicit Residual Scaling")) {
799  std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
800 
801  // Only update the scaling if it's different.
802  if (expResScale_ != tempExpResScale) {
803  Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
804  expResScale_ = tempExpResScale;
805 
806  // Update parameter in our list and residual tests
807  params_->set("Explicit Residual Scaling", expResScale_);
808  if (expConvTest_ != Teuchos::null) {
809  try {
810  expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
811  }
812  catch (std::exception& e) {
813  // Make sure the convergence test gets constructed again.
814  isSTSet_ = false;
815  }
816  }
817  }
818  }
819 
820 
821  if (params->isParameter("Show Maximum Residual Norm Only")) {
822  showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
823 
824  // Update parameter in our list and residual tests
825  params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
826  if (impConvTest_ != Teuchos::null)
827  impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
828  if (expConvTest_ != Teuchos::null)
829  expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
830  }
831 
832  // Create orthogonalization manager if we need to.
833  if (ortho_ == Teuchos::null) {
834  params_->set("Orthogonalization", orthoType_);
835  if (orthoType_=="DGKS") {
836  if (orthoKappa_ <= 0) {
837  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
838  }
839  else {
840  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
841  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
842  }
843  }
844  else if (orthoType_=="ICGS") {
845  ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
846  }
847  else if (orthoType_=="IMGS") {
848  ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
849  }
850  else {
851  TEUCHOS_TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::logic_error,
852  "Belos::GmresPolySolMgr(): Invalid orthogonalization type.");
853  }
854  }
855 
856  // Create the timers if we need to.
857  if (timerSolve_ == Teuchos::null) {
858  std::string solveLabel = label_ + ": GmresPolySolMgr total solve time";
859 #ifdef BELOS_TEUCHOS_TIME_MONITOR
860  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
861 #endif
862  }
863 
864  if (timerPoly_ == Teuchos::null) {
865  std::string polyLabel = label_ + ": GmresPolySolMgr polynomial creation time";
866 #ifdef BELOS_TEUCHOS_TIME_MONITOR
867  timerPoly_ = Teuchos::TimeMonitor::getNewCounter(polyLabel);
868 #endif
869  }
870 
871  // Inform the solver manager that the current parameters were set.
872  isSet_ = true;
873 }
874 
875 // Check the status test versus the defined linear problem
876 template<class ScalarType, class MV, class OP>
878 
879  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
880  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestGenResNorm_t;
881  typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP> StatusTestImpResNorm_t;
882 
883  // Basic test checks maximum iterations and native residual.
884  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
885 
886  // If there is a left preconditioner, we create a combined status test that checks the implicit
887  // and then explicit residual norm to see if we have convergence.
888  if (!Teuchos::is_null(problem_->getLeftPrec())) {
889  expResTest_ = true;
890  }
891 
892  if (expResTest_) {
893 
894  // Implicit residual test, using the native residual to determine if convergence was achieved.
895  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
896  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
897  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
898  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
899  impConvTest_ = tmpImpConvTest;
900 
901  // Explicit residual test once the native residual is below the tolerance
902  Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
903  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
904  tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
905  tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
906  tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
907  expConvTest_ = tmpExpConvTest;
908 
909  // The convergence test is a combination of the "cheap" implicit test and explicit test.
910  convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
911  }
912  else {
913 
914  // Implicit residual test, using the native residual to determine if convergence was achieved.
915  // Use test that checks for loss of accuracy.
916  Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
917  Teuchos::rcp( new StatusTestImpResNorm_t( convtol_ ) );
918  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
919  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
920  impConvTest_ = tmpImpConvTest;
921 
922  // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
923  expConvTest_ = impConvTest_;
924  convTest_ = impConvTest_;
925  }
926 
927  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
928 
929  // Create the status test output class.
930  // This class manages and formats the output from the status test.
931  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
932  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
933 
934  // Set the solver string for the output test
935  std::string solverDesc = " Gmres Polynomial ";
936  outputTest_->setSolverDesc( solverDesc );
937 
938 
939  // The status test is now set.
940  isSTSet_ = true;
941 
942  return false;
943 }
944 
945 template<class ScalarType, class MV, class OP>
947 {
948  // Create a copy of the linear problem that has a zero initial guess and random RHS.
949  Teuchos::RCP<MV> newX = MVT::Clone( *(problem_->getLHS()), 1 );
950  Teuchos::RCP<MV> newB = MVT::Clone( *(problem_->getRHS()), 1 );
951  MVT::MvInit( *newX, STS::zero() );
952  MVT::MvRandom( *newB );
953  Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > newProblem =
954  Teuchos::rcp( new LinearProblem<ScalarType,MV,OP>( problem_->getOperator(), newX, newB ) );
955  newProblem->setLeftPrec( problem_->getLeftPrec() );
956  newProblem->setRightPrec( problem_->getRightPrec() );
957  newProblem->setLabel("Belos GMRES Poly Generation");
958  newProblem->setProblem();
959  std::vector<int> idx(1,0); // Must set the index to be the first vector (0)!
960  newProblem->setLSIndex( idx );
961 
962  // Create a parameter list for the GMRES iteration.
963  Teuchos::ParameterList polyList;
964 
965  // Tell the block solver that the block size is one.
966  polyList.set("Num Blocks",maxDegree_);
967  polyList.set("Block Size",1);
968  polyList.set("Keep Hessenberg", true);
969 
970  // Create a simple status test that either reaches the relative residual tolerance or maximum polynomial size.
971  Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxItrTst =
972  Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxDegree_ ) );
973 
974  // Implicit residual test, using the native residual to determine if convergence was achieved.
975  Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTst =
976  Teuchos::rcp( new StatusTestGenResNorm<ScalarType,MV,OP>( polytol_ ) );
977  convTst->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
978 
979  // Convergence test that stops the iteration when either are satisfied.
980  Teuchos::RCP<StatusTestCombo<ScalarType,MV,OP> > polyTest =
981  Teuchos::rcp( new StatusTestCombo<ScalarType,MV,OP>( StatusTestCombo<ScalarType,MV,OP>::OR, maxItrTst, convTst ) );
982 
983  // Create Gmres iteration object to perform one cycle of Gmres.
984  Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > gmres_iter;
985  gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(newProblem,printer_,polyTest,ortho_,polyList) );
986 
987  // Create the first block in the current Krylov basis (residual).
988  Teuchos::RCP<MV> V_0 = MVT::Clone( *(newProblem->getRHS()), 1 );
989  newProblem->computeCurrPrecResVec( &*V_0 );
990 
991  // Get a matrix to hold the orthonormalization coefficients.
992  poly_r0_ = Teuchos::rcp( new Teuchos::SerialDenseVector<int,ScalarType>(1) );
993 
994  // Orthonormalize the new V_0
995  int rank = ortho_->normalize( *V_0, poly_r0_ );
996  TEUCHOS_TEST_FOR_EXCEPTION(rank != 1,GmresPolySolMgrOrthoFailure,
997  "Belos::GmresPolySolMgr::generatePoly(): Failed to compute initial block of orthonormal vectors for polynomial generation.");
998 
999  // Set the new state and initialize the solver.
1001  newstate.V = V_0;
1002  newstate.z = poly_r0_;
1003  newstate.curDim = 0;
1004  gmres_iter->initializeGmres(newstate);
1005 
1006  // Perform Gmres iteration
1007  bool polyConverged = false;
1008  try {
1009  gmres_iter->iterate();
1010 
1011  // Check convergence first
1012  if ( convTst->getStatus() == Passed ) {
1013  // we have convergence
1014  polyConverged = true;
1015  }
1016  }
1017  catch (GmresIterationOrthoFailure e) {
1018  // Try to recover the most recent least-squares solution
1019  gmres_iter->updateLSQR( gmres_iter->getCurSubspaceDim() );
1020 
1021  // Check to see if the most recent least-squares solution yielded convergence.
1022  polyTest->checkStatus( &*gmres_iter );
1023  if (convTst->getStatus() == Passed)
1024  polyConverged = true;
1025  }
1026  catch (std::exception e) {
1027  printer_->stream(Errors) << "Error! Caught exception in BlockGmresIter::iterate() at iteration "
1028  << gmres_iter->getNumIters() << std::endl
1029  << e.what() << std::endl;
1030  throw;
1031  }
1032 
1033  // FIXME (mfh 27 Apr 2013) Why aren't we using polyConverged after
1034  // setting it? I'm tired of the compile warning so I'll silence it
1035  // here, but I'm curious why we aren't using the variable.
1036  (void) polyConverged;
1037 
1038  // Get the solution for this polynomial, use in comparison below
1039  Teuchos::RCP<MV> currX = gmres_iter->getCurrentUpdate();
1040 
1041  // Record polynomial info, get current GMRES state
1042  GmresIterationState<ScalarType,MV> gmresState = gmres_iter->getState();
1043 
1044  // If the polynomial has no dimension, the tolerance is too low, return false
1045  poly_dim_ = gmresState.curDim;
1046  if (poly_dim_ == 0) {
1047  return false;
1048  }
1049  //
1050  // Make a view and then copy the RHS of the least squares problem.
1051  //
1052  poly_y_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::Copy, *gmresState.z, poly_dim_, 1 ) );
1053  poly_H_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( *gmresState.H ) );
1054  //
1055  // Solve the least squares problem.
1056  //
1057  const ScalarType one = STS::one ();
1058  Teuchos::BLAS<int,ScalarType> blas;
1059  blas.TRSM( Teuchos::LEFT_SIDE, Teuchos::UPPER_TRI, Teuchos::NO_TRANS,
1060  Teuchos::NON_UNIT_DIAG, poly_dim_, 1, one,
1061  gmresState.R->values(), gmresState.R->stride(),
1062  poly_y_->values(), poly_y_->stride() );
1063  //
1064  // Generate the polynomial operator
1065  //
1066  poly_Op_ = Teuchos::rcp(
1067  new Belos::GmresPolyOp<ScalarType,MV,OP>( problem_, poly_H_, poly_y_, poly_r0_ ) );
1068 
1069  return true;
1070 }
1071 
1072 
1073 template<class ScalarType, class MV, class OP>
1075 {
1076  using Teuchos::RCP;
1077  using Teuchos::rcp;
1078  typedef Teuchos::SerialDenseMatrix<int, ScalarType> SDM;
1079 
1080  // Set the current parameters if they were not set before. NOTE:
1081  // This may occur if the user generated the solver manager with the
1082  // default constructor and then didn't set any parameters using
1083  // setParameters().
1084  if (! isSet_) {
1085  setParameters (Teuchos::parameterList (*getValidParameters ()));
1086  }
1087 
1088  TEUCHOS_TEST_FOR_EXCEPTION(
1089  problem_.is_null (), GmresPolySolMgrLinearProblemFailure,
1090  "Belos::GmresPolySolMgr::solve: The linear problem has not been set yet, "
1091  "or was set to null. Please call setProblem() with a nonnull input before "
1092  "calling solve().");
1093 
1094  TEUCHOS_TEST_FOR_EXCEPTION(
1095  ! problem_->isProblemSet (), GmresPolySolMgrLinearProblemFailure,
1096  "Belos::GmresPolySolMgr::solve: The linear problem is not ready. Please "
1097  "call setProblem() on the LinearProblem object before calling solve().");
1098 
1099  if (! isSTSet_ || (! expResTest_ && ! problem_->getLeftPrec ().is_null ())) {
1100  // checkStatusTest() shouldn't have side effects, but it's still
1101  // not such a good idea to put possibly side-effect-y function
1102  // calls in a macro invocation. It could be expensive if the
1103  // macro evaluates the argument more than once, for example.
1104  const bool check = checkStatusTest ();
1105  TEUCHOS_TEST_FOR_EXCEPTION(
1107  "Belos::GmresPolySolMgr::solve: Linear problem and requested status "
1108  "tests are incompatible.");
1109  }
1110 
1111  // If the GMRES polynomial has not been constructed for this
1112  // (nmatrix, preconditioner) pair, generate it.
1113  if (! isPolyBuilt_) {
1114 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1115  Teuchos::TimeMonitor slvtimer(*timerPoly_);
1116 #endif
1117  isPolyBuilt_ = generatePoly();
1118  TEUCHOS_TEST_FOR_EXCEPTION( !isPolyBuilt_ && poly_dim_==0, GmresPolySolMgrPolynomialFailure,
1119  "Belos::GmresPolySolMgr::generatePoly(): Failed to generate a non-trivial polynomial, reduce polynomial tolerance.");
1120  TEUCHOS_TEST_FOR_EXCEPTION( !isPolyBuilt_, GmresPolySolMgrPolynomialFailure,
1121  "Belos::GmresPolySolMgr::generatePoly(): Failed to generate polynomial that satisfied requirements.");
1122  }
1123 
1124  // Assume convergence is achieved if user does not require strict convergence.
1125  bool isConverged = true;
1126 
1127  // Solve the linear system using the polynomial
1128  {
1129 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1130  Teuchos::TimeMonitor slvtimer(*timerSolve_);
1131 #endif
1132 
1133  // Apply the polynomial to the current linear system
1134  poly_Op_->Apply( *problem_->getRHS(), *problem_->getLHS() );
1135 
1136  // Reset the problem to acknowledge the updated solution
1137  problem_->setProblem ();
1138 
1139  // If we have to strictly adhere to the requested convergence tolerance, set up a standard GMRES solver.
1140  if (strictConvTol_) {
1141 
1142  // Create indices for the linear systems to be solved.
1143  int startPtr = 0;
1144  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
1145  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1146 
1147 
1148  // If an adaptive block size is allowed then only the linear
1149  // systems that need to be solved are solved. Otherwise, the
1150  // index set is generated that informs the linear problem that
1151  // some linear systems are augmented.
1152 
1153  std::vector<int> currIdx (blockSize_);
1154  for (int i = 0; i < numCurrRHS; ++i) {
1155  currIdx[i] = startPtr+i;
1156  }
1157  for (int i = numCurrRHS; i < blockSize_; ++i) {
1158  currIdx[i] = -1;
1159  }
1160 
1161  // Inform the linear problem of the current linear system to solve.
1162  problem_->setLSIndex (currIdx);
1163 
1165  // Parameter list
1166  Teuchos::ParameterList plist;
1167  plist.set ("Block Size", blockSize_);
1168 
1169  ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
1170  if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {
1171  ptrdiff_t tmpNumBlocks = 0;
1172  if (blockSize_ == 1) {
1173  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
1174  } else {
1175  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
1176  }
1177  printer_->stream(Warnings)
1178  << "Warning! Requested Krylov subspace dimension is larger than "
1179  << "operator dimension!" << std::endl << "The maximum number of "
1180  << "blocks allowed for the Krylov subspace will be adjusted to "
1181  << tmpNumBlocks << std::endl;
1182  plist.set ("Num Blocks", Teuchos::asSafe<int> (tmpNumBlocks));
1183  }
1184  else {
1185  plist.set ("Num Blocks", numBlocks_);
1186  }
1187 
1188  // Reset the status test.
1189  outputTest_->reset ();
1190  loaDetected_ = false;
1191 
1192  // Assume convergence is achieved, then let any failed
1193  // convergence set this to false.
1194  isConverged = true;
1195 
1196  //
1197  // Solve using BlockGmres
1198  //
1199 
1200  RCP<GmresIteration<ScalarType,MV,OP> > block_gmres_iter =
1201  rcp (new BlockGmresIter<ScalarType,MV,OP> (problem_, printer_,
1202  outputTest_, ortho_, plist));
1203 
1204  // Enter solve() iterations
1205  while (numRHS2Solve > 0) {
1206  // Set the current number of blocks and block size with the
1207  // Gmres iteration.
1208  if (blockSize_*numBlocks_ > dim) {
1209  int tmpNumBlocks = 0;
1210  if (blockSize_ == 1) {
1211  // Leave space for happy breakdown.
1212  tmpNumBlocks = dim / blockSize_;
1213  } else {
1214  // Leave space for restarting.
1215  tmpNumBlocks = (dim - blockSize_) / blockSize_;
1216  }
1217  block_gmres_iter->setSize (blockSize_, tmpNumBlocks);
1218  }
1219  else {
1220  block_gmres_iter->setSize (blockSize_, numBlocks_);
1221  }
1222 
1223  // Reset the number of iterations.
1224  block_gmres_iter->resetNumIters ();
1225 
1226  // Reset the number of calls that the status test output knows about.
1227  outputTest_->resetNumCalls ();
1228 
1229  // Create the first block in the current Krylov basis.
1230  RCP<MV> V_0 = MVT::CloneCopy (*(problem_->getInitPrecResVec ()), currIdx);
1231 
1232  // Get a matrix to hold the orthonormalization coefficients.
1233  RCP<SDM> z_0 = rcp (new SDM (blockSize_, blockSize_));
1234 
1235  // Orthonormalize the new V_0
1236  int rank = ortho_->normalize (*V_0, z_0);
1237  TEUCHOS_TEST_FOR_EXCEPTION(
1238  rank != blockSize_, GmresPolySolMgrOrthoFailure,
1239  "Belos::GmresPolySolMgr::solve: Failed to compute initial block of "
1240  "orthonormal vectors.");
1241 
1242  // Set the new state and initialize the solver.
1244  newstate.V = V_0;
1245  newstate.z = z_0;
1246  newstate.curDim = 0;
1247  block_gmres_iter->initializeGmres(newstate);
1248  int numRestarts = 0;
1249 
1250  while(1) {
1251  try {
1252  block_gmres_iter->iterate();
1253 
1254  // check convergence first
1255  if ( convTest_->getStatus() == Passed ) {
1256  if ( expConvTest_->getLOADetected() ) {
1257  // we don't have convergence
1258  loaDetected_ = true;
1259  isConverged = false;
1260  }
1261  break; // break from while(1){block_gmres_iter->iterate()}
1262  }
1263 
1264  // check for maximum iterations
1265  else if ( maxIterTest_->getStatus() == Passed ) {
1266  // we don't have convergence
1267  isConverged = false;
1268  break; // break from while(1){block_gmres_iter->iterate()}
1269  }
1270  // check for restarting, i.e. the subspace is full
1271  else if (block_gmres_iter->getCurSubspaceDim () ==
1272  block_gmres_iter->getMaxSubspaceDim ()) {
1273  if (numRestarts >= maxRestarts_) {
1274  isConverged = false;
1275  break; // break from while(1){block_gmres_iter->iterate()}
1276  }
1277  numRestarts++;
1278 
1279  printer_->stream(Debug)
1280  << " Performing restart number " << numRestarts << " of "
1281  << maxRestarts_ << std::endl << std::endl;
1282 
1283  // Update the linear problem.
1284  RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1285  problem_->updateSolution (update, true);
1286 
1287  // Get the state.
1288  GmresIterationState<ScalarType,MV> oldState = block_gmres_iter->getState();
1289 
1290  // Compute the restart vector.
1291  // Get a view of the current Krylov basis.
1292  //
1293  // We call this VV_0 to avoid shadowing the previously
1294  // defined V_0 above.
1295  RCP<MV> VV_0 = MVT::Clone (*(oldState.V), blockSize_);
1296  problem_->computeCurrPrecResVec (&*VV_0);
1297 
1298  // Get a view of the first block of the Krylov basis.
1299  //
1300  // We call this zz_0 to avoid shadowing the previously
1301  // defined z_0 above.
1302  RCP<SDM> zz_0 = rcp (new SDM (blockSize_, blockSize_));
1303 
1304  // Orthonormalize the new VV_0
1305  const int theRank = ortho_->normalize( *VV_0, zz_0 );
1306  TEUCHOS_TEST_FOR_EXCEPTION(
1307  theRank != blockSize_, GmresPolySolMgrOrthoFailure,
1308  "Belos::GmresPolySolMgr::solve: Failed to compute initial "
1309  "block of orthonormal vectors after restart.");
1310 
1311  // Set the new state and initialize the solver.
1313  theNewState.V = VV_0;
1314  theNewState.z = zz_0;
1315  theNewState.curDim = 0;
1316  block_gmres_iter->initializeGmres (theNewState);
1317  } // end of restarting
1318  //
1319  // We returned from iterate(), but none of our status
1320  // tests Passed. Something is wrong, and it is probably
1321  // our fault.
1322  //
1323  else {
1324  TEUCHOS_TEST_FOR_EXCEPTION(
1325  true, std::logic_error,
1326  "Belos::GmresPolySolMgr::solve: Invalid return from "
1327  "BlockGmresIter::iterate(). Please report this bug "
1328  "to the Belos developers.");
1329  }
1330  }
1331  catch (const GmresIterationOrthoFailure& e) {
1332  // If the block size is not one, it's not considered a lucky breakdown.
1333  if (blockSize_ != 1) {
1334  printer_->stream(Errors)
1335  << "Error! Caught std::exception in BlockGmresIter::iterate() "
1336  << "at iteration " << block_gmres_iter->getNumIters()
1337  << std::endl << e.what() << std::endl;
1338  if (convTest_->getStatus() != Passed) {
1339  isConverged = false;
1340  }
1341  break;
1342  }
1343  else {
1344  // If the block size is one, try to recover the most
1345  // recent least-squares solution
1346  block_gmres_iter->updateLSQR (block_gmres_iter->getCurSubspaceDim ());
1347 
1348  // Check to see if the most recent least-squares
1349  // solution yielded convergence.
1350  sTest_->checkStatus (&*block_gmres_iter);
1351  if (convTest_->getStatus() != Passed) {
1352  isConverged = false;
1353  }
1354  break;
1355  }
1356  }
1357  catch (const std::exception &e) {
1358  printer_->stream(Errors)
1359  << "Error! Caught std::exception in BlockGmresIter::iterate() "
1360  << "at iteration " << block_gmres_iter->getNumIters() << std::endl
1361  << e.what() << std::endl;
1362  throw;
1363  }
1364  }
1365 
1366  // Compute the current solution. Update the linear problem.
1367  // Attempt to get the current solution from the residual
1368  // status test, if it has one.
1369  if (! Teuchos::is_null (expConvTest_->getSolution ()) ) {
1370  RCP<MV> newX = expConvTest_->getSolution ();
1371  RCP<MV> curX = problem_->getCurrLHSVec ();
1372  MVT::MvAddMv (STS::zero (), *newX, STS::one (), *newX, *curX);
1373  }
1374  else {
1375  RCP<MV> update = block_gmres_iter->getCurrentUpdate ();
1376  problem_->updateSolution (update, true);
1377  }
1378 
1379  // Inform the linear problem that we are finished with this block linear system.
1380  problem_->setCurrLS ();
1381 
1382  // Update indices for the linear systems to be solved.
1383  startPtr += numCurrRHS;
1384  numRHS2Solve -= numCurrRHS;
1385  if (numRHS2Solve > 0) {
1386  numCurrRHS = (numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1387 
1388  currIdx.resize (blockSize_);
1389  for (int i=0; i<numCurrRHS; ++i) {
1390  currIdx[i] = startPtr+i;
1391  }
1392  for (int i=numCurrRHS; i<blockSize_; ++i) {
1393  currIdx[i] = -1;
1394  }
1395 
1396  // Set the next indices.
1397  problem_->setLSIndex( currIdx );
1398  }
1399  else {
1400  currIdx.resize( numRHS2Solve );
1401  }
1402 
1403  }// while ( numRHS2Solve > 0 )
1404 
1405  // print final summary
1406  sTest_->print( printer_->stream(FinalSummary) );
1407 
1408  } // if (strictConvTol_)
1409  } // timing block
1410 
1411  // print timing information
1412 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1413  // Calling summarize() can be expensive, so don't call unless the
1414  // user wants to print out timing details. summarize() will do all
1415  // the work even if it's passed a "black hole" output stream.
1416  if (verbosity_ & TimingDetails)
1417  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1418 #endif
1419 
1420  if (!isConverged || loaDetected_) {
1421  return Unconverged; // return from GmresPolySolMgr::solve()
1422  }
1423  return Converged; // return from GmresPolySolMgr::solve()
1424 }
1425 
1426 
1427 template<class ScalarType, class MV, class OP>
1429 {
1430  std::ostringstream out;
1431 
1432  out << "\"Belos::GmresPolySolMgr\": {"
1433  << "ScalarType: " << Teuchos::TypeNameTraits<ScalarType>::name ()
1434  << ", Ortho Type: " << orthoType_
1435  << ", Block Size: " << blockSize_
1436  << ", Num Blocks: " << numBlocks_
1437  << ", Max Restarts: " << maxRestarts_;
1438  out << "}";
1439  return out.str ();
1440 }
1441 
1442 } // namespace Belos
1443 
1444 #endif // BELOS_GMRES_POLY_SOLMGR_HPP
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:87
Collection of types and exceptions used within the Belos solvers.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
virtual ~GmresPolySolMgr()
Destructor.
Defines the GMRES polynomial operator hybrid-GMRES iterative linear solver.
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:119
static void MvRandom(MV &mv)
Replace the vectors in mv with random vectors.
This class implements the block GMRES iteration, where a block Krylov subspace is constructed...
Convergence test using the implicit residual norm(s), with an explicit residual norm(s) check for los...
Teuchos::RCP< const MV > V
The current Krylov basis.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
A factory class for generating StatusTestOutput objects.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
An implementation of StatusTestResNorm using a family of residual norms.
Structure to contain pointers to GmresIteration state variables.
Belos::StatusTest class for specifying a maximum number of iterations.
static void MvInit(MV &mv, const ScalarType alpha=Teuchos::ScalarTraits< ScalarType >::zero())
Replace each element of the vectors in mv with alpha.
static int GetNumberVecs(const MV &mv)
Obtain the number of vectors in mv.
A factory class for generating StatusTestOutput objects.
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
void reset(const ResetType type)
Reset the solver.
static void MvAddMv(const ScalarType alpha, const MV &A, const ScalarType beta, const MV &B, MV &mv)
Replace mv with .
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
std::string description() const
Method to return description of the hybrid block GMRES solver manager.
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class...
A Belos::StatusTest class for specifying a maximum number of iterations.
static Teuchos::RCP< MV > Clone(const MV &mv, const int numvecs)
Creates a new empty MV containing numvecs columns.
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:205
GmresPolySolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
Pure virtual base class which describes the basic interface for a solver manager. ...
int curDim
The current dimension of the reduction.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a parameter list containing the valid parameters for this object.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > R
The current upper-triangular matrix from the QR reduction of H.
A linear system to solve, and its associated information.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Class which describes the linear problem to be solved by the iterative solver.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params)
Set the parameters the solver manager should use to solve the linear problem.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
Hybrid block GMRES iterative linear solver.
static Teuchos::RCP< MV > CloneCopy(const MV &mv)
Creates a new MV and copies contents of mv into the new vector (deep copy).
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:154
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
ReturnType solve()
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
Belos&#39;s class for applying the GMRES polynomial operator that is used by the hybrid-GMRES linear solv...
Belos concrete class for performing the block GMRES iteration.
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
static ptrdiff_t GetGlobalLength(const MV &mv)
Return the number of rows in the given multivector mv.
GmresPolySolMgrPolynomialFailure is thrown when their is a problem generating the GMRES polynomial fo...
GmresPolySolMgrPolynomialFailure(const std::string &what_arg)
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > create(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, int mod, int printStates)
Create the StatusTestOutput object specified by the outputStyle.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
A class for extending the status testing capabilities of Belos via logical combinations.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
GmresPolySolMgr()
Empty constructor for GmresPolySolMgr. This constructor takes no arguments and sets the default value...
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration...
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
int getNumIters() const
Get the iteration count for the most recent call to solve().
Belos header file which uses auto-configuration information to include necessary C++ headers...
bool isLOADetected() const
Return whether a loss of accuracy was detected by this solver during the most current solve...
GmresPolySolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthon...
GmresPolySolMgrOrthoFailure(const std::string &what_arg)
GmresPolySolMgrLinearProblemFailure(const std::string &what_arg)
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.

Generated on Mon Feb 5 2018 15:01:51 for Belos by doxygen 1.8.13