Belos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosBlockGmresSolMgr.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #ifndef BELOS_BLOCK_GMRES_SOLMGR_HPP
11 #define BELOS_BLOCK_GMRES_SOLMGR_HPP
12 
17 #include "BelosConfigDefs.hpp"
18 #include "BelosTypes.hpp"
19 
20 #include "BelosLinearProblem.hpp"
21 #include "BelosSolverManager.hpp"
22 
23 #include "BelosGmresIteration.hpp"
24 #include "BelosBlockGmresIter.hpp"
25 #include "BelosBlockFGmresIter.hpp"
30 #include "BelosStatusTestCombo.hpp"
33 #include "BelosOutputManager.hpp"
34 #ifdef BELOS_TEUCHOS_TIME_MONITOR
35 #include "Teuchos_TimeMonitor.hpp"
36 #endif
37 
54 namespace Belos {
55 
57 
58 
66  BlockGmresSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
67  {}};
68 
76  BlockGmresSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
77  {}};
78 
95 template<class ScalarType, class MV, class OP>
96 class BlockGmresSolMgr : public SolverManager<ScalarType,MV,OP> {
97 
98 private:
104 
105 public:
106 
108 
109 
116 
139 
141  virtual ~BlockGmresSolMgr() {};
142 
146  }
148 
150 
151 
154  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
155  return *problem_;
156  }
157 
161 
165 
172  return Teuchos::tuple(timerSolve_);
173  }
174 
185  MagnitudeType achievedTol() const override {
186  return achievedTol_;
187  }
188 
190  int getNumIters() const override {
191  return numIters_;
192  }
193 
197  bool isLOADetected() const override { return loaDetected_; }
198 
200 
202 
203 
205  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; isSTSet_ = false; }
206 
208  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
209 
211  void setDebugStatusTest( const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest ) override;
212 
214 
216 
217 
221  void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
223 
225 
226 
244  ReturnType solve() override;
245 
247 
250 
257  void
259  const Teuchos::EVerbosityLevel verbLevel =
261 
263  std::string description () const override;
264 
266 
267 private:
268 
269  // Method for checking current status test against defined linear problem.
270  bool checkStatusTest();
271 
272  // Linear problem.
274 
275  // Output manager.
278 
279  // Status test.
286 
287  // Orthogonalization manager.
289 
290  // Current parameter list.
292 
293  // Default solver values.
294  static constexpr int maxRestarts_default_ = 20;
295  static constexpr int maxIters_default_ = 1000;
296  static constexpr bool adaptiveBlockSize_default_ = true;
297  static constexpr bool showMaxResNormOnly_default_ = false;
298  static constexpr bool flexibleGmres_default_ = false;
299  static constexpr bool expResTest_default_ = false;
300  static constexpr int blockSize_default_ = 1;
301  static constexpr int numBlocks_default_ = 300;
302  static constexpr int verbosity_default_ = Belos::Errors;
303  static constexpr int outputStyle_default_ = Belos::General;
304  static constexpr int outputFreq_default_ = -1;
305  static constexpr const char * impResScale_default_ = "Norm of Preconditioned Initial Residual";
306  static constexpr const char * expResScale_default_ = "Norm of Initial Residual";
307  static constexpr const char * label_default_ = "Belos";
308  static constexpr const char * orthoType_default_ = "ICGS";
309 
310  // Current solver values.
315  std::string orthoType_;
317 
318  // Timers.
319  std::string label_;
321 
322  // Internal state variables.
325 };
326 
327 
328 // Empty Constructor
329 template<class ScalarType, class MV, class OP>
331  outputStream_(Teuchos::rcpFromRef(std::cout)),
332  convtol_(DefaultSolverParameters::convTol),
333  orthoKappa_(DefaultSolverParameters::orthoKappa),
334  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
335  maxRestarts_(maxRestarts_default_),
336  maxIters_(maxIters_default_),
337  numIters_(0),
338  blockSize_(blockSize_default_),
339  numBlocks_(numBlocks_default_),
340  verbosity_(verbosity_default_),
341  outputStyle_(outputStyle_default_),
342  outputFreq_(outputFreq_default_),
343  adaptiveBlockSize_(adaptiveBlockSize_default_),
344  showMaxResNormOnly_(showMaxResNormOnly_default_),
345  isFlexible_(flexibleGmres_default_),
346  expResTest_(expResTest_default_),
347  orthoType_(orthoType_default_),
348  impResScale_(impResScale_default_),
349  expResScale_(expResScale_default_),
350  label_(label_default_),
351  isSet_(false),
352  isSTSet_(false),
353  loaDetected_(false)
354 {}
355 
356 
357 // Basic Constructor
358 template<class ScalarType, class MV, class OP>
362  problem_(problem),
363  outputStream_(Teuchos::rcpFromRef(std::cout)),
364  convtol_(DefaultSolverParameters::convTol),
365  orthoKappa_(DefaultSolverParameters::orthoKappa),
366  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
367  maxRestarts_(maxRestarts_default_),
368  maxIters_(maxIters_default_),
369  numIters_(0),
370  blockSize_(blockSize_default_),
371  numBlocks_(numBlocks_default_),
372  verbosity_(verbosity_default_),
373  outputStyle_(outputStyle_default_),
374  outputFreq_(outputFreq_default_),
375  adaptiveBlockSize_(adaptiveBlockSize_default_),
376  showMaxResNormOnly_(showMaxResNormOnly_default_),
377  isFlexible_(flexibleGmres_default_),
378  expResTest_(expResTest_default_),
379  orthoType_(orthoType_default_),
380  impResScale_(impResScale_default_),
381  expResScale_(expResScale_default_),
382  label_(label_default_),
383  isSet_(false),
384  isSTSet_(false),
385  loaDetected_(false)
386 {
387 
388  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
389 
390  // If the parameter list pointer is null, then set the current parameters to the default parameter list.
391  if ( !is_null(pl) ) {
392  setParameters( pl );
393  }
394 
395 }
396 
397 
398 template<class ScalarType, class MV, class OP>
401 {
403  if (is_null(validPL)) {
404  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
405 
406  // The static_cast is to resolve an issue with older clang versions which
407  // would cause the constexpr to link fail. With c++17 the problem is resolved.
408  pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
409  "The relative residual tolerance that needs to be achieved by the\n"
410  "iterative solver in order for the linear system to be declared converged." );
411  pl->set("Maximum Restarts", static_cast<int>(maxRestarts_default_),
412  "The maximum number of restarts allowed for each\n"
413  "set of RHS solved.");
414  pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
415  "The maximum number of block iterations allowed for each\n"
416  "set of RHS solved.");
417  pl->set("Num Blocks", static_cast<int>(numBlocks_default_),
418  "The maximum number of blocks allowed in the Krylov subspace\n"
419  "for each set of RHS solved.");
420  pl->set("Block Size", static_cast<int>(blockSize_default_),
421  "The number of vectors in each block. This number times the\n"
422  "number of blocks is the total Krylov subspace dimension.");
423  pl->set("Adaptive Block Size", static_cast<bool>(adaptiveBlockSize_default_),
424  "Whether the solver manager should adapt the block size\n"
425  "based on the number of RHS to solve.");
426  pl->set("Verbosity", static_cast<int>(verbosity_default_),
427  "What type(s) of solver information should be outputted\n"
428  "to the output stream.");
429  pl->set("Output Style", static_cast<int>(outputStyle_default_),
430  "What style is used for the solver information outputted\n"
431  "to the output stream.");
432  pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
433  "How often convergence information should be outputted\n"
434  "to the output stream.");
435  pl->set("Output Stream", Teuchos::rcpFromRef(std::cout),
436  "A reference-counted pointer to the output stream where all\n"
437  "solver output is sent.");
438  pl->set("Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
439  "When convergence information is printed, only show the maximum\n"
440  "relative residual norm when the block size is greater than one.");
441  pl->set("Flexible Gmres", static_cast<bool>(flexibleGmres_default_),
442  "Whether the solver manager should use the flexible variant\n"
443  "of GMRES.");
444  pl->set("Explicit Residual Test", static_cast<bool>(expResTest_default_),
445  "Whether the explicitly computed residual should be used in the convergence test.");
446  pl->set("Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
447  "The type of scaling used in the implicit residual convergence test.");
448  pl->set("Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
449  "The type of scaling used in the explicit residual convergence test.");
450  pl->set("Timer Label", static_cast<const char *>(label_default_),
451  "The string to use as a prefix for the timer labels.");
452  pl->set("Orthogonalization", static_cast<const char *>(orthoType_default_),
453  "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
454  pl->set("Orthogonalization Constant",static_cast<MagnitudeType>(DefaultSolverParameters::orthoKappa),
455  "The constant used by DGKS orthogonalization to determine\n"
456  "whether another step of classical Gram-Schmidt is necessary.");
457  validPL = pl;
458  }
459  return validPL;
460 }
461 
462 
463 template<class ScalarType, class MV, class OP>
465 {
466 
467  // Create the internal parameter list if ones doesn't already exist.
468  if (params_ == Teuchos::null) {
469  params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
470  }
471  else {
472  params->validateParameters(*getValidParameters());
473  }
474 
475  // Check for maximum number of restarts
476  if (params->isParameter("Maximum Restarts")) {
477  maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
478 
479  // Update parameter in our list.
480  params_->set("Maximum Restarts", maxRestarts_);
481  }
482 
483  // Check for maximum number of iterations
484  if (params->isParameter("Maximum Iterations")) {
485  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
486 
487  // Update parameter in our list and in status test.
488  params_->set("Maximum Iterations", maxIters_);
489  if (maxIterTest_!=Teuchos::null)
490  maxIterTest_->setMaxIters( maxIters_ );
491  }
492 
493  // Check for blocksize
494  if (params->isParameter("Block Size")) {
495  blockSize_ = params->get("Block Size",blockSize_default_);
496  TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
497  "Belos::BlockGmresSolMgr: \"Block Size\" must be strictly positive.");
498 
499  // Update parameter in our list.
500  params_->set("Block Size", blockSize_);
501  }
502 
503  // Check if the blocksize should be adaptive
504  if (params->isParameter("Adaptive Block Size")) {
505  adaptiveBlockSize_ = params->get("Adaptive Block Size",adaptiveBlockSize_default_);
506 
507  // Update parameter in our list.
508  params_->set("Adaptive Block Size", adaptiveBlockSize_);
509  }
510 
511  // Check for the maximum number of blocks.
512  if (params->isParameter("Num Blocks")) {
513  numBlocks_ = params->get("Num Blocks",numBlocks_default_);
514  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
515  "Belos::BlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
516 
517  // Update parameter in our list.
518  params_->set("Num Blocks", numBlocks_);
519  }
520 
521  // Check to see if the timer label changed.
522  if (params->isParameter("Timer Label")) {
523  std::string tempLabel = params->get("Timer Label", label_default_);
524 
525  // Update parameter in our list, solver timer, and orthogonalization label
526  if (tempLabel != label_) {
527  label_ = tempLabel;
528  params_->set("Timer Label", label_);
529  std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
530 #ifdef BELOS_TEUCHOS_TIME_MONITOR
531  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
532 #endif
533  if (ortho_ != Teuchos::null) {
534  ortho_->setLabel( label_ );
535  }
536  }
537  }
538 
539  // Determine whether this solver should be "flexible".
540  if (params->isParameter("Flexible Gmres")) {
541  isFlexible_ = Teuchos::getParameter<bool>(*params,"Flexible Gmres");
542  params_->set("Flexible Gmres", isFlexible_);
543  if (isFlexible_ && expResTest_) {
544  // Use an implicit convergence test if the Gmres solver is flexible
545  isSTSet_ = false;
546  }
547  }
548 
549  // Check for a change in verbosity level
550  if (params->isParameter("Verbosity")) {
551  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
552  verbosity_ = params->get("Verbosity", verbosity_default_);
553  } else {
554  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
555  }
556 
557  // Update parameter in our list.
558  params_->set("Verbosity", verbosity_);
559  if (printer_ != Teuchos::null)
560  printer_->setVerbosity(verbosity_);
561  }
562 
563  // Check for a change in output style
564  if (params->isParameter("Output Style")) {
565  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
566  outputStyle_ = params->get("Output Style", outputStyle_default_);
567  } else {
568  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
569  }
570 
571  // Reconstruct the convergence test if the explicit residual test is not being used.
572  params_->set("Output Style", outputStyle_);
573  if (outputTest_ != Teuchos::null) {
574  isSTSet_ = false;
575  }
576  }
577 
578  // output stream
579  if (params->isParameter("Output Stream")) {
580  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
581 
582  // Update parameter in our list.
583  params_->set("Output Stream", outputStream_);
584  if (printer_ != Teuchos::null)
585  printer_->setOStream( outputStream_ );
586  }
587 
588  // frequency level
589  if (verbosity_ & Belos::StatusTestDetails) {
590  if (params->isParameter("Output Frequency")) {
591  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
592  }
593 
594  // Update parameter in out list and output status test.
595  params_->set("Output Frequency", outputFreq_);
596  if (outputTest_ != Teuchos::null)
597  outputTest_->setOutputFrequency( outputFreq_ );
598  }
599 
600  // Create output manager if we need to.
601  if (printer_ == Teuchos::null) {
602  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
603  }
604 
605  // Check if the orthogonalization changed.
606  bool changedOrthoType = false;
607  if (params->isParameter("Orthogonalization")) {
608  std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
609  if (tempOrthoType != orthoType_) {
610  orthoType_ = tempOrthoType;
611  changedOrthoType = true;
612  }
613  }
614  params_->set("Orthogonalization", orthoType_);
615 
616  // Check which orthogonalization constant to use.
617  if (params->isParameter("Orthogonalization Constant")) {
618  if (params->isType<MagnitudeType> ("Orthogonalization Constant")) {
619  orthoKappa_ = params->get ("Orthogonalization Constant",
620  static_cast<MagnitudeType> (DefaultSolverParameters::orthoKappa));
621  }
622  else {
623  orthoKappa_ = params->get ("Orthogonalization Constant",
625  }
626 
627  // Update parameter in our list.
628  params_->set("Orthogonalization Constant",orthoKappa_);
629  if (orthoType_=="DGKS") {
630  if (orthoKappa_ > 0 && ortho_ != Teuchos::null && !changedOrthoType) {
631  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
632  }
633  }
634  }
635 
636  // Create orthogonalization manager if we need to.
637  if (ortho_ == Teuchos::null || changedOrthoType) {
640  if (orthoType_=="DGKS" && orthoKappa_ > 0) {
641  paramsOrtho = Teuchos::rcp(new Teuchos::ParameterList());
642  paramsOrtho->set ("depTol", orthoKappa_ );
643  }
644 
645  ortho_ = factory.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, label_, paramsOrtho);
646  }
647 
648  // Check for convergence tolerance
649  if (params->isParameter("Convergence Tolerance")) {
650  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
651  convtol_ = params->get ("Convergence Tolerance",
652  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
653  }
654  else {
655  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
656  }
657 
658  // Update parameter in our list and residual tests.
659  params_->set("Convergence Tolerance", convtol_);
660  if (impConvTest_ != Teuchos::null)
661  impConvTest_->setTolerance( convtol_ );
662  if (expConvTest_ != Teuchos::null)
663  expConvTest_->setTolerance( convtol_ );
664  }
665 
666  // Check for a change in scaling, if so we need to build new residual tests.
667  if (params->isParameter("Implicit Residual Scaling")) {
668  std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
669 
670  // Only update the scaling if it's different.
671  if (impResScale_ != tempImpResScale) {
672  Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
673  impResScale_ = tempImpResScale;
674 
675  // Update parameter in our list and residual tests
676  params_->set("Implicit Residual Scaling", impResScale_);
677  if (impConvTest_ != Teuchos::null) {
678  try {
679  impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
680  }
681  catch (std::exception& e) {
682  // Make sure the convergence test gets constructed again.
683  isSTSet_ = false;
684  }
685  }
686  }
687  }
688 
689  if (params->isParameter("Explicit Residual Scaling")) {
690  std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
691 
692  // Only update the scaling if it's different.
693  if (expResScale_ != tempExpResScale) {
694  Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
695  expResScale_ = tempExpResScale;
696 
697  // Update parameter in our list and residual tests
698  params_->set("Explicit Residual Scaling", expResScale_);
699  if (expConvTest_ != Teuchos::null) {
700  try {
701  expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
702  }
703  catch (std::exception& e) {
704  // Make sure the convergence test gets constructed again.
705  isSTSet_ = false;
706  }
707  }
708  }
709  }
710 
711  if (params->isParameter("Explicit Residual Test")) {
712  expResTest_ = Teuchos::getParameter<bool>( *params,"Explicit Residual Test" );
713 
714  // Reconstruct the convergence test if the explicit residual test is not being used.
715  params_->set("Explicit Residual Test", expResTest_);
716  if (expConvTest_ == Teuchos::null) {
717  isSTSet_ = false;
718  }
719  }
720 
721  if (params->isParameter("Show Maximum Residual Norm Only")) {
722  showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
723 
724  // Update parameter in our list and residual tests
725  params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
726  if (impConvTest_ != Teuchos::null)
727  impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
728  if (expConvTest_ != Teuchos::null)
729  expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
730  }
731 
732 
733  // Create the timer if we need to.
734  if (timerSolve_ == Teuchos::null) {
735  std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
736 #ifdef BELOS_TEUCHOS_TIME_MONITOR
737  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
738 #endif
739  }
740 
741  // Inform the solver manager that the current parameters were set.
742  isSet_ = true;
743 }
744 
745 // Check the status test versus the defined linear problem
746 template<class ScalarType, class MV, class OP>
748 
749  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
750  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestGenResNorm_t;
751  typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP> StatusTestImpResNorm_t;
752 
753  // Basic test checks maximum iterations and native residual.
754  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
755 
756  // Perform sanity checking for flexible Gmres here.
757  // NOTE: If the user requests that the solver manager use flexible GMRES, but there is no right preconditioner, don't use flexible GMRES.
758  // Throw an error is the user provided a left preconditioner, as that is inconsistent with flexible GMRES.
759  if (isFlexible_ && Teuchos::is_null(problem_->getRightPrec())) {
760  isFlexible_ = false;
761  params_->set("Flexible Gmres", isFlexible_);
762 
763  // If the user specified the preconditioner as a left preconditioner, throw an error.
765  "Belos::BlockGmresSolMgr::solve(): Linear problem has a left preconditioner, not a right preconditioner, which is required for flexible GMRES.");
766  }
767 
768  // If there is a left preconditioner, we create a combined status test that checks the implicit
769  // and then explicit residual norm to see if we have convergence.
770  if (!Teuchos::is_null(problem_->getLeftPrec()) && !isFlexible_) {
771  expResTest_ = true;
772  }
773 
774  if (expResTest_) {
775 
776  // Implicit residual test, using the native residual to determine if convergence was achieved.
777  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
778  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
779  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
780  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
781  impConvTest_ = tmpImpConvTest;
782 
783  // Explicit residual test once the native residual is below the tolerance
784  Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
785  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
786  tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
787  tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
788  tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
789  expConvTest_ = tmpExpConvTest;
790 
791  // The convergence test is a combination of the "cheap" implicit test and explicit test.
792  convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
793  }
794  else {
795 
796  if (isFlexible_) {
797  // Implicit residual test, using the native residual to determine if convergence was achieved.
798  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
799  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
800  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
801  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
802  impConvTest_ = tmpImpConvTest;
803  }
804  else {
805  // Implicit residual test, using the native residual to determine if convergence was achieved.
806  // Use test that checks for loss of accuracy.
807  Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
808  Teuchos::rcp( new StatusTestImpResNorm_t( convtol_ ) );
809  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
810  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
811  impConvTest_ = tmpImpConvTest;
812  }
813 
814  // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
815  expConvTest_ = impConvTest_;
816  convTest_ = impConvTest_;
817  }
818 
819  // Create the status test.
820  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
821 
822  // Add debug status test, if one is provided by the user
823  if (nonnull(debugStatusTest_) ) {
824  // Add debug convergence test
825  Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(sTest_)->addStatusTest( debugStatusTest_ );
826  }
827 
828  // Create the status test output class.
829  // This class manages and formats the output from the status test.
830  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
831  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
832 
833  // Set the solver string for the output test
834  std::string solverDesc = " Block Gmres ";
835  if (isFlexible_)
836  solverDesc = "Flexible" + solverDesc;
837  outputTest_->setSolverDesc( solverDesc );
838 
839  // The status test is now set.
840  isSTSet_ = true;
841 
842  return false;
843 }
844 
845 template<class ScalarType, class MV, class OP>
847  const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest
848  )
849 {
850  debugStatusTest_ = debugStatusTest;
851 }
852 
853 
854 // solve()
855 template<class ScalarType, class MV, class OP>
857 
858  // Set the current parameters if they were not set before.
859  // NOTE: This may occur if the user generated the solver manager with the default constructor and
860  // then didn't set any parameters using setParameters().
861  if (!isSet_) {
862  setParameters(Teuchos::parameterList(*getValidParameters()));
863  }
864 
866  "Belos::BlockGmresSolMgr::solve(): Linear problem is not a valid object.");
867 
869  "Belos::BlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
870 
871  if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
873  "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
874  }
875 
876  // Create indices for the linear systems to be solved.
877  int startPtr = 0;
878  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
879  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
880 
881  std::vector<int> currIdx;
882  // If an adaptive block size is allowed then only the linear systems that need to be solved are solved.
883  // Otherwise, the index set is generated that informs the linear problem that some linear systems are augmented.
884  if ( adaptiveBlockSize_ ) {
885  blockSize_ = numCurrRHS;
886  currIdx.resize( numCurrRHS );
887  for (int i=0; i<numCurrRHS; ++i)
888  { currIdx[i] = startPtr+i; }
889 
890  }
891  else {
892  currIdx.resize( blockSize_ );
893  for (int i=0; i<numCurrRHS; ++i)
894  { currIdx[i] = startPtr+i; }
895  for (int i=numCurrRHS; i<blockSize_; ++i)
896  { currIdx[i] = -1; }
897  }
898 
899  // Inform the linear problem of the current linear system to solve.
900  problem_->setLSIndex( currIdx );
901 
903  // Parameter list
905  plist.set("Block Size",blockSize_);
906 
907  ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
908  if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {
909  int tmpNumBlocks = 0;
910  if (blockSize_ == 1)
911  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
912  else
913  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
914  printer_->stream(Warnings) <<
915  "Belos::BlockGmresSolMgr::solve(): Warning! Requested Krylov subspace dimension is larger than operator dimension!"
916  << std::endl << " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
917  plist.set("Num Blocks",tmpNumBlocks);
918  }
919  else
920  plist.set("Num Blocks",numBlocks_);
921 
922  // Reset the status test.
923  outputTest_->reset();
924  loaDetected_ = false;
925 
926  // Assume convergence is achieved, then let any failed convergence set this to false.
927  bool isConverged = true;
928 
930  // BlockGmres solver
931 
933 
934  if (isFlexible_)
935  block_gmres_iter = Teuchos::rcp( new BlockFGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
936  else
937  block_gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
938 
939  // Enter solve() iterations
940  {
941 #ifdef BELOS_TEUCHOS_TIME_MONITOR
942  Teuchos::TimeMonitor slvtimer(*timerSolve_);
943 #endif
944 
945  while ( numRHS2Solve > 0 ) {
946 
947  // Set the current number of blocks and blocksize with the Gmres iteration.
948  if (blockSize_*numBlocks_ > dim) {
949  int tmpNumBlocks = 0;
950  if (blockSize_ == 1)
951  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
952  else
953  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
954  block_gmres_iter->setSize( blockSize_, tmpNumBlocks );
955  }
956  else
957  block_gmres_iter->setSize( blockSize_, numBlocks_ );
958 
959  // Reset the number of iterations.
960  block_gmres_iter->resetNumIters();
961 
962  // Reset the number of calls that the status test output knows about.
963  outputTest_->resetNumCalls();
964 
965  // Create the first block in the current Krylov basis.
966  Teuchos::RCP<MV> V_0;
967  if (isFlexible_) {
968  // Load the correct residual if the system is augmented
969  if (currIdx[blockSize_-1] == -1) {
970  V_0 = MVT::Clone( *(problem_->getInitResVec()), blockSize_ );
971  problem_->computeCurrResVec( &*V_0 );
972  }
973  else {
974  V_0 = MVT::CloneCopy( *(problem_->getInitResVec()), currIdx );
975  }
976  }
977  else {
978  // Load the correct residual if the system is augmented
979  if (currIdx[blockSize_-1] == -1) {
980  V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
981  problem_->computeCurrPrecResVec( &*V_0 );
982  }
983  else {
984  V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
985  }
986  }
987 
988  // Get a matrix to hold the orthonormalization coefficients.
990  Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
991 
992  // Orthonormalize the new V_0
993  int rank = ortho_->normalize( *V_0, z_0 );
995  "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
996 
997  // Set the new state and initialize the solver.
999  newstate.V = V_0;
1000  newstate.z = z_0;
1001  newstate.curDim = 0;
1002  block_gmres_iter->initializeGmres(newstate);
1003  int numRestarts = 0;
1004 
1005  while(1) {
1006  // tell block_gmres_iter to iterate
1007  try {
1008  block_gmres_iter->iterate();
1009 
1011  //
1012  // check convergence first
1013  //
1015  if ( convTest_->getStatus() == Passed ) {
1016  if ( expConvTest_->getLOADetected() ) {
1017  // we don't have convergence
1018  loaDetected_ = true;
1019  printer_->stream(Warnings) <<
1020  "Belos::BlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1021  isConverged = false;
1022  }
1023  break; // break from while(1){block_gmres_iter->iterate()}
1024  }
1026  //
1027  // check for maximum iterations
1028  //
1030  else if ( maxIterTest_->getStatus() == Passed ) {
1031  // we don't have convergence
1032  isConverged = false;
1033  break; // break from while(1){block_gmres_iter->iterate()}
1034  }
1036  //
1037  // check for restarting, i.e. the subspace is full
1038  //
1040  else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
1041 
1042  if ( numRestarts >= maxRestarts_ ) {
1043  isConverged = false;
1044  break; // break from while(1){block_gmres_iter->iterate()}
1045  }
1046  numRestarts++;
1047 
1048  printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
1049 
1050  // Update the linear problem.
1051  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1052  if (isFlexible_) {
1053  // Update the solution manually, since the preconditioning doesn't need to be undone.
1054  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1055  MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1056  }
1057  else
1058  problem_->updateSolution( update, true );
1059 
1060  // Get the state.
1061  GmresIterationState<ScalarType,MV> oldState = block_gmres_iter->getState();
1062 
1063  // Compute the restart std::vector.
1064  // Get a view of the current Krylov basis.
1065  V_0 = MVT::Clone( *(oldState.V), blockSize_ );
1066  if (isFlexible_)
1067  problem_->computeCurrResVec( &*V_0 );
1068  else
1069  problem_->computeCurrPrecResVec( &*V_0 );
1070 
1071  // Get a view of the first block of the Krylov basis.
1072  z_0 = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1073 
1074  // Orthonormalize the new V_0
1075  rank = ortho_->normalize( *V_0, z_0 );
1077  "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
1078 
1079  // Set the new state and initialize the solver.
1080  newstate.V = V_0;
1081  newstate.z = z_0;
1082  newstate.curDim = 0;
1083  block_gmres_iter->initializeGmres(newstate);
1084 
1085  } // end of restarting
1086 
1088  //
1089  // we returned from iterate(), but none of our status tests Passed.
1090  // something is wrong, and it is probably our fault.
1091  //
1093 
1094  else {
1095  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
1096  "Belos::BlockGmresSolMgr::solve(): Invalid return from BlockGmresIter::iterate().");
1097  }
1098  }
1099  catch (const GmresIterationOrthoFailure &e) {
1100  // If the block size is not one, it's not considered a lucky breakdown.
1101  if (blockSize_ != 1) {
1102  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1103  << block_gmres_iter->getNumIters() << std::endl
1104  << e.what() << std::endl;
1105  if (convTest_->getStatus() != Passed)
1106  isConverged = false;
1107  break;
1108  }
1109  else {
1110  // If the block size is one, try to recover the most recent least-squares solution
1111  block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
1112 
1113  // Check to see if the most recent least-squares solution yielded convergence.
1114  sTest_->checkStatus( &*block_gmres_iter );
1115  if (convTest_->getStatus() != Passed)
1116  isConverged = false;
1117  break;
1118  }
1119  }
1120  catch (const StatusTestNaNError& e) {
1121  // A NaN was detected in the solver. Set the solution to zero and return unconverged.
1122  achievedTol_ = MT::one();
1123  Teuchos::RCP<MV> X = problem_->getLHS();
1124  MVT::MvInit( *X, SCT::zero() );
1125  printer_->stream(Warnings) << "Belos::BlockGmresSolMgr::solve(): Warning! NaN has been detected!"
1126  << std::endl;
1127  return Unconverged;
1128  }
1129  catch (const std::exception &e) {
1130  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1131  << block_gmres_iter->getNumIters() << std::endl
1132  << e.what() << std::endl;
1133  throw;
1134  }
1135  }
1136 
1137  // Compute the current solution.
1138  // Update the linear problem.
1139  if (isFlexible_) {
1140  // Update the solution manually, since the preconditioning doesn't need to be undone.
1141  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1142  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1143  // Update the solution only if there is a valid update from the iteration
1144  if (update != Teuchos::null)
1145  MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1146  }
1147  else {
1148  // Attempt to get the current solution from the residual status test, if it has one.
1149  if ( !Teuchos::is_null(expConvTest_->getSolution()) ) {
1150  Teuchos::RCP<MV> newX = expConvTest_->getSolution();
1151  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1152  MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
1153  }
1154  else {
1155  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1156  problem_->updateSolution( update, true );
1157  }
1158  }
1159 
1160  // Inform the linear problem that we are finished with this block linear system.
1161  problem_->setCurrLS();
1162 
1163  // Update indices for the linear systems to be solved.
1164  startPtr += numCurrRHS;
1165  numRHS2Solve -= numCurrRHS;
1166  if ( numRHS2Solve > 0 ) {
1167  numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1168 
1169  if ( adaptiveBlockSize_ ) {
1170  blockSize_ = numCurrRHS;
1171  currIdx.resize( numCurrRHS );
1172  for (int i=0; i<numCurrRHS; ++i)
1173  { currIdx[i] = startPtr+i; }
1174  }
1175  else {
1176  currIdx.resize( blockSize_ );
1177  for (int i=0; i<numCurrRHS; ++i)
1178  { currIdx[i] = startPtr+i; }
1179  for (int i=numCurrRHS; i<blockSize_; ++i)
1180  { currIdx[i] = -1; }
1181  }
1182  // Set the next indices.
1183  problem_->setLSIndex( currIdx );
1184  }
1185  else {
1186  currIdx.resize( numRHS2Solve );
1187  }
1188 
1189  }// while ( numRHS2Solve > 0 )
1190 
1191  }
1192 
1193  // print final summary
1194  sTest_->print( printer_->stream(FinalSummary) );
1195 
1196  // print timing information
1197 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1198  // Calling summarize() can be expensive, so don't call unless the
1199  // user wants to print out timing details. summarize() will do all
1200  // the work even if it's passed a "black hole" output stream.
1201  if (verbosity_ & TimingDetails)
1202  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1203 #endif
1204 
1205  // get iteration information for this solve
1206  numIters_ = maxIterTest_->getNumIters();
1207 
1208  // Save the convergence test value ("achieved tolerance") for this
1209  // solve. This requires a bit more work than for BlockCGSolMgr,
1210  // since for this solver, convTest_ may either be a single residual
1211  // norm test, or a combination of two residual norm tests. In the
1212  // latter case, the master convergence test convTest_ is a SEQ combo
1213  // of the implicit resp. explicit tests. If the implicit test never
1214  // passes, then the explicit test won't ever be executed. This
1215  // manifests as expConvTest_->getTestValue()->size() < 1. We deal
1216  // with this case by using the values returned by
1217  // impConvTest_->getTestValue().
1218  {
1219  // We'll fetch the vector of residual norms one way or the other.
1220  const std::vector<MagnitudeType>* pTestValues = NULL;
1221  if (expResTest_) {
1222  pTestValues = expConvTest_->getTestValue();
1223  if (pTestValues == NULL || pTestValues->size() < 1) {
1224  pTestValues = impConvTest_->getTestValue();
1225  }
1226  }
1227  else {
1228  // Only the implicit residual norm test is being used.
1229  pTestValues = impConvTest_->getTestValue();
1230  }
1231  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1232  "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1233  "getTestValue() method returned NULL. Please report this bug to the "
1234  "Belos developers.");
1235  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1236  "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1237  "getTestValue() method returned a vector of length zero. Please report "
1238  "this bug to the Belos developers.");
1239 
1240  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1241  // achieved tolerances for all vectors in the current solve(), or
1242  // just for the vectors from the last deflation?
1243  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1244  }
1245 
1246  if (!isConverged || loaDetected_) {
1247  return Unconverged; // return from BlockGmresSolMgr::solve()
1248  }
1249  return Converged; // return from BlockGmresSolMgr::solve()
1250 }
1251 
1252 
1253 template<class ScalarType, class MV, class OP>
1255 {
1256  std::ostringstream out;
1257  out << "\"Belos::BlockGmresSolMgr\": {";
1258  if (this->getObjectLabel () != "") {
1259  out << "Label: " << this->getObjectLabel () << ", ";
1260  }
1261  out << "Flexible: " << (isFlexible_ ? "true" : "false")
1262  << ", Num Blocks: " << numBlocks_
1263  << ", Maximum Iterations: " << maxIters_
1264  << ", Maximum Restarts: " << maxRestarts_
1265  << ", Convergence Tolerance: " << convtol_
1266  << "}";
1267  return out.str ();
1268 }
1269 
1270 
1271 template<class ScalarType, class MV, class OP>
1272 void
1275  const Teuchos::EVerbosityLevel verbLevel) const
1276 {
1278  using Teuchos::VERB_DEFAULT;
1279  using Teuchos::VERB_NONE;
1280  using Teuchos::VERB_LOW;
1281  // using Teuchos::VERB_MEDIUM;
1282  // using Teuchos::VERB_HIGH;
1283  // using Teuchos::VERB_EXTREME;
1284  using std::endl;
1285 
1286  // Set default verbosity if applicable.
1287  const Teuchos::EVerbosityLevel vl =
1288  (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
1289 
1290  if (vl != VERB_NONE) {
1291  Teuchos::OSTab tab0 (out);
1292 
1293  out << "\"Belos::BlockGmresSolMgr\":" << endl;
1294  Teuchos::OSTab tab1 (out);
1295  out << "Template parameters:" << endl;
1296  {
1297  Teuchos::OSTab tab2 (out);
1298  out << "ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
1299  << "MV: " << TypeNameTraits<MV>::name () << endl
1300  << "OP: " << TypeNameTraits<OP>::name () << endl;
1301  }
1302  if (this->getObjectLabel () != "") {
1303  out << "Label: " << this->getObjectLabel () << endl;
1304  }
1305  out << "Flexible: " << (isFlexible_ ? "true" : "false") << endl
1306  << "Num Blocks: " << numBlocks_ << endl
1307  << "Maximum Iterations: " << maxIters_ << endl
1308  << "Maximum Restarts: " << maxRestarts_ << endl
1309  << "Convergence Tolerance: " << convtol_ << endl;
1310  }
1311 }
1312 
1313 
1314 } // end Belos namespace
1315 
1316 #endif /* BELOS_BLOCK_GMRES_SOLMGR_HPP */
static const double orthoKappa
DGKS orthogonalization constant.
Definition: BelosTypes.hpp:267
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:74
Teuchos::RCP< OutputManager< ScalarType > > printer_
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Collection of types and exceptions used within the Belos solvers.
Teuchos::ScalarTraits< MagnitudeType > MT
Teuchos::ScalarTraits< ScalarType >::magnitudeType MagnitudeType
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
static constexpr int maxIters_default_
Class which manages the output and verbosity of the Belos solvers.
bool is_null(const boost::shared_ptr< T > &p)
void setDebugStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debugStatusTest) override
Set a debug status test that will be checked at the same time as the top-level status test...
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > convTest_
static constexpr int blockSize_default_
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > outputTest_
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:88
Teuchos::RCP< Teuchos::Time > timerSolve_
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.
T & get(ParameterList &l, const std::string &name)
Virtual base class for StatusTest that printing status tests.
static RCP< Time > getNewCounter(const std::string &name)
bool is_null(const std::shared_ptr< T > &p)
MagnitudeType achievedTol() const override
Tolerance achieved by the last solve() invocation.
Teuchos::RCP< MatOrthoManager< ScalarType, MV, OP > > ortho_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Get a parameter list containing the valid parameters for this object.
static constexpr int numBlocks_default_
static constexpr bool expResTest_default_
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.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
BlockGmresSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
Structure to contain pointers to GmresIteration state variables.
static const double convTol
Default convergence tolerance.
Definition: BelosTypes.hpp:261
Belos::StatusTest class for specifying a maximum number of iterations.
Teuchos::ScalarTraits< ScalarType > SCT
Interface to Block GMRES and Flexible GMRES.
A pure virtual class for defining the status tests for the Belos iterative solvers.
virtual ~BlockGmresSolMgr()
Destructor.
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > debugStatusTest_
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem) override
Set the linear problem that needs to be solved.
A factory class for generating StatusTestOutput objects.
ReturnType solve() override
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
Traits class which defines basic operations on multivectors.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params) override
Set the parameters the solver manager should use to solve the linear problem.
Belos::StatusTest for logically combining several status tests.
bool isParameter(const std::string &name) const
A Belos::StatusTest class for specifying a maximum number of iterations.
Teuchos::RCP< StatusTestMaxIters< ScalarType, MV, OP > > maxIterTest_
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:174
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Pure virtual base class which describes the basic interface for a solver manager. ...
int curDim
The current dimension of the reduction.
static constexpr const char * impResScale_default_
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > sTest_
BlockGmresSolMgrLinearProblemFailure(const std::string &what_arg)
void reset(const ResetType type) override
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
static constexpr bool flexibleGmres_default_
static constexpr int outputFreq_default_
Teuchos::RCP< std::ostream > outputStream_
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
This class implements the block flexible GMRES iteration, where a block Krylov subspace is constructe...
int getNumIters() const override
Get the iteration count for the most recent call to solve().
static constexpr const char * label_default_
Teuchos::RCP< SolverManager< ScalarType, MV, OP > > clone() const override
clone for Inverted Injection (DII)
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:123
static constexpr int maxRestarts_default_
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
static constexpr int verbosity_default_
Belos concrete class for performing the block GMRES iteration.
void validateParameters(ParameterList const &validParamList, int const depth=1000, EValidateUsed const validateUsed=VALIDATE_USED_ENABLED, EValidateDefaults const validateDefaults=VALIDATE_DEFAULTS_ENABLED) const
static const EVerbosityLevel verbLevel_default
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const override
Print the object with the given verbosity level to a FancyOStream.
static constexpr bool adaptiveBlockSize_default_
static constexpr const char * expResScale_default_
static constexpr const char * orthoType_default_
static constexpr int outputStyle_default_
bool nonnull(const boost::shared_ptr< T > &p)
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.
bool isLOADetected() const override
Return whether a loss of accuracy was detected by this solver during the most current solve...
static constexpr bool showMaxResNormOnly_default_
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
BlockGmresSolMgr()
Empty constructor for BlockGmresSolMgr. This constructor takes no arguments and sets the default valu...
bool isType(const std::string &name) const
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
BlockGmresSolMgrOrthoFailure(const std::string &what_arg)
A class for extending the status testing capabilities of Belos via logical combinations.
Teuchos::RCP< StatusTestResNorm< ScalarType, MV, OP > > expConvTest_
std::string description() const override
Return a one-line description of this object.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:28
BlockGmresSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate ortho...
Default parameters common to most Belos solvers.
Definition: BelosTypes.hpp:251
Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > problem_
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 ...
OperatorTraits< ScalarType, MV, OP > OPT
Belos header file which uses auto-configuration information to include necessary C++ headers...
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const override
Get a parameter list containing the current parameters for this object.
Teuchos::RCP< Teuchos::ParameterList > params_
Teuchos::RCP< StatusTestResNorm< ScalarType, MV, OP > > impConvTest_
const LinearProblem< ScalarType, MV, OP > & getProblem() const override
Get current linear problem being solved for in this object.
Belos concrete class for performing the block, flexible GMRES iteration.
MultiVecTraits< ScalarType, MV > MVT
Teuchos::RCP< Belos::MatOrthoManager< Scalar, MV, OP > > makeMatOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > &params)
Return an instance of the specified MatOrthoManager subclass.