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