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 //
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_BLOCK_GMRES_SOLMGR_HPP
43 #define BELOS_BLOCK_GMRES_SOLMGR_HPP
44 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosLinearProblem.hpp"
53 #include "BelosSolverManager.hpp"
54 
55 #include "BelosGmresIteration.hpp"
56 #include "BelosBlockGmresIter.hpp"
57 #include "BelosBlockFGmresIter.hpp"
62 #include "BelosStatusTestCombo.hpp"
65 #include "BelosOutputManager.hpp"
66 #ifdef BELOS_TEUCHOS_TIME_MONITOR
67 #include "Teuchos_TimeMonitor.hpp"
68 #endif
69 
80 namespace Belos {
81 
83 
84 
92  BlockGmresSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
93  {}};
94 
102  BlockGmresSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
103  {}};
104 
121 template<class ScalarType, class MV, class OP>
122 class BlockGmresSolMgr : public SolverManager<ScalarType,MV,OP> {
123 
124 private:
130 
131 public:
132 
134 
135 
142 
165 
167  virtual ~BlockGmresSolMgr() {};
168 
172  }
174 
176 
177 
180  const LinearProblem<ScalarType,MV,OP>& getProblem() const override {
181  return *problem_;
182  }
183 
187 
191 
198  return Teuchos::tuple(timerSolve_);
199  }
200 
211  MagnitudeType achievedTol() const override {
212  return achievedTol_;
213  }
214 
216  int getNumIters() const override {
217  return numIters_;
218  }
219 
223  bool isLOADetected() const override { return loaDetected_; }
224 
226 
228 
229 
231  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) override { problem_ = problem; isSTSet_ = false; }
232 
234  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params ) override;
235 
237  void setDebugStatusTest( const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest ) override;
238 
240 
242 
243 
247  void reset( const ResetType type ) override { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
249 
251 
252 
270  ReturnType solve() override;
271 
273 
276 
283  void
285  const Teuchos::EVerbosityLevel verbLevel =
287 
289  std::string description () const override;
290 
292 
293 private:
294 
295  // Method for checking current status test against defined linear problem.
296  bool checkStatusTest();
297 
298  // Linear problem.
300 
301  // Output manager.
304 
305  // Status test.
312 
313  // Orthogonalization manager.
315 
316  // Current parameter list.
318 
319  // Default solver values.
320  static constexpr int maxRestarts_default_ = 20;
321  static constexpr int maxIters_default_ = 1000;
322  static constexpr bool adaptiveBlockSize_default_ = true;
323  static constexpr bool showMaxResNormOnly_default_ = false;
324  static constexpr bool flexibleGmres_default_ = false;
325  static constexpr bool expResTest_default_ = false;
326  static constexpr int blockSize_default_ = 1;
327  static constexpr int numBlocks_default_ = 300;
328  static constexpr int verbosity_default_ = Belos::Errors;
329  static constexpr int outputStyle_default_ = Belos::General;
330  static constexpr int outputFreq_default_ = -1;
331  static constexpr const char * impResScale_default_ = "Norm of Preconditioned Initial Residual";
332  static constexpr const char * expResScale_default_ = "Norm of Initial Residual";
333  static constexpr const char * label_default_ = "Belos";
334  static constexpr const char * orthoType_default_ = "ICGS";
335 
336  // Current solver values.
341  std::string orthoType_;
343 
344  // Timers.
345  std::string label_;
347 
348  // Internal state variables.
351 };
352 
353 
354 // Empty Constructor
355 template<class ScalarType, class MV, class OP>
357  outputStream_(Teuchos::rcpFromRef(std::cout)),
358  convtol_(DefaultSolverParameters::convTol),
359  orthoKappa_(DefaultSolverParameters::orthoKappa),
360  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
361  maxRestarts_(maxRestarts_default_),
362  maxIters_(maxIters_default_),
363  numIters_(0),
364  blockSize_(blockSize_default_),
365  numBlocks_(numBlocks_default_),
366  verbosity_(verbosity_default_),
367  outputStyle_(outputStyle_default_),
368  outputFreq_(outputFreq_default_),
369  adaptiveBlockSize_(adaptiveBlockSize_default_),
370  showMaxResNormOnly_(showMaxResNormOnly_default_),
371  isFlexible_(flexibleGmres_default_),
372  expResTest_(expResTest_default_),
373  orthoType_(orthoType_default_),
374  impResScale_(impResScale_default_),
375  expResScale_(expResScale_default_),
376  label_(label_default_),
377  isSet_(false),
378  isSTSet_(false),
379  loaDetected_(false)
380 {}
381 
382 
383 // Basic Constructor
384 template<class ScalarType, class MV, class OP>
388  problem_(problem),
389  outputStream_(Teuchos::rcpFromRef(std::cout)),
390  convtol_(DefaultSolverParameters::convTol),
391  orthoKappa_(DefaultSolverParameters::orthoKappa),
392  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
393  maxRestarts_(maxRestarts_default_),
394  maxIters_(maxIters_default_),
395  numIters_(0),
396  blockSize_(blockSize_default_),
397  numBlocks_(numBlocks_default_),
398  verbosity_(verbosity_default_),
399  outputStyle_(outputStyle_default_),
400  outputFreq_(outputFreq_default_),
401  adaptiveBlockSize_(adaptiveBlockSize_default_),
402  showMaxResNormOnly_(showMaxResNormOnly_default_),
403  isFlexible_(flexibleGmres_default_),
404  expResTest_(expResTest_default_),
405  orthoType_(orthoType_default_),
406  impResScale_(impResScale_default_),
407  expResScale_(expResScale_default_),
408  label_(label_default_),
409  isSet_(false),
410  isSTSet_(false),
411  loaDetected_(false)
412 {
413 
414  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
415 
416  // If the parameter list pointer is null, then set the current parameters to the default parameter list.
417  if ( !is_null(pl) ) {
418  setParameters( pl );
419  }
420 
421 }
422 
423 
424 template<class ScalarType, class MV, class OP>
427 {
429  if (is_null(validPL)) {
430  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
431 
432  // The static_cast is to resolve an issue with older clang versions which
433  // would cause the constexpr to link fail. With c++17 the problem is resolved.
434  pl->set("Convergence Tolerance", static_cast<MagnitudeType>(DefaultSolverParameters::convTol),
435  "The relative residual tolerance that needs to be achieved by the\n"
436  "iterative solver in order for the linear system to be declared converged." );
437  pl->set("Maximum Restarts", static_cast<int>(maxRestarts_default_),
438  "The maximum number of restarts allowed for each\n"
439  "set of RHS solved.");
440  pl->set("Maximum Iterations", static_cast<int>(maxIters_default_),
441  "The maximum number of block iterations allowed for each\n"
442  "set of RHS solved.");
443  pl->set("Num Blocks", static_cast<int>(numBlocks_default_),
444  "The maximum number of blocks allowed in the Krylov subspace\n"
445  "for each set of RHS solved.");
446  pl->set("Block Size", static_cast<int>(blockSize_default_),
447  "The number of vectors in each block. This number times the\n"
448  "number of blocks is the total Krylov subspace dimension.");
449  pl->set("Adaptive Block Size", static_cast<bool>(adaptiveBlockSize_default_),
450  "Whether the solver manager should adapt the block size\n"
451  "based on the number of RHS to solve.");
452  pl->set("Verbosity", static_cast<int>(verbosity_default_),
453  "What type(s) of solver information should be outputted\n"
454  "to the output stream.");
455  pl->set("Output Style", static_cast<int>(outputStyle_default_),
456  "What style is used for the solver information outputted\n"
457  "to the output stream.");
458  pl->set("Output Frequency", static_cast<int>(outputFreq_default_),
459  "How often convergence information should be outputted\n"
460  "to the output stream.");
461  pl->set("Output Stream", Teuchos::rcpFromRef(std::cout),
462  "A reference-counted pointer to the output stream where all\n"
463  "solver output is sent.");
464  pl->set("Show Maximum Residual Norm Only", static_cast<bool>(showMaxResNormOnly_default_),
465  "When convergence information is printed, only show the maximum\n"
466  "relative residual norm when the block size is greater than one.");
467  pl->set("Flexible Gmres", static_cast<bool>(flexibleGmres_default_),
468  "Whether the solver manager should use the flexible variant\n"
469  "of GMRES.");
470  pl->set("Explicit Residual Test", static_cast<bool>(expResTest_default_),
471  "Whether the explicitly computed residual should be used in the convergence test.");
472  pl->set("Implicit Residual Scaling", static_cast<const char *>(impResScale_default_),
473  "The type of scaling used in the implicit residual convergence test.");
474  pl->set("Explicit Residual Scaling", static_cast<const char *>(expResScale_default_),
475  "The type of scaling used in the explicit residual convergence test.");
476  pl->set("Timer Label", static_cast<const char *>(label_default_),
477  "The string to use as a prefix for the timer labels.");
478  pl->set("Orthogonalization", static_cast<const char *>(orthoType_default_),
479  "The type of orthogonalization to use: DGKS, ICGS, or IMGS.");
480  pl->set("Orthogonalization Constant",static_cast<MagnitudeType>(DefaultSolverParameters::orthoKappa),
481  "The constant used by DGKS orthogonalization to determine\n"
482  "whether another step of classical Gram-Schmidt is necessary.");
483  validPL = pl;
484  }
485  return validPL;
486 }
487 
488 
489 template<class ScalarType, class MV, class OP>
491 {
492 
493  // Create the internal parameter list if ones doesn't already exist.
494  if (params_ == Teuchos::null) {
495  params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
496  }
497  else {
498  params->validateParameters(*getValidParameters());
499  }
500 
501  // Check for maximum number of restarts
502  if (params->isParameter("Maximum Restarts")) {
503  maxRestarts_ = params->get("Maximum Restarts",maxRestarts_default_);
504 
505  // Update parameter in our list.
506  params_->set("Maximum Restarts", maxRestarts_);
507  }
508 
509  // Check for maximum number of iterations
510  if (params->isParameter("Maximum Iterations")) {
511  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
512 
513  // Update parameter in our list and in status test.
514  params_->set("Maximum Iterations", maxIters_);
515  if (maxIterTest_!=Teuchos::null)
516  maxIterTest_->setMaxIters( maxIters_ );
517  }
518 
519  // Check for blocksize
520  if (params->isParameter("Block Size")) {
521  blockSize_ = params->get("Block Size",blockSize_default_);
522  TEUCHOS_TEST_FOR_EXCEPTION(blockSize_ <= 0, std::invalid_argument,
523  "Belos::BlockGmresSolMgr: \"Block Size\" must be strictly positive.");
524 
525  // Update parameter in our list.
526  params_->set("Block Size", blockSize_);
527  }
528 
529  // Check if the blocksize should be adaptive
530  if (params->isParameter("Adaptive Block Size")) {
531  adaptiveBlockSize_ = params->get("Adaptive Block Size",adaptiveBlockSize_default_);
532 
533  // Update parameter in our list.
534  params_->set("Adaptive Block Size", adaptiveBlockSize_);
535  }
536 
537  // Check for the maximum number of blocks.
538  if (params->isParameter("Num Blocks")) {
539  numBlocks_ = params->get("Num Blocks",numBlocks_default_);
540  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
541  "Belos::BlockGmresSolMgr: \"Num Blocks\" must be strictly positive.");
542 
543  // Update parameter in our list.
544  params_->set("Num Blocks", numBlocks_);
545  }
546 
547  // Check to see if the timer label changed.
548  if (params->isParameter("Timer Label")) {
549  std::string tempLabel = params->get("Timer Label", label_default_);
550 
551  // Update parameter in our list, solver timer, and orthogonalization label
552  if (tempLabel != label_) {
553  label_ = tempLabel;
554  params_->set("Timer Label", label_);
555  std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
556 #ifdef BELOS_TEUCHOS_TIME_MONITOR
557  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
558 #endif
559  if (ortho_ != Teuchos::null) {
560  ortho_->setLabel( label_ );
561  }
562  }
563  }
564 
565  // Determine whether this solver should be "flexible".
566  if (params->isParameter("Flexible Gmres")) {
567  isFlexible_ = Teuchos::getParameter<bool>(*params,"Flexible Gmres");
568  params_->set("Flexible Gmres", isFlexible_);
569  if (isFlexible_ && expResTest_) {
570  // Use an implicit convergence test if the Gmres solver is flexible
571  isSTSet_ = false;
572  }
573  }
574 
575  // Check for a change in verbosity level
576  if (params->isParameter("Verbosity")) {
577  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
578  verbosity_ = params->get("Verbosity", verbosity_default_);
579  } else {
580  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
581  }
582 
583  // Update parameter in our list.
584  params_->set("Verbosity", verbosity_);
585  if (printer_ != Teuchos::null)
586  printer_->setVerbosity(verbosity_);
587  }
588 
589  // Check for a change in output style
590  if (params->isParameter("Output Style")) {
591  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
592  outputStyle_ = params->get("Output Style", outputStyle_default_);
593  } else {
594  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
595  }
596 
597  // Reconstruct the convergence test if the explicit residual test is not being used.
598  params_->set("Output Style", outputStyle_);
599  if (outputTest_ != Teuchos::null) {
600  isSTSet_ = false;
601  }
602  }
603 
604  // output stream
605  if (params->isParameter("Output Stream")) {
606  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
607 
608  // Update parameter in our list.
609  params_->set("Output Stream", outputStream_);
610  if (printer_ != Teuchos::null)
611  printer_->setOStream( outputStream_ );
612  }
613 
614  // frequency level
615  if (verbosity_ & Belos::StatusTestDetails) {
616  if (params->isParameter("Output Frequency")) {
617  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
618  }
619 
620  // Update parameter in out list and output status test.
621  params_->set("Output Frequency", outputFreq_);
622  if (outputTest_ != Teuchos::null)
623  outputTest_->setOutputFrequency( outputFreq_ );
624  }
625 
626  // Create output manager if we need to.
627  if (printer_ == Teuchos::null) {
628  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
629  }
630 
631  // Check if the orthogonalization changed.
632  bool changedOrthoType = false;
633  if (params->isParameter("Orthogonalization")) {
634  std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
635  if (tempOrthoType != orthoType_) {
636  orthoType_ = tempOrthoType;
637  changedOrthoType = true;
638  }
639  }
640  params_->set("Orthogonalization", orthoType_);
641 
642  // Check which orthogonalization constant to use.
643  if (params->isParameter("Orthogonalization Constant")) {
644  if (params->isType<MagnitudeType> ("Orthogonalization Constant")) {
645  orthoKappa_ = params->get ("Orthogonalization Constant",
646  static_cast<MagnitudeType> (DefaultSolverParameters::orthoKappa));
647  }
648  else {
649  orthoKappa_ = params->get ("Orthogonalization Constant",
651  }
652 
653  // Update parameter in our list.
654  params_->set("Orthogonalization Constant",orthoKappa_);
655  if (orthoType_=="DGKS") {
656  if (orthoKappa_ > 0 && ortho_ != Teuchos::null && !changedOrthoType) {
657  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
658  }
659  }
660  }
661 
662  // Create orthogonalization manager if we need to.
663  if (ortho_ == Teuchos::null || changedOrthoType) {
666  if (orthoType_=="DGKS" && orthoKappa_ > 0) {
667  paramsOrtho = Teuchos::rcp(new Teuchos::ParameterList());
668  paramsOrtho->set ("depTol", orthoKappa_ );
669  }
670 
671  ortho_ = factory.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, label_, paramsOrtho);
672  }
673 
674  // Check for convergence tolerance
675  if (params->isParameter("Convergence Tolerance")) {
676  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
677  convtol_ = params->get ("Convergence Tolerance",
678  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
679  }
680  else {
681  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
682  }
683 
684  // Update parameter in our list and residual tests.
685  params_->set("Convergence Tolerance", convtol_);
686  if (impConvTest_ != Teuchos::null)
687  impConvTest_->setTolerance( convtol_ );
688  if (expConvTest_ != Teuchos::null)
689  expConvTest_->setTolerance( convtol_ );
690  }
691 
692  // Check for a change in scaling, if so we need to build new residual tests.
693  if (params->isParameter("Implicit Residual Scaling")) {
694  std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
695 
696  // Only update the scaling if it's different.
697  if (impResScale_ != tempImpResScale) {
698  Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
699  impResScale_ = tempImpResScale;
700 
701  // Update parameter in our list and residual tests
702  params_->set("Implicit Residual Scaling", impResScale_);
703  if (impConvTest_ != Teuchos::null) {
704  try {
705  impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
706  }
707  catch (std::exception& e) {
708  // Make sure the convergence test gets constructed again.
709  isSTSet_ = false;
710  }
711  }
712  }
713  }
714 
715  if (params->isParameter("Explicit Residual Scaling")) {
716  std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
717 
718  // Only update the scaling if it's different.
719  if (expResScale_ != tempExpResScale) {
720  Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
721  expResScale_ = tempExpResScale;
722 
723  // Update parameter in our list and residual tests
724  params_->set("Explicit Residual Scaling", expResScale_);
725  if (expConvTest_ != Teuchos::null) {
726  try {
727  expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
728  }
729  catch (std::exception& e) {
730  // Make sure the convergence test gets constructed again.
731  isSTSet_ = false;
732  }
733  }
734  }
735  }
736 
737  if (params->isParameter("Explicit Residual Test")) {
738  expResTest_ = Teuchos::getParameter<bool>( *params,"Explicit Residual Test" );
739 
740  // Reconstruct the convergence test if the explicit residual test is not being used.
741  params_->set("Explicit Residual Test", expResTest_);
742  if (expConvTest_ == Teuchos::null) {
743  isSTSet_ = false;
744  }
745  }
746 
747  if (params->isParameter("Show Maximum Residual Norm Only")) {
748  showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
749 
750  // Update parameter in our list and residual tests
751  params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
752  if (impConvTest_ != Teuchos::null)
753  impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
754  if (expConvTest_ != Teuchos::null)
755  expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
756  }
757 
758 
759  // Create the timer if we need to.
760  if (timerSolve_ == Teuchos::null) {
761  std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
762 #ifdef BELOS_TEUCHOS_TIME_MONITOR
763  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
764 #endif
765  }
766 
767  // Inform the solver manager that the current parameters were set.
768  isSet_ = true;
769 }
770 
771 // Check the status test versus the defined linear problem
772 template<class ScalarType, class MV, class OP>
774 
775  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
776  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestGenResNorm_t;
777  typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP> StatusTestImpResNorm_t;
778 
779  // Basic test checks maximum iterations and native residual.
780  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
781 
782  // Perform sanity checking for flexible Gmres here.
783  // NOTE: If the user requests that the solver manager use flexible GMRES, but there is no right preconditioner, don't use flexible GMRES.
784  // Throw an error is the user provided a left preconditioner, as that is inconsistent with flexible GMRES.
785  if (isFlexible_ && Teuchos::is_null(problem_->getRightPrec())) {
786  isFlexible_ = false;
787  params_->set("Flexible Gmres", isFlexible_);
788 
789  // If the user specified the preconditioner as a left preconditioner, throw an error.
791  "Belos::BlockGmresSolMgr::solve(): Linear problem has a left preconditioner, not a right preconditioner, which is required for flexible GMRES.");
792  }
793 
794  // If there is a left preconditioner, we create a combined status test that checks the implicit
795  // and then explicit residual norm to see if we have convergence.
796  if (!Teuchos::is_null(problem_->getLeftPrec()) && !isFlexible_) {
797  expResTest_ = true;
798  }
799 
800  if (expResTest_) {
801 
802  // Implicit residual test, using the native residual to determine if convergence was achieved.
803  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
804  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
805  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
806  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
807  impConvTest_ = tmpImpConvTest;
808 
809  // Explicit residual test once the native residual is below the tolerance
810  Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
811  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
812  tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
813  tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
814  tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
815  expConvTest_ = tmpExpConvTest;
816 
817  // The convergence test is a combination of the "cheap" implicit test and explicit test.
818  convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
819  }
820  else {
821 
822  if (isFlexible_) {
823  // Implicit residual test, using the native residual to determine if convergence was achieved.
824  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
825  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
826  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
827  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
828  impConvTest_ = tmpImpConvTest;
829  }
830  else {
831  // Implicit residual test, using the native residual to determine if convergence was achieved.
832  // Use test that checks for loss of accuracy.
833  Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
834  Teuchos::rcp( new StatusTestImpResNorm_t( convtol_ ) );
835  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
836  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
837  impConvTest_ = tmpImpConvTest;
838  }
839 
840  // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
841  expConvTest_ = impConvTest_;
842  convTest_ = impConvTest_;
843  }
844 
845  // Create the status test.
846  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
847 
848  // Add debug status test, if one is provided by the user
849  if (nonnull(debugStatusTest_) ) {
850  // Add debug convergence test
851  Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(sTest_)->addStatusTest( debugStatusTest_ );
852  }
853 
854  // Create the status test output class.
855  // This class manages and formats the output from the status test.
856  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
857  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
858 
859  // Set the solver string for the output test
860  std::string solverDesc = " Block Gmres ";
861  if (isFlexible_)
862  solverDesc = "Flexible" + solverDesc;
863  outputTest_->setSolverDesc( solverDesc );
864 
865  // The status test is now set.
866  isSTSet_ = true;
867 
868  return false;
869 }
870 
871 template<class ScalarType, class MV, class OP>
873  const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest
874  )
875 {
876  debugStatusTest_ = debugStatusTest;
877 }
878 
879 
880 // solve()
881 template<class ScalarType, class MV, class OP>
883 
884  // Set the current parameters if they were not set before.
885  // NOTE: This may occur if the user generated the solver manager with the default constructor and
886  // then didn't set any parameters using setParameters().
887  if (!isSet_) {
888  setParameters(Teuchos::parameterList(*getValidParameters()));
889  }
890 
892  "Belos::BlockGmresSolMgr::solve(): Linear problem is not a valid object.");
893 
895  "Belos::BlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
896 
897  if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
899  "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
900  }
901 
902  // Create indices for the linear systems to be solved.
903  int startPtr = 0;
904  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
905  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
906 
907  std::vector<int> currIdx;
908  // If an adaptive block size is allowed then only the linear systems that need to be solved are solved.
909  // Otherwise, the index set is generated that informs the linear problem that some linear systems are augmented.
910  if ( adaptiveBlockSize_ ) {
911  blockSize_ = numCurrRHS;
912  currIdx.resize( numCurrRHS );
913  for (int i=0; i<numCurrRHS; ++i)
914  { currIdx[i] = startPtr+i; }
915 
916  }
917  else {
918  currIdx.resize( blockSize_ );
919  for (int i=0; i<numCurrRHS; ++i)
920  { currIdx[i] = startPtr+i; }
921  for (int i=numCurrRHS; i<blockSize_; ++i)
922  { currIdx[i] = -1; }
923  }
924 
925  // Inform the linear problem of the current linear system to solve.
926  problem_->setLSIndex( currIdx );
927 
929  // Parameter list
931  plist.set("Block Size",blockSize_);
932 
933  ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
934  if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {
935  int tmpNumBlocks = 0;
936  if (blockSize_ == 1)
937  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
938  else
939  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
940  printer_->stream(Warnings) <<
941  "Belos::BlockGmresSolMgr::solve(): Warning! Requested Krylov subspace dimension is larger than operator dimension!"
942  << std::endl << " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
943  plist.set("Num Blocks",tmpNumBlocks);
944  }
945  else
946  plist.set("Num Blocks",numBlocks_);
947 
948  // Reset the status test.
949  outputTest_->reset();
950  loaDetected_ = false;
951 
952  // Assume convergence is achieved, then let any failed convergence set this to false.
953  bool isConverged = true;
954 
956  // BlockGmres solver
957 
959 
960  if (isFlexible_)
961  block_gmres_iter = Teuchos::rcp( new BlockFGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
962  else
963  block_gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
964 
965  // Enter solve() iterations
966  {
967 #ifdef BELOS_TEUCHOS_TIME_MONITOR
968  Teuchos::TimeMonitor slvtimer(*timerSolve_);
969 #endif
970 
971  while ( numRHS2Solve > 0 ) {
972 
973  // Set the current number of blocks and blocksize with the Gmres iteration.
974  if (blockSize_*numBlocks_ > dim) {
975  int tmpNumBlocks = 0;
976  if (blockSize_ == 1)
977  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
978  else
979  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
980  block_gmres_iter->setSize( blockSize_, tmpNumBlocks );
981  }
982  else
983  block_gmres_iter->setSize( blockSize_, numBlocks_ );
984 
985  // Reset the number of iterations.
986  block_gmres_iter->resetNumIters();
987 
988  // Reset the number of calls that the status test output knows about.
989  outputTest_->resetNumCalls();
990 
991  // Create the first block in the current Krylov basis.
992  Teuchos::RCP<MV> V_0;
993  if (isFlexible_) {
994  // Load the correct residual if the system is augmented
995  if (currIdx[blockSize_-1] == -1) {
996  V_0 = MVT::Clone( *(problem_->getInitResVec()), blockSize_ );
997  problem_->computeCurrResVec( &*V_0 );
998  }
999  else {
1000  V_0 = MVT::CloneCopy( *(problem_->getInitResVec()), currIdx );
1001  }
1002  }
1003  else {
1004  // Load the correct residual if the system is augmented
1005  if (currIdx[blockSize_-1] == -1) {
1006  V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
1007  problem_->computeCurrPrecResVec( &*V_0 );
1008  }
1009  else {
1010  V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
1011  }
1012  }
1013 
1014  // Get a matrix to hold the orthonormalization coefficients.
1016  Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1017 
1018  // Orthonormalize the new V_0
1019  int rank = ortho_->normalize( *V_0, z_0 );
1021  "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
1022 
1023  // Set the new state and initialize the solver.
1025  newstate.V = V_0;
1026  newstate.z = z_0;
1027  newstate.curDim = 0;
1028  block_gmres_iter->initializeGmres(newstate);
1029  int numRestarts = 0;
1030 
1031  while(1) {
1032  // tell block_gmres_iter to iterate
1033  try {
1034  block_gmres_iter->iterate();
1035 
1037  //
1038  // check convergence first
1039  //
1041  if ( convTest_->getStatus() == Passed ) {
1042  if ( expConvTest_->getLOADetected() ) {
1043  // we don't have convergence
1044  loaDetected_ = true;
1045  printer_->stream(Warnings) <<
1046  "Belos::BlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1047  isConverged = false;
1048  }
1049  break; // break from while(1){block_gmres_iter->iterate()}
1050  }
1052  //
1053  // check for maximum iterations
1054  //
1056  else if ( maxIterTest_->getStatus() == Passed ) {
1057  // we don't have convergence
1058  isConverged = false;
1059  break; // break from while(1){block_gmres_iter->iterate()}
1060  }
1062  //
1063  // check for restarting, i.e. the subspace is full
1064  //
1066  else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
1067 
1068  if ( numRestarts >= maxRestarts_ ) {
1069  isConverged = false;
1070  break; // break from while(1){block_gmres_iter->iterate()}
1071  }
1072  numRestarts++;
1073 
1074  printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
1075 
1076  // Update the linear problem.
1077  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1078  if (isFlexible_) {
1079  // Update the solution manually, since the preconditioning doesn't need to be undone.
1080  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1081  MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1082  }
1083  else
1084  problem_->updateSolution( update, true );
1085 
1086  // Get the state.
1087  GmresIterationState<ScalarType,MV> oldState = block_gmres_iter->getState();
1088 
1089  // Compute the restart std::vector.
1090  // Get a view of the current Krylov basis.
1091  V_0 = MVT::Clone( *(oldState.V), blockSize_ );
1092  if (isFlexible_)
1093  problem_->computeCurrResVec( &*V_0 );
1094  else
1095  problem_->computeCurrPrecResVec( &*V_0 );
1096 
1097  // Get a view of the first block of the Krylov basis.
1098  z_0 = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1099 
1100  // Orthonormalize the new V_0
1101  rank = ortho_->normalize( *V_0, z_0 );
1103  "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
1104 
1105  // Set the new state and initialize the solver.
1106  newstate.V = V_0;
1107  newstate.z = z_0;
1108  newstate.curDim = 0;
1109  block_gmres_iter->initializeGmres(newstate);
1110 
1111  } // end of restarting
1112 
1114  //
1115  // we returned from iterate(), but none of our status tests Passed.
1116  // something is wrong, and it is probably our fault.
1117  //
1119 
1120  else {
1121  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
1122  "Belos::BlockGmresSolMgr::solve(): Invalid return from BlockGmresIter::iterate().");
1123  }
1124  }
1125  catch (const GmresIterationOrthoFailure &e) {
1126  // If the block size is not one, it's not considered a lucky breakdown.
1127  if (blockSize_ != 1) {
1128  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1129  << block_gmres_iter->getNumIters() << std::endl
1130  << e.what() << std::endl;
1131  if (convTest_->getStatus() != Passed)
1132  isConverged = false;
1133  break;
1134  }
1135  else {
1136  // If the block size is one, try to recover the most recent least-squares solution
1137  block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
1138 
1139  // Check to see if the most recent least-squares solution yielded convergence.
1140  sTest_->checkStatus( &*block_gmres_iter );
1141  if (convTest_->getStatus() != Passed)
1142  isConverged = false;
1143  break;
1144  }
1145  }
1146  catch (const StatusTestNaNError& e) {
1147  // A NaN was detected in the solver. Set the solution to zero and return unconverged.
1148  achievedTol_ = MT::one();
1149  Teuchos::RCP<MV> X = problem_->getLHS();
1150  MVT::MvInit( *X, SCT::zero() );
1151  printer_->stream(Warnings) << "Belos::BlockGmresSolMgr::solve(): Warning! NaN has been detected!"
1152  << std::endl;
1153  return Unconverged;
1154  }
1155  catch (const std::exception &e) {
1156  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1157  << block_gmres_iter->getNumIters() << std::endl
1158  << e.what() << std::endl;
1159  throw;
1160  }
1161  }
1162 
1163  // Compute the current solution.
1164  // Update the linear problem.
1165  if (isFlexible_) {
1166  // Update the solution manually, since the preconditioning doesn't need to be undone.
1167  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1168  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1169  // Update the solution only if there is a valid update from the iteration
1170  if (update != Teuchos::null)
1171  MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1172  }
1173  else {
1174  // Attempt to get the current solution from the residual status test, if it has one.
1175  if ( !Teuchos::is_null(expConvTest_->getSolution()) ) {
1176  Teuchos::RCP<MV> newX = expConvTest_->getSolution();
1177  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1178  MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
1179  }
1180  else {
1181  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1182  problem_->updateSolution( update, true );
1183  }
1184  }
1185 
1186  // Inform the linear problem that we are finished with this block linear system.
1187  problem_->setCurrLS();
1188 
1189  // Update indices for the linear systems to be solved.
1190  startPtr += numCurrRHS;
1191  numRHS2Solve -= numCurrRHS;
1192  if ( numRHS2Solve > 0 ) {
1193  numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1194 
1195  if ( adaptiveBlockSize_ ) {
1196  blockSize_ = numCurrRHS;
1197  currIdx.resize( numCurrRHS );
1198  for (int i=0; i<numCurrRHS; ++i)
1199  { currIdx[i] = startPtr+i; }
1200  }
1201  else {
1202  currIdx.resize( blockSize_ );
1203  for (int i=0; i<numCurrRHS; ++i)
1204  { currIdx[i] = startPtr+i; }
1205  for (int i=numCurrRHS; i<blockSize_; ++i)
1206  { currIdx[i] = -1; }
1207  }
1208  // Set the next indices.
1209  problem_->setLSIndex( currIdx );
1210  }
1211  else {
1212  currIdx.resize( numRHS2Solve );
1213  }
1214 
1215  }// while ( numRHS2Solve > 0 )
1216 
1217  }
1218 
1219  // print final summary
1220  sTest_->print( printer_->stream(FinalSummary) );
1221 
1222  // print timing information
1223 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1224  // Calling summarize() can be expensive, so don't call unless the
1225  // user wants to print out timing details. summarize() will do all
1226  // the work even if it's passed a "black hole" output stream.
1227  if (verbosity_ & TimingDetails)
1228  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1229 #endif
1230 
1231  // get iteration information for this solve
1232  numIters_ = maxIterTest_->getNumIters();
1233 
1234  // Save the convergence test value ("achieved tolerance") for this
1235  // solve. This requires a bit more work than for BlockCGSolMgr,
1236  // since for this solver, convTest_ may either be a single residual
1237  // norm test, or a combination of two residual norm tests. In the
1238  // latter case, the master convergence test convTest_ is a SEQ combo
1239  // of the implicit resp. explicit tests. If the implicit test never
1240  // passes, then the explicit test won't ever be executed. This
1241  // manifests as expConvTest_->getTestValue()->size() < 1. We deal
1242  // with this case by using the values returned by
1243  // impConvTest_->getTestValue().
1244  {
1245  // We'll fetch the vector of residual norms one way or the other.
1246  const std::vector<MagnitudeType>* pTestValues = NULL;
1247  if (expResTest_) {
1248  pTestValues = expConvTest_->getTestValue();
1249  if (pTestValues == NULL || pTestValues->size() < 1) {
1250  pTestValues = impConvTest_->getTestValue();
1251  }
1252  }
1253  else {
1254  // Only the implicit residual norm test is being used.
1255  pTestValues = impConvTest_->getTestValue();
1256  }
1257  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1258  "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1259  "getTestValue() method returned NULL. Please report this bug to the "
1260  "Belos developers.");
1261  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1262  "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1263  "getTestValue() method returned a vector of length zero. Please report "
1264  "this bug to the Belos developers.");
1265 
1266  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1267  // achieved tolerances for all vectors in the current solve(), or
1268  // just for the vectors from the last deflation?
1269  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1270  }
1271 
1272  if (!isConverged || loaDetected_) {
1273  return Unconverged; // return from BlockGmresSolMgr::solve()
1274  }
1275  return Converged; // return from BlockGmresSolMgr::solve()
1276 }
1277 
1278 
1279 template<class ScalarType, class MV, class OP>
1281 {
1282  std::ostringstream out;
1283  out << "\"Belos::BlockGmresSolMgr\": {";
1284  if (this->getObjectLabel () != "") {
1285  out << "Label: " << this->getObjectLabel () << ", ";
1286  }
1287  out << "Flexible: " << (isFlexible_ ? "true" : "false")
1288  << ", Num Blocks: " << numBlocks_
1289  << ", Maximum Iterations: " << maxIters_
1290  << ", Maximum Restarts: " << maxRestarts_
1291  << ", Convergence Tolerance: " << convtol_
1292  << "}";
1293  return out.str ();
1294 }
1295 
1296 
1297 template<class ScalarType, class MV, class OP>
1298 void
1301  const Teuchos::EVerbosityLevel verbLevel) const
1302 {
1304  using Teuchos::VERB_DEFAULT;
1305  using Teuchos::VERB_NONE;
1306  using Teuchos::VERB_LOW;
1307  // using Teuchos::VERB_MEDIUM;
1308  // using Teuchos::VERB_HIGH;
1309  // using Teuchos::VERB_EXTREME;
1310  using std::endl;
1311 
1312  // Set default verbosity if applicable.
1313  const Teuchos::EVerbosityLevel vl =
1314  (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
1315 
1316  if (vl != VERB_NONE) {
1317  Teuchos::OSTab tab0 (out);
1318 
1319  out << "\"Belos::BlockGmresSolMgr\":" << endl;
1320  Teuchos::OSTab tab1 (out);
1321  out << "Template parameters:" << endl;
1322  {
1323  Teuchos::OSTab tab2 (out);
1324  out << "ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
1325  << "MV: " << TypeNameTraits<MV>::name () << endl
1326  << "OP: " << TypeNameTraits<OP>::name () << endl;
1327  }
1328  if (this->getObjectLabel () != "") {
1329  out << "Label: " << this->getObjectLabel () << endl;
1330  }
1331  out << "Flexible: " << (isFlexible_ ? "true" : "false") << endl
1332  << "Num Blocks: " << numBlocks_ << endl
1333  << "Maximum Iterations: " << maxIters_ << endl
1334  << "Maximum Restarts: " << maxRestarts_ << endl
1335  << "Convergence Tolerance: " << convtol_ << endl;
1336  }
1337 }
1338 
1339 
1340 } // end Belos namespace
1341 
1342 #endif /* BELOS_BLOCK_GMRES_SOLMGR_HPP */
static const double orthoKappa
DGKS orthogonalization constant.
Definition: BelosTypes.hpp:299
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:106
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:120
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.
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
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.
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:293
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:206
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:155
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:60
BlockGmresSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate ortho...
Default parameters common to most Belos solvers.
Definition: BelosTypes.hpp:283
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.