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 #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:
128  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
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.
303  Teuchos::RCP<std::ostream> outputStream_;
304 
305  // Status test.
310  Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
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.
337  MagnitudeType convtol_, orthoKappa_, achievedTol_;
338  int maxRestarts_, maxIters_, numIters_;
339  int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_;
340  bool adaptiveBlockSize_, showMaxResNormOnly_, isFlexible_, expResTest_;
341  std::string orthoType_;
342  std::string impResScale_, expResScale_;
343 
344  // Timers.
345  std::string label_;
346  Teuchos::RCP<Teuchos::Time> timerSolve_;
347 
348  // Internal state variables.
349  bool isSet_, isSTSet_;
350  bool loaDetected_;
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) {
665  Teuchos::RCP<Teuchos::ParameterList> paramsOrtho; // can be null
666  if (orthoType_=="DGKS" && orthoKappa_ > 0) {
667  paramsOrtho->set ("depTol", orthoKappa_ );
668  }
669 
670  ortho_ = factory.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, label_, paramsOrtho);
671  }
672 
673  // Check for convergence tolerance
674  if (params->isParameter("Convergence Tolerance")) {
675  if (params->isType<MagnitudeType> ("Convergence Tolerance")) {
676  convtol_ = params->get ("Convergence Tolerance",
677  static_cast<MagnitudeType> (DefaultSolverParameters::convTol));
678  }
679  else {
680  convtol_ = params->get ("Convergence Tolerance", DefaultSolverParameters::convTol);
681  }
682 
683  // Update parameter in our list and residual tests.
684  params_->set("Convergence Tolerance", convtol_);
685  if (impConvTest_ != Teuchos::null)
686  impConvTest_->setTolerance( convtol_ );
687  if (expConvTest_ != Teuchos::null)
688  expConvTest_->setTolerance( convtol_ );
689  }
690 
691  // Check for a change in scaling, if so we need to build new residual tests.
692  if (params->isParameter("Implicit Residual Scaling")) {
693  std::string tempImpResScale = Teuchos::getParameter<std::string>( *params, "Implicit Residual Scaling" );
694 
695  // Only update the scaling if it's different.
696  if (impResScale_ != tempImpResScale) {
697  Belos::ScaleType impResScaleType = convertStringToScaleType( tempImpResScale );
698  impResScale_ = tempImpResScale;
699 
700  // Update parameter in our list and residual tests
701  params_->set("Implicit Residual Scaling", impResScale_);
702  if (impConvTest_ != Teuchos::null) {
703  try {
704  impConvTest_->defineScaleForm( impResScaleType, Belos::TwoNorm );
705  }
706  catch (std::exception& e) {
707  // Make sure the convergence test gets constructed again.
708  isSTSet_ = false;
709  }
710  }
711  }
712  }
713 
714  if (params->isParameter("Explicit Residual Scaling")) {
715  std::string tempExpResScale = Teuchos::getParameter<std::string>( *params, "Explicit Residual Scaling" );
716 
717  // Only update the scaling if it's different.
718  if (expResScale_ != tempExpResScale) {
719  Belos::ScaleType expResScaleType = convertStringToScaleType( tempExpResScale );
720  expResScale_ = tempExpResScale;
721 
722  // Update parameter in our list and residual tests
723  params_->set("Explicit Residual Scaling", expResScale_);
724  if (expConvTest_ != Teuchos::null) {
725  try {
726  expConvTest_->defineScaleForm( expResScaleType, Belos::TwoNorm );
727  }
728  catch (std::exception& e) {
729  // Make sure the convergence test gets constructed again.
730  isSTSet_ = false;
731  }
732  }
733  }
734  }
735 
736  if (params->isParameter("Explicit Residual Test")) {
737  expResTest_ = Teuchos::getParameter<bool>( *params,"Explicit Residual Test" );
738 
739  // Reconstruct the convergence test if the explicit residual test is not being used.
740  params_->set("Explicit Residual Test", expResTest_);
741  if (expConvTest_ == Teuchos::null) {
742  isSTSet_ = false;
743  }
744  }
745 
746  if (params->isParameter("Show Maximum Residual Norm Only")) {
747  showMaxResNormOnly_ = Teuchos::getParameter<bool>(*params,"Show Maximum Residual Norm Only");
748 
749  // Update parameter in our list and residual tests
750  params_->set("Show Maximum Residual Norm Only", showMaxResNormOnly_);
751  if (impConvTest_ != Teuchos::null)
752  impConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
753  if (expConvTest_ != Teuchos::null)
754  expConvTest_->setShowMaxResNormOnly( showMaxResNormOnly_ );
755  }
756 
757 
758  // Create the timer if we need to.
759  if (timerSolve_ == Teuchos::null) {
760  std::string solveLabel = label_ + ": BlockGmresSolMgr total solve time";
761 #ifdef BELOS_TEUCHOS_TIME_MONITOR
762  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
763 #endif
764  }
765 
766  // Inform the solver manager that the current parameters were set.
767  isSet_ = true;
768 }
769 
770 // Check the status test versus the defined linear problem
771 template<class ScalarType, class MV, class OP>
773 
774  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
775  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestGenResNorm_t;
776  typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP> StatusTestImpResNorm_t;
777 
778  // Basic test checks maximum iterations and native residual.
779  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
780 
781  // Perform sanity checking for flexible Gmres here.
782  // NOTE: If the user requests that the solver manager use flexible GMRES, but there is no right preconditioner, don't use flexible GMRES.
783  // Throw an error is the user provided a left preconditioner, as that is inconsistent with flexible GMRES.
784  if (isFlexible_ && Teuchos::is_null(problem_->getRightPrec())) {
785  isFlexible_ = false;
786  params_->set("Flexible Gmres", isFlexible_);
787 
788  // If the user specified the preconditioner as a left preconditioner, throw an error.
790  "Belos::BlockGmresSolMgr::solve(): Linear problem has a left preconditioner, not a right preconditioner, which is required for flexible GMRES.");
791  }
792 
793  // If there is a left preconditioner, we create a combined status test that checks the implicit
794  // and then explicit residual norm to see if we have convergence.
795  if (!Teuchos::is_null(problem_->getLeftPrec()) && !isFlexible_) {
796  expResTest_ = true;
797  }
798 
799  if (expResTest_) {
800 
801  // Implicit residual test, using the native residual to determine if convergence was achieved.
802  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
803  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
804  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
805  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
806  impConvTest_ = tmpImpConvTest;
807 
808  // Explicit residual test once the native residual is below the tolerance
809  Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
810  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
811  tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
812  tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
813  tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
814  expConvTest_ = tmpExpConvTest;
815 
816  // The convergence test is a combination of the "cheap" implicit test and explicit test.
817  convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
818  }
819  else {
820 
821  if (isFlexible_) {
822  // Implicit residual test, using the native residual to determine if convergence was achieved.
823  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
824  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_ ) );
825  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
826  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
827  impConvTest_ = tmpImpConvTest;
828  }
829  else {
830  // Implicit residual test, using the native residual to determine if convergence was achieved.
831  // Use test that checks for loss of accuracy.
832  Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
833  Teuchos::rcp( new StatusTestImpResNorm_t( convtol_ ) );
834  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
835  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
836  impConvTest_ = tmpImpConvTest;
837  }
838 
839  // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
840  expConvTest_ = impConvTest_;
841  convTest_ = impConvTest_;
842  }
843 
844  // Create the status test.
845  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
846 
847  // Add debug status test, if one is provided by the user
848  if (nonnull(debugStatusTest_) ) {
849  // Add debug convergence test
850  Teuchos::rcp_dynamic_cast<StatusTestCombo_t>(sTest_)->addStatusTest( debugStatusTest_ );
851  }
852 
853  // Create the status test output class.
854  // This class manages and formats the output from the status test.
855  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
856  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
857 
858  // Set the solver string for the output test
859  std::string solverDesc = " Block Gmres ";
860  if (isFlexible_)
861  solverDesc = "Flexible" + solverDesc;
862  outputTest_->setSolverDesc( solverDesc );
863 
864  // The status test is now set.
865  isSTSet_ = true;
866 
867  return false;
868 }
869 
870 template<class ScalarType, class MV, class OP>
872  const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest
873  )
874 {
875  debugStatusTest_ = debugStatusTest;
876 }
877 
878 
879 // solve()
880 template<class ScalarType, class MV, class OP>
882 
883  // Set the current parameters if they were not set before.
884  // NOTE: This may occur if the user generated the solver manager with the default constructor and
885  // then didn't set any parameters using setParameters().
886  if (!isSet_) {
887  setParameters(Teuchos::parameterList(*getValidParameters()));
888  }
889 
891  "Belos::BlockGmresSolMgr::solve(): Linear problem is not a valid object.");
892 
894  "Belos::BlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
895 
896  if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
898  "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
899  }
900 
901  // Create indices for the linear systems to be solved.
902  int startPtr = 0;
903  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
904  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
905 
906  std::vector<int> currIdx;
907  // If an adaptive block size is allowed then only the linear systems that need to be solved are solved.
908  // Otherwise, the index set is generated that informs the linear problem that some linear systems are augmented.
909  if ( adaptiveBlockSize_ ) {
910  blockSize_ = numCurrRHS;
911  currIdx.resize( numCurrRHS );
912  for (int i=0; i<numCurrRHS; ++i)
913  { currIdx[i] = startPtr+i; }
914 
915  }
916  else {
917  currIdx.resize( blockSize_ );
918  for (int i=0; i<numCurrRHS; ++i)
919  { currIdx[i] = startPtr+i; }
920  for (int i=numCurrRHS; i<blockSize_; ++i)
921  { currIdx[i] = -1; }
922  }
923 
924  // Inform the linear problem of the current linear system to solve.
925  problem_->setLSIndex( currIdx );
926 
928  // Parameter list
930  plist.set("Block Size",blockSize_);
931 
932  ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
933  if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {
934  int tmpNumBlocks = 0;
935  if (blockSize_ == 1)
936  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
937  else
938  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
939  printer_->stream(Warnings) <<
940  "Belos::BlockGmresSolMgr::solve(): Warning! Requested Krylov subspace dimension is larger than operator dimension!"
941  << std::endl << " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
942  plist.set("Num Blocks",tmpNumBlocks);
943  }
944  else
945  plist.set("Num Blocks",numBlocks_);
946 
947  // Reset the status test.
948  outputTest_->reset();
949  loaDetected_ = false;
950 
951  // Assume convergence is achieved, then let any failed convergence set this to false.
952  bool isConverged = true;
953 
955  // BlockGmres solver
956 
958 
959  if (isFlexible_)
960  block_gmres_iter = Teuchos::rcp( new BlockFGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
961  else
962  block_gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
963 
964  // Enter solve() iterations
965  {
966 #ifdef BELOS_TEUCHOS_TIME_MONITOR
967  Teuchos::TimeMonitor slvtimer(*timerSolve_);
968 #endif
969 
970  while ( numRHS2Solve > 0 ) {
971 
972  // Set the current number of blocks and blocksize with the Gmres iteration.
973  if (blockSize_*numBlocks_ > dim) {
974  int tmpNumBlocks = 0;
975  if (blockSize_ == 1)
976  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
977  else
978  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
979  block_gmres_iter->setSize( blockSize_, tmpNumBlocks );
980  }
981  else
982  block_gmres_iter->setSize( blockSize_, numBlocks_ );
983 
984  // Reset the number of iterations.
985  block_gmres_iter->resetNumIters();
986 
987  // Reset the number of calls that the status test output knows about.
988  outputTest_->resetNumCalls();
989 
990  // Create the first block in the current Krylov basis.
991  Teuchos::RCP<MV> V_0;
992  if (isFlexible_) {
993  // Load the correct residual if the system is augmented
994  if (currIdx[blockSize_-1] == -1) {
995  V_0 = MVT::Clone( *(problem_->getInitResVec()), blockSize_ );
996  problem_->computeCurrResVec( &*V_0 );
997  }
998  else {
999  V_0 = MVT::CloneCopy( *(problem_->getInitResVec()), currIdx );
1000  }
1001  }
1002  else {
1003  // Load the correct residual if the system is augmented
1004  if (currIdx[blockSize_-1] == -1) {
1005  V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
1006  problem_->computeCurrPrecResVec( &*V_0 );
1007  }
1008  else {
1009  V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
1010  }
1011  }
1012 
1013  // Get a matrix to hold the orthonormalization coefficients.
1015  Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1016 
1017  // Orthonormalize the new V_0
1018  int rank = ortho_->normalize( *V_0, z_0 );
1020  "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
1021 
1022  // Set the new state and initialize the solver.
1024  newstate.V = V_0;
1025  newstate.z = z_0;
1026  newstate.curDim = 0;
1027  block_gmres_iter->initializeGmres(newstate);
1028  int numRestarts = 0;
1029 
1030  while(1) {
1031  // tell block_gmres_iter to iterate
1032  try {
1033  block_gmres_iter->iterate();
1034 
1036  //
1037  // check convergence first
1038  //
1040  if ( convTest_->getStatus() == Passed ) {
1041  if ( expConvTest_->getLOADetected() ) {
1042  // we don't have convergence
1043  loaDetected_ = true;
1044  printer_->stream(Warnings) <<
1045  "Belos::BlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1046  isConverged = false;
1047  }
1048  break; // break from while(1){block_gmres_iter->iterate()}
1049  }
1051  //
1052  // check for maximum iterations
1053  //
1055  else if ( maxIterTest_->getStatus() == Passed ) {
1056  // we don't have convergence
1057  isConverged = false;
1058  break; // break from while(1){block_gmres_iter->iterate()}
1059  }
1061  //
1062  // check for restarting, i.e. the subspace is full
1063  //
1065  else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
1066 
1067  if ( numRestarts >= maxRestarts_ ) {
1068  isConverged = false;
1069  break; // break from while(1){block_gmres_iter->iterate()}
1070  }
1071  numRestarts++;
1072 
1073  printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
1074 
1075  // Update the linear problem.
1076  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1077  if (isFlexible_) {
1078  // Update the solution manually, since the preconditioning doesn't need to be undone.
1079  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1080  MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1081  }
1082  else
1083  problem_->updateSolution( update, true );
1084 
1085  // Get the state.
1086  GmresIterationState<ScalarType,MV> oldState = block_gmres_iter->getState();
1087 
1088  // Compute the restart std::vector.
1089  // Get a view of the current Krylov basis.
1090  V_0 = MVT::Clone( *(oldState.V), blockSize_ );
1091  if (isFlexible_)
1092  problem_->computeCurrResVec( &*V_0 );
1093  else
1094  problem_->computeCurrPrecResVec( &*V_0 );
1095 
1096  // Get a view of the first block of the Krylov basis.
1097  z_0 = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( blockSize_, blockSize_ ) );
1098 
1099  // Orthonormalize the new V_0
1100  rank = ortho_->normalize( *V_0, z_0 );
1102  "Belos::BlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after restart.");
1103 
1104  // Set the new state and initialize the solver.
1105  newstate.V = V_0;
1106  newstate.z = z_0;
1107  newstate.curDim = 0;
1108  block_gmres_iter->initializeGmres(newstate);
1109 
1110  } // end of restarting
1111 
1113  //
1114  // we returned from iterate(), but none of our status tests Passed.
1115  // something is wrong, and it is probably our fault.
1116  //
1118 
1119  else {
1120  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
1121  "Belos::BlockGmresSolMgr::solve(): Invalid return from BlockGmresIter::iterate().");
1122  }
1123  }
1124  catch (const GmresIterationOrthoFailure &e) {
1125  // If the block size is not one, it's not considered a lucky breakdown.
1126  if (blockSize_ != 1) {
1127  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1128  << block_gmres_iter->getNumIters() << std::endl
1129  << e.what() << std::endl;
1130  if (convTest_->getStatus() != Passed)
1131  isConverged = false;
1132  break;
1133  }
1134  else {
1135  // If the block size is one, try to recover the most recent least-squares solution
1136  block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
1137 
1138  // Check to see if the most recent least-squares solution yielded convergence.
1139  sTest_->checkStatus( &*block_gmres_iter );
1140  if (convTest_->getStatus() != Passed)
1141  isConverged = false;
1142  break;
1143  }
1144  }
1145  catch (const StatusTestNaNError& e) {
1146  // A NaN was detected in the solver. Set the solution to zero and return unconverged.
1147  achievedTol_ = MT::one();
1148  Teuchos::RCP<MV> X = problem_->getLHS();
1149  MVT::MvInit( *X, SCT::zero() );
1150  printer_->stream(Warnings) << "Belos::BlockGmresSolMgr::solve(): Warning! NaN has been detected!"
1151  << std::endl;
1152  return Unconverged;
1153  }
1154  catch (const std::exception &e) {
1155  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
1156  << block_gmres_iter->getNumIters() << std::endl
1157  << e.what() << std::endl;
1158  throw;
1159  }
1160  }
1161 
1162  // Compute the current solution.
1163  // Update the linear problem.
1164  if (isFlexible_) {
1165  // Update the solution manually, since the preconditioning doesn't need to be undone.
1166  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1167  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1168  // Update the solution only if there is a valid update from the iteration
1169  if (update != Teuchos::null)
1170  MVT::MvAddMv( 1.0, *curX, 1.0, *update, *curX );
1171  }
1172  else {
1173  // Attempt to get the current solution from the residual status test, if it has one.
1174  if ( !Teuchos::is_null(expConvTest_->getSolution()) ) {
1175  Teuchos::RCP<MV> newX = expConvTest_->getSolution();
1176  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1177  MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
1178  }
1179  else {
1180  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1181  problem_->updateSolution( update, true );
1182  }
1183  }
1184 
1185  // Inform the linear problem that we are finished with this block linear system.
1186  problem_->setCurrLS();
1187 
1188  // Update indices for the linear systems to be solved.
1189  startPtr += numCurrRHS;
1190  numRHS2Solve -= numCurrRHS;
1191  if ( numRHS2Solve > 0 ) {
1192  numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1193 
1194  if ( adaptiveBlockSize_ ) {
1195  blockSize_ = numCurrRHS;
1196  currIdx.resize( numCurrRHS );
1197  for (int i=0; i<numCurrRHS; ++i)
1198  { currIdx[i] = startPtr+i; }
1199  }
1200  else {
1201  currIdx.resize( blockSize_ );
1202  for (int i=0; i<numCurrRHS; ++i)
1203  { currIdx[i] = startPtr+i; }
1204  for (int i=numCurrRHS; i<blockSize_; ++i)
1205  { currIdx[i] = -1; }
1206  }
1207  // Set the next indices.
1208  problem_->setLSIndex( currIdx );
1209  }
1210  else {
1211  currIdx.resize( numRHS2Solve );
1212  }
1213 
1214  }// while ( numRHS2Solve > 0 )
1215 
1216  }
1217 
1218  // print final summary
1219  sTest_->print( printer_->stream(FinalSummary) );
1220 
1221  // print timing information
1222 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1223  // Calling summarize() can be expensive, so don't call unless the
1224  // user wants to print out timing details. summarize() will do all
1225  // the work even if it's passed a "black hole" output stream.
1226  if (verbosity_ & TimingDetails)
1227  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1228 #endif
1229 
1230  // get iteration information for this solve
1231  numIters_ = maxIterTest_->getNumIters();
1232 
1233  // Save the convergence test value ("achieved tolerance") for this
1234  // solve. This requires a bit more work than for BlockCGSolMgr,
1235  // since for this solver, convTest_ may either be a single residual
1236  // norm test, or a combination of two residual norm tests. In the
1237  // latter case, the master convergence test convTest_ is a SEQ combo
1238  // of the implicit resp. explicit tests. If the implicit test never
1239  // passes, then the explicit test won't ever be executed. This
1240  // manifests as expConvTest_->getTestValue()->size() < 1. We deal
1241  // with this case by using the values returned by
1242  // impConvTest_->getTestValue().
1243  {
1244  // We'll fetch the vector of residual norms one way or the other.
1245  const std::vector<MagnitudeType>* pTestValues = NULL;
1246  if (expResTest_) {
1247  pTestValues = expConvTest_->getTestValue();
1248  if (pTestValues == NULL || pTestValues->size() < 1) {
1249  pTestValues = impConvTest_->getTestValue();
1250  }
1251  }
1252  else {
1253  // Only the implicit residual norm test is being used.
1254  pTestValues = impConvTest_->getTestValue();
1255  }
1256  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1257  "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1258  "getTestValue() method returned NULL. Please report this bug to the "
1259  "Belos developers.");
1260  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1261  "Belos::BlockGmresSolMgr::solve(): The implicit convergence test's "
1262  "getTestValue() method returned a vector of length zero. Please report "
1263  "this bug to the Belos developers.");
1264 
1265  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1266  // achieved tolerances for all vectors in the current solve(), or
1267  // just for the vectors from the last deflation?
1268  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1269  }
1270 
1271  if (!isConverged || loaDetected_) {
1272  return Unconverged; // return from BlockGmresSolMgr::solve()
1273  }
1274  return Converged; // return from BlockGmresSolMgr::solve()
1275 }
1276 
1277 
1278 template<class ScalarType, class MV, class OP>
1280 {
1281  std::ostringstream out;
1282  out << "\"Belos::BlockGmresSolMgr\": {";
1283  if (this->getObjectLabel () != "") {
1284  out << "Label: " << this->getObjectLabel () << ", ";
1285  }
1286  out << "Flexible: " << (isFlexible_ ? "true" : "false")
1287  << ", Num Blocks: " << numBlocks_
1288  << ", Maximum Iterations: " << maxIters_
1289  << ", Maximum Restarts: " << maxRestarts_
1290  << ", Convergence Tolerance: " << convtol_
1291  << "}";
1292  return out.str ();
1293 }
1294 
1295 
1296 template<class ScalarType, class MV, class OP>
1297 void
1300  const Teuchos::EVerbosityLevel verbLevel) const
1301 {
1303  using Teuchos::VERB_DEFAULT;
1304  using Teuchos::VERB_NONE;
1305  using Teuchos::VERB_LOW;
1306  // using Teuchos::VERB_MEDIUM;
1307  // using Teuchos::VERB_HIGH;
1308  // using Teuchos::VERB_EXTREME;
1309  using std::endl;
1310 
1311  // Set default verbosity if applicable.
1312  const Teuchos::EVerbosityLevel vl =
1313  (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
1314 
1315  if (vl != VERB_NONE) {
1316  Teuchos::OSTab tab0 (out);
1317 
1318  out << "\"Belos::BlockGmresSolMgr\":" << endl;
1319  Teuchos::OSTab tab1 (out);
1320  out << "Template parameters:" << endl;
1321  {
1322  Teuchos::OSTab tab2 (out);
1323  out << "ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
1324  << "MV: " << TypeNameTraits<MV>::name () << endl
1325  << "OP: " << TypeNameTraits<OP>::name () << endl;
1326  }
1327  if (this->getObjectLabel () != "") {
1328  out << "Label: " << this->getObjectLabel () << endl;
1329  }
1330  out << "Flexible: " << (isFlexible_ ? "true" : "false") << endl
1331  << "Num Blocks: " << numBlocks_ << endl
1332  << "Maximum Iterations: " << maxIters_ << endl
1333  << "Maximum Restarts: " << maxRestarts_ << endl
1334  << "Convergence Tolerance: " << convtol_ << endl;
1335  }
1336 }
1337 
1338 
1339 } // end Belos namespace
1340 
1341 #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 Apr 19 2024 09:25:11 for Belos by doxygen 1.8.5