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

Generated on Mon Feb 5 2018 15:01:51 for Belos by doxygen 1.8.13