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

Generated on Fri Nov 22 2024 09:26:14 for Belos by doxygen 1.8.5