Belos  Version of the Day
BelosPseudoBlockGmresSolMgr.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_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
43 #define BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP
44 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 
52 #include "BelosLinearProblem.hpp"
53 #include "BelosSolverManager.hpp"
54 
59 #ifdef HAVE_BELOS_TSQR
60 # include "BelosTsqrOrthoManager.hpp"
61 #endif // HAVE_BELOS_TSQR
65 #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 
80 namespace Belos {
81 
83 
84 
92  PseudoBlockGmresSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
93  {}};
94 
102  PseudoBlockGmresSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
103  {}};
104 
128  template<class ScalarType, class MV, class OP>
129  class PseudoBlockGmresSolMgr : public SolverManager<ScalarType,MV,OP> {
130 
131  private:
134  typedef Teuchos::ScalarTraits<ScalarType> SCT;
135  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
136  typedef Teuchos::ScalarTraits<MagnitudeType> MT;
137 
138  public:
139 
141 
142 
151 
255  PseudoBlockGmresSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
256  const Teuchos::RCP<Teuchos::ParameterList> &pl );
257 
261 
263 
264 
266  return *problem_;
267  }
268 
270  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
271 
273  Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
274 
280  Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
281  return Teuchos::tuple(timerSolve_);
282  }
283 
294  MagnitudeType achievedTol() const {
295  return achievedTol_;
296  }
297 
299  int getNumIters() const {
300  return numIters_;
301  }
302 
358  bool isLOADetected() const { return loaDetected_; }
359 
361 
363 
364 
366  void setProblem (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem) {
367  problem_ = problem;
368  }
369 
371  void setParameters (const Teuchos::RCP<Teuchos::ParameterList> &params);
372 
379  virtual void setUserConvStatusTest(
380  const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &userConvStatusTest
381  );
382 
388  virtual void setDebugStatusTest(
389  const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest
390  );
391 
393 
395 
396 
400  void reset( const ResetType type ) { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
402 
404 
405 
423  ReturnType solve();
424 
426 
429 
436  void
437  describe (Teuchos::FancyOStream& out,
438  const Teuchos::EVerbosityLevel verbLevel =
439  Teuchos::Describable::verbLevel_default) const;
440 
442  std::string description () const;
443 
445 
446  private:
447 
462  bool checkStatusTest();
463 
465  Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
466 
467  // Output manager.
468  Teuchos::RCP<OutputManager<ScalarType> > printer_;
469  Teuchos::RCP<std::ostream> outputStream_;
470 
471  // Status tests.
472  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > userConvStatusTest_;
473  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > debugStatusTest_;
474  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
475  Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
476  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
477  Teuchos::RCP<StatusTestResNorm<ScalarType,MV,OP> > impConvTest_, expConvTest_;
478  Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
479 
480  // Orthogonalization manager.
481  Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
482 
483  // Current parameter list.
484  Teuchos::RCP<Teuchos::ParameterList> params_;
485 
486  // Default solver values.
487  static const MagnitudeType convtol_default_;
488  static const MagnitudeType orthoKappa_default_;
489  static const int maxRestarts_default_;
490  static const int maxIters_default_;
491  static const bool showMaxResNormOnly_default_;
492  static const int blockSize_default_;
493  static const int numBlocks_default_;
494  static const int verbosity_default_;
495  static const int outputStyle_default_;
496  static const int outputFreq_default_;
497  static const int defQuorum_default_;
498  static const std::string impResScale_default_;
499  static const std::string expResScale_default_;
500  static const std::string label_default_;
501  static const std::string orthoType_default_;
502  static const Teuchos::RCP<std::ostream> outputStream_default_;
503 
504  // Current solver values.
505  MagnitudeType convtol_, orthoKappa_, achievedTol_;
506  int maxRestarts_, maxIters_, numIters_;
507  int blockSize_, numBlocks_, verbosity_, outputStyle_, outputFreq_, defQuorum_;
508  bool showMaxResNormOnly_;
509  std::string orthoType_;
510  std::string impResScale_, expResScale_;
511 
512  // Timers.
513  std::string label_;
514  Teuchos::RCP<Teuchos::Time> timerSolve_;
515 
516  // Internal state variables.
517  bool isSet_, isSTSet_, expResTest_;
518  bool loaDetected_;
519  };
520 
521 
522 // Default solver values.
523 template<class ScalarType, class MV, class OP>
524 const typename PseudoBlockGmresSolMgr<ScalarType,MV,OP>::MagnitudeType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::convtol_default_ = 1e-8;
525 
526 template<class ScalarType, class MV, class OP>
527 const typename PseudoBlockGmresSolMgr<ScalarType,MV,OP>::MagnitudeType PseudoBlockGmresSolMgr<ScalarType,MV,OP>::orthoKappa_default_ = -1.0;
528 
529 template<class ScalarType, class MV, class OP>
531 
532 template<class ScalarType, class MV, class OP>
534 
535 template<class ScalarType, class MV, class OP>
537 
538 template<class ScalarType, class MV, class OP>
540 
541 template<class ScalarType, class MV, class OP>
543 
544 template<class ScalarType, class MV, class OP>
546 
547 template<class ScalarType, class MV, class OP>
549 
550 template<class ScalarType, class MV, class OP>
552 
553 template<class ScalarType, class MV, class OP>
555 
556 template<class ScalarType, class MV, class OP>
557 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::impResScale_default_ = "Norm of Preconditioned Initial Residual";
558 
559 template<class ScalarType, class MV, class OP>
560 const std::string PseudoBlockGmresSolMgr<ScalarType,MV,OP>::expResScale_default_ = "Norm of Initial Residual";
561 
562 template<class ScalarType, class MV, class OP>
564 
565 template<class ScalarType, class MV, class OP>
567 
568 template<class ScalarType, class MV, class OP>
569 const Teuchos::RCP<std::ostream> PseudoBlockGmresSolMgr<ScalarType,MV,OP>::outputStream_default_ = Teuchos::rcp(&std::cout,false);
570 
571 
572 // Empty Constructor
573 template<class ScalarType, class MV, class OP>
575  outputStream_(outputStream_default_),
576  convtol_(convtol_default_),
577  orthoKappa_(orthoKappa_default_),
578  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
579  maxRestarts_(maxRestarts_default_),
580  maxIters_(maxIters_default_),
581  numIters_(0),
582  blockSize_(blockSize_default_),
583  numBlocks_(numBlocks_default_),
584  verbosity_(verbosity_default_),
585  outputStyle_(outputStyle_default_),
586  outputFreq_(outputFreq_default_),
587  defQuorum_(defQuorum_default_),
588  showMaxResNormOnly_(showMaxResNormOnly_default_),
589  orthoType_(orthoType_default_),
590  impResScale_(impResScale_default_),
591  expResScale_(expResScale_default_),
592  label_(label_default_),
593  isSet_(false),
594  isSTSet_(false),
595  expResTest_(false),
596  loaDetected_(false)
597 {}
598 
599 // Basic Constructor
600 template<class ScalarType, class MV, class OP>
603  const Teuchos::RCP<Teuchos::ParameterList> &pl) :
604  problem_(problem),
605  outputStream_(outputStream_default_),
606  convtol_(convtol_default_),
607  orthoKappa_(orthoKappa_default_),
608  achievedTol_(Teuchos::ScalarTraits<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>::zero()),
609  maxRestarts_(maxRestarts_default_),
610  maxIters_(maxIters_default_),
611  numIters_(0),
612  blockSize_(blockSize_default_),
613  numBlocks_(numBlocks_default_),
614  verbosity_(verbosity_default_),
615  outputStyle_(outputStyle_default_),
616  outputFreq_(outputFreq_default_),
617  defQuorum_(defQuorum_default_),
618  showMaxResNormOnly_(showMaxResNormOnly_default_),
619  orthoType_(orthoType_default_),
620  impResScale_(impResScale_default_),
621  expResScale_(expResScale_default_),
622  label_(label_default_),
623  isSet_(false),
624  isSTSet_(false),
625  expResTest_(false),
626  loaDetected_(false)
627 {
628  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null, std::invalid_argument, "Problem not given to solver manager.");
629 
630  // If the parameter list pointer is null, then set the current parameters to the default parameter list.
631  if (!is_null(pl)) {
632  // Set the parameters using the list that was passed in.
633  setParameters( pl );
634  }
635 }
636 
637 template<class ScalarType, class MV, class OP>
638 void
640 setParameters (const Teuchos::RCP<Teuchos::ParameterList>& params)
641 {
642  using Teuchos::ParameterList;
643  using Teuchos::parameterList;
644  using Teuchos::rcp;
645  using Teuchos::rcp_dynamic_cast;
646 
647  // Create the internal parameter list if one doesn't already exist.
648  if (params_ == Teuchos::null) {
649  params_ = parameterList (*getValidParameters ());
650  } else {
651  params->validateParameters (*getValidParameters ());
652  }
653 
654  // Check for maximum number of restarts
655  if (params->isParameter ("Maximum Restarts")) {
656  maxRestarts_ = params->get ("Maximum Restarts", maxRestarts_default_);
657 
658  // Update parameter in our list.
659  params_->set ("Maximum Restarts", maxRestarts_);
660  }
661 
662  // Check for maximum number of iterations
663  if (params->isParameter ("Maximum Iterations")) {
664  maxIters_ = params->get ("Maximum Iterations", maxIters_default_);
665 
666  // Update parameter in our list and in status test.
667  params_->set ("Maximum Iterations", maxIters_);
668  if (! maxIterTest_.is_null ()) {
669  maxIterTest_->setMaxIters (maxIters_);
670  }
671  }
672 
673  // Check for blocksize
674  if (params->isParameter ("Block Size")) {
675  blockSize_ = params->get ("Block Size", blockSize_default_);
676  TEUCHOS_TEST_FOR_EXCEPTION(
677  blockSize_ <= 0, std::invalid_argument,
678  "Belos::PseudoBlockGmresSolMgr::setParameters: "
679  "The \"Block Size\" parameter must be strictly positive, "
680  "but you specified a value of " << blockSize_ << ".");
681 
682  // Update parameter in our list.
683  params_->set ("Block Size", blockSize_);
684  }
685 
686  // Check for the maximum number of blocks.
687  if (params->isParameter ("Num Blocks")) {
688  numBlocks_ = params->get ("Num Blocks", numBlocks_default_);
689  TEUCHOS_TEST_FOR_EXCEPTION(
690  numBlocks_ <= 0, std::invalid_argument,
691  "Belos::PseudoBlockGmresSolMgr::setParameters: "
692  "The \"Num Blocks\" parameter must be strictly positive, "
693  "but you specified a value of " << numBlocks_ << ".");
694 
695  // Update parameter in our list.
696  params_->set ("Num Blocks", numBlocks_);
697  }
698 
699  // Check to see if the timer label changed.
700  if (params->isParameter ("Timer Label")) {
701  const std::string tempLabel = params->get ("Timer Label", label_default_);
702 
703  // Update parameter in our list and solver timer
704  if (tempLabel != label_) {
705  label_ = tempLabel;
706  params_->set ("Timer Label", label_);
707  const std::string solveLabel =
708  label_ + ": PseudoBlockGmresSolMgr total solve time";
709 #ifdef BELOS_TEUCHOS_TIME_MONITOR
710  timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
711 #endif // BELOS_TEUCHOS_TIME_MONITOR
712  if (ortho_ != Teuchos::null) {
713  ortho_->setLabel( label_ );
714  }
715  }
716  }
717 
718  // Check if the orthogonalization changed.
719  if (params->isParameter ("Orthogonalization")) {
720  std::string tempOrthoType = params->get ("Orthogonalization", orthoType_default_);
721 #ifdef HAVE_BELOS_TSQR
722  TEUCHOS_TEST_FOR_EXCEPTION(
723  tempOrthoType != "DGKS" && tempOrthoType != "ICGS" &&
724  tempOrthoType != "IMGS" && tempOrthoType != "TSQR",
725  std::invalid_argument,
726  "Belos::PseudoBlockGmresSolMgr::setParameters: "
727  "The \"Orthogonalization\" parameter must be one of \"DGKS\", \"ICGS\", "
728  "\"IMGS\", or \"TSQR\".");
729 #else
730  TEUCHOS_TEST_FOR_EXCEPTION(
731  tempOrthoType != "DGKS" && tempOrthoType != "ICGS" &&
732  tempOrthoType != "IMGS",
733  std::invalid_argument,
734  "Belos::PseudoBlockGmresSolMgr::setParameters: "
735  "The \"Orthogonalization\" parameter must be one of \"DGKS\", \"ICGS\", "
736  "or \"IMGS\".");
737 #endif // HAVE_BELOS_TSQR
738 
739  if (tempOrthoType != orthoType_) {
740  orthoType_ = tempOrthoType;
741  // Create orthogonalization manager
742  if (orthoType_ == "DGKS") {
743  typedef DGKSOrthoManager<ScalarType, MV, OP> ortho_type;
744  if (orthoKappa_ <= 0) {
745  ortho_ = rcp (new ortho_type (label_));
746  }
747  else {
748  ortho_ = rcp (new ortho_type (label_));
749  rcp_dynamic_cast<ortho_type> (ortho_)->setDepTol (orthoKappa_);
750  }
751  }
752  else if (orthoType_ == "ICGS") {
753  typedef ICGSOrthoManager<ScalarType, MV, OP> ortho_type;
754  ortho_ = rcp (new ortho_type (label_));
755  }
756  else if (orthoType_ == "IMGS") {
757  typedef IMGSOrthoManager<ScalarType, MV, OP> ortho_type;
758  ortho_ = rcp (new ortho_type (label_));
759  }
760 #ifdef HAVE_BELOS_TSQR
761  else if (orthoType_ == "TSQR") {
762  typedef TsqrMatOrthoManager<ScalarType, MV, OP> ortho_type;
763  ortho_ = rcp (new ortho_type (label_));
764  }
765 #endif // HAVE_BELOS_TSQR
766  }
767  }
768 
769  // Check which orthogonalization constant to use.
770  if (params->isParameter ("Orthogonalization Constant")) {
771  orthoKappa_ = params->get ("Orthogonalization Constant", orthoKappa_default_);
772 
773  // Update parameter in our list.
774  params_->set ("Orthogonalization Constant", orthoKappa_);
775  if (orthoType_ == "DGKS") {
776  if (orthoKappa_ > 0 && ! ortho_.is_null ()) {
777  typedef DGKSOrthoManager<ScalarType, MV, OP> ortho_type;
778  rcp_dynamic_cast<ortho_type> (ortho_)->setDepTol (orthoKappa_);
779  }
780  }
781  }
782 
783  // Check for a change in verbosity level
784  if (params->isParameter ("Verbosity")) {
785  if (Teuchos::isParameterType<int> (*params, "Verbosity")) {
786  verbosity_ = params->get ("Verbosity", verbosity_default_);
787  } else {
788  verbosity_ = (int) Teuchos::getParameter<Belos::MsgType> (*params, "Verbosity");
789  }
790 
791  // Update parameter in our list.
792  params_->set ("Verbosity", verbosity_);
793  if (! printer_.is_null ()) {
794  printer_->setVerbosity (verbosity_);
795  }
796  }
797 
798  // Check for a change in output style.
799  if (params->isParameter ("Output Style")) {
800  if (Teuchos::isParameterType<int> (*params, "Output Style")) {
801  outputStyle_ = params->get ("Output Style", outputStyle_default_);
802  } else {
803  outputStyle_ = (int) Teuchos::getParameter<Belos::OutputType> (*params, "Output Style");
804  }
805 
806  // Update parameter in our list.
807  params_->set ("Output Style", verbosity_);
808  if (! outputTest_.is_null ()) {
809  isSTSet_ = false;
810  }
811  }
812 
813  // output stream
814  if (params->isParameter ("Output Stream")) {
815  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> > (*params, "Output Stream");
816 
817  // Update parameter in our list.
818  params_->set("Output Stream", outputStream_);
819  if (! printer_.is_null ()) {
820  printer_->setOStream (outputStream_);
821  }
822  }
823 
824  // frequency level
825  if (verbosity_ & Belos::StatusTestDetails) {
826  if (params->isParameter ("Output Frequency")) {
827  outputFreq_ = params->get ("Output Frequency", outputFreq_default_);
828  }
829 
830  // Update parameter in out list and output status test.
831  params_->set ("Output Frequency", outputFreq_);
832  if (! outputTest_.is_null ()) {
833  outputTest_->setOutputFrequency (outputFreq_);
834  }
835  }
836 
837  // Create output manager if we need to.
838  if (printer_.is_null ()) {
839  printer_ = rcp (new OutputManager<ScalarType> (verbosity_, outputStream_));
840  }
841 
842  // Convergence
843 
844  // Check for convergence tolerance
845  if (params->isParameter ("Convergence Tolerance")) {
846  convtol_ = params->get ("Convergence Tolerance", convtol_default_);
847 
848  // Update parameter in our list and residual tests.
849  params_->set ("Convergence Tolerance", convtol_);
850  if (! impConvTest_.is_null ()) {
851  impConvTest_->setTolerance (convtol_);
852  }
853  if (! expConvTest_.is_null ()) {
854  expConvTest_->setTolerance (convtol_);
855  }
856  }
857 
858  // Check for a change in scaling, if so we need to build new residual tests.
859  if (params->isParameter ("Implicit Residual Scaling")) {
860  const std::string tempImpResScale =
861  Teuchos::getParameter<std::string> (*params, "Implicit Residual Scaling");
862 
863  // Only update the scaling if it's different.
864  if (impResScale_ != tempImpResScale) {
865  Belos::ScaleType impResScaleType = convertStringToScaleType (tempImpResScale);
866  impResScale_ = tempImpResScale;
867 
868  // Update parameter in our list and residual tests
869  params_->set ("Implicit Residual Scaling", impResScale_);
870  if (! impConvTest_.is_null ()) {
871  try {
872  impConvTest_->defineScaleForm (impResScaleType, Belos::TwoNorm);
873  }
874  catch (std::exception& e) {
875  // Make sure the convergence test gets constructed again.
876  isSTSet_ = false;
877  }
878  }
879  }
880  }
881 
882  if (params->isParameter ("Explicit Residual Scaling")) {
883  const std::string tempExpResScale =
884  Teuchos::getParameter<std::string> (*params, "Explicit Residual Scaling");
885 
886  // Only update the scaling if it's different.
887  if (expResScale_ != tempExpResScale) {
888  Belos::ScaleType expResScaleType = convertStringToScaleType (tempExpResScale);
889  expResScale_ = tempExpResScale;
890 
891  // Update parameter in our list and residual tests
892  params_->set ("Explicit Residual Scaling", expResScale_);
893  if (! expConvTest_.is_null ()) {
894  try {
895  expConvTest_->defineScaleForm (expResScaleType, Belos::TwoNorm);
896  }
897  catch (std::exception& e) {
898  // Make sure the convergence test gets constructed again.
899  isSTSet_ = false;
900  }
901  }
902  }
903  }
904 
905 
906  if (params->isParameter ("Show Maximum Residual Norm Only")) {
907  showMaxResNormOnly_ =
908  Teuchos::getParameter<bool> (*params, "Show Maximum Residual Norm Only");
909 
910  // Update parameter in our list and residual tests
911  params_->set ("Show Maximum Residual Norm Only", showMaxResNormOnly_);
912  if (! impConvTest_.is_null ()) {
913  impConvTest_->setShowMaxResNormOnly (showMaxResNormOnly_);
914  }
915  if (! expConvTest_.is_null ()) {
916  expConvTest_->setShowMaxResNormOnly (showMaxResNormOnly_);
917  }
918  }
919 
920  // Create status tests if we need to.
921 
922  // Get the deflation quorum, or number of converged systems before deflation is allowed
923  if (params->isParameter("Deflation Quorum")) {
924  defQuorum_ = params->get("Deflation Quorum", defQuorum_);
925  TEUCHOS_TEST_FOR_EXCEPTION(
926  defQuorum_ > blockSize_, std::invalid_argument,
927  "Belos::PseudoBlockGmresSolMgr::setParameters: "
928  "The \"Deflation Quorum\" parameter (= " << defQuorum_ << ") must not be "
929  "larger than \"Block Size\" (= " << blockSize_ << ").");
930  params_->set ("Deflation Quorum", defQuorum_);
931  if (! impConvTest_.is_null ()) {
932  impConvTest_->setQuorum (defQuorum_);
933  }
934  if (! expConvTest_.is_null ()) {
935  expConvTest_->setQuorum (defQuorum_);
936  }
937  }
938 
939  // Create orthogonalization manager if we need to.
940  if (ortho_.is_null ()) {
941  if (orthoType_ == "DGKS") {
942  typedef DGKSOrthoManager<ScalarType, MV, OP> ortho_type;
943  if (orthoKappa_ <= 0) {
944  ortho_ = rcp (new ortho_type (label_));
945  }
946  else {
947  ortho_ = rcp (new ortho_type (label_));
948  rcp_dynamic_cast<ortho_type> (ortho_)->setDepTol (orthoKappa_);
949  }
950  }
951  else if (orthoType_ == "ICGS") {
952  typedef ICGSOrthoManager<ScalarType, MV, OP> ortho_type;
953  ortho_ = rcp (new ortho_type (label_));
954  }
955  else if (orthoType_ == "IMGS") {
956  typedef IMGSOrthoManager<ScalarType, MV, OP> ortho_type;
957  ortho_ = rcp (new ortho_type (label_));
958  }
959 #ifdef HAVE_BELOS_TSQR
960  else if (orthoType_ == "TSQR") {
961  typedef TsqrMatOrthoManager<ScalarType, MV, OP> ortho_type;
962  ortho_ = rcp (new ortho_type (label_));
963  }
964 #endif // HAVE_BELOS_TSQR
965  else {
966 #ifdef HAVE_BELOS_TSQR
967  TEUCHOS_TEST_FOR_EXCEPTION(
968  orthoType_ != "ICGS" && orthoType_ != "DGKS" &&
969  orthoType_ != "IMGS" && orthoType_ != "TSQR",
970  std::logic_error,
971  "Belos::PseudoBlockGmresSolMgr::setParameters(): "
972  "Invalid orthogonalization type \"" << orthoType_ << "\".");
973 #else
974  TEUCHOS_TEST_FOR_EXCEPTION(
975  orthoType_ != "ICGS" && orthoType_ != "DGKS" &&
976  orthoType_ != "IMGS",
977  std::logic_error,
978  "Belos::PseudoBlockGmresSolMgr::setParameters(): "
979  "Invalid orthogonalization type \"" << orthoType_ << "\".");
980 #endif // HAVE_BELOS_TSQR
981  }
982  }
983 
984  // Create the timer if we need to.
985  if (timerSolve_ == Teuchos::null) {
986  std::string solveLabel = label_ + ": PseudoBlockGmresSolMgr total solve time";
987 #ifdef BELOS_TEUCHOS_TIME_MONITOR
988  timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
989 #endif
990  }
991 
992  // Inform the solver manager that the current parameters were set.
993  isSet_ = true;
994 }
995 
996 
997 template<class ScalarType, class MV, class OP>
998 void
1000  const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &userConvStatusTest
1001  )
1002 {
1003  userConvStatusTest_ = userConvStatusTest;
1004 }
1005 
1006 template<class ScalarType, class MV, class OP>
1007 void
1009  const Teuchos::RCP<StatusTest<ScalarType,MV,OP> > &debugStatusTest
1010  )
1011 {
1012  debugStatusTest_ = debugStatusTest;
1013 }
1014 
1015 
1016 
1017 template<class ScalarType, class MV, class OP>
1018 Teuchos::RCP<const Teuchos::ParameterList>
1020 {
1021  static Teuchos::RCP<const Teuchos::ParameterList> validPL;
1022  if (is_null(validPL)) {
1023  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
1024  // Set all the valid parameters and their default values.
1025  pl= Teuchos::rcp( new Teuchos::ParameterList() );
1026  pl->set("Convergence Tolerance", convtol_default_,
1027  "The relative residual tolerance that needs to be achieved by the\n"
1028  "iterative solver in order for the linear system to be declared converged.");
1029  pl->set("Maximum Restarts", maxRestarts_default_,
1030  "The maximum number of restarts allowed for each\n"
1031  "set of RHS solved.");
1032  pl->set("Maximum Iterations", maxIters_default_,
1033  "The maximum number of block iterations allowed for each\n"
1034  "set of RHS solved.");
1035  pl->set("Num Blocks", numBlocks_default_,
1036  "The maximum number of vectors allowed in the Krylov subspace\n"
1037  "for each set of RHS solved.");
1038  pl->set("Block Size", blockSize_default_,
1039  "The number of RHS solved simultaneously.");
1040  pl->set("Verbosity", verbosity_default_,
1041  "What type(s) of solver information should be outputted\n"
1042  "to the output stream.");
1043  pl->set("Output Style", outputStyle_default_,
1044  "What style is used for the solver information outputted\n"
1045  "to the output stream.");
1046  pl->set("Output Frequency", outputFreq_default_,
1047  "How often convergence information should be outputted\n"
1048  "to the output stream.");
1049  pl->set("Deflation Quorum", defQuorum_default_,
1050  "The number of linear systems that need to converge before\n"
1051  "they are deflated. This number should be <= block size.");
1052  pl->set("Output Stream", outputStream_default_,
1053  "A reference-counted pointer to the output stream where all\n"
1054  "solver output is sent.");
1055  pl->set("Show Maximum Residual Norm Only", showMaxResNormOnly_default_,
1056  "When convergence information is printed, only show the maximum\n"
1057  "relative residual norm when the block size is greater than one.");
1058  pl->set("Implicit Residual Scaling", impResScale_default_,
1059  "The type of scaling used in the implicit residual convergence test.");
1060  pl->set("Explicit Residual Scaling", expResScale_default_,
1061  "The type of scaling used in the explicit residual convergence test.");
1062  pl->set("Timer Label", label_default_,
1063  "The string to use as a prefix for the timer labels.");
1064  // defaultParams_->set("Restart Timers", restartTimers_);
1065  pl->set("Orthogonalization", orthoType_default_,
1066  "The type of orthogonalization to use: DGKS, ICGS, IMGS.");
1067  pl->set("Orthogonalization Constant",orthoKappa_default_,
1068  "The constant used by DGKS orthogonalization to determine\n"
1069  "whether another step of classical Gram-Schmidt is necessary.");
1070  validPL = pl;
1071  }
1072  return validPL;
1073 }
1074 
1075 // Check the status test versus the defined linear problem
1076 template<class ScalarType, class MV, class OP>
1078 
1079  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
1080  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestGenResNorm_t;
1081  typedef Belos::StatusTestImpResNorm<ScalarType,MV,OP> StatusTestImpResNorm_t;
1082 
1083  // Basic test checks maximum iterations and native residual.
1084  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
1085 
1086  // If there is a left preconditioner, we create a combined status test that checks the implicit
1087  // and then explicit residual norm to see if we have convergence.
1088  if ( !Teuchos::is_null(problem_->getLeftPrec()) ) {
1089  expResTest_ = true;
1090  }
1091 
1092  if (expResTest_) {
1093 
1094  // Implicit residual test, using the native residual to determine if convergence was achieved.
1095  Teuchos::RCP<StatusTestGenResNorm_t> tmpImpConvTest =
1096  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
1097  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
1098  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
1099  impConvTest_ = tmpImpConvTest;
1100 
1101  // Explicit residual test once the native residual is below the tolerance
1102  Teuchos::RCP<StatusTestGenResNorm_t> tmpExpConvTest =
1103  Teuchos::rcp( new StatusTestGenResNorm_t( convtol_, defQuorum_ ) );
1104  tmpExpConvTest->defineResForm( StatusTestGenResNorm_t::Explicit, Belos::TwoNorm );
1105  tmpExpConvTest->defineScaleForm( convertStringToScaleType(expResScale_), Belos::TwoNorm );
1106  tmpExpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
1107  expConvTest_ = tmpExpConvTest;
1108 
1109  // The convergence test is a combination of the "cheap" implicit test and explicit test.
1110  convTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::SEQ, impConvTest_, expConvTest_ ) );
1111  }
1112  else {
1113 
1114  // Implicit residual test, using the native residual to determine if convergence was achieved.
1115  // Use test that checks for loss of accuracy.
1116  Teuchos::RCP<StatusTestImpResNorm_t> tmpImpConvTest =
1117  Teuchos::rcp( new StatusTestImpResNorm_t( convtol_, defQuorum_ ) );
1118  tmpImpConvTest->defineScaleForm( convertStringToScaleType(impResScale_), Belos::TwoNorm );
1119  tmpImpConvTest->setShowMaxResNormOnly( showMaxResNormOnly_ );
1120  impConvTest_ = tmpImpConvTest;
1121 
1122  // Set the explicit and total convergence test to this implicit test that checks for accuracy loss.
1123  expConvTest_ = impConvTest_;
1124  convTest_ = impConvTest_;
1125  }
1126 
1127  if (nonnull(debugStatusTest_) ) {
1128  // Add debug convergence test
1129  convTest_ = Teuchos::rcp(
1130  new StatusTestCombo_t( StatusTestCombo_t::AND, convTest_, debugStatusTest_ ) );
1131  }
1132 
1133  if (nonnull(userConvStatusTest_) ) {
1134  // Override the overall convergence test with the users convergence test
1135  convTest_ = Teuchos::rcp(
1136  new StatusTestCombo_t( StatusTestCombo_t::SEQ, convTest_, userConvStatusTest_ ) );
1137  // NOTE: Above, you have to run the other convergence tests also because
1138  // the logic in this class depends on it. This is very unfortunate.
1139  }
1140 
1141  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
1142 
1143  // Create the status test output class.
1144  // This class manages and formats the output from the status test.
1145  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
1146  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
1147 
1148  // Set the solver string for the output test
1149  std::string solverDesc = " Pseudo Block Gmres ";
1150  outputTest_->setSolverDesc( solverDesc );
1151 
1152 
1153  // The status test is now set.
1154  isSTSet_ = true;
1155 
1156  return false;
1157 }
1158 
1159 
1160 // solve()
1161 template<class ScalarType, class MV, class OP>
1163 
1164  // Set the current parameters if they were not set before.
1165  // NOTE: This may occur if the user generated the solver manager with the default constructor and
1166  // then didn't set any parameters using setParameters().
1167  if (!isSet_) { setParameters( params_ ); }
1168 
1169  Teuchos::BLAS<int,ScalarType> blas;
1170 
1171  TEUCHOS_TEST_FOR_EXCEPTION(!problem_->isProblemSet(),PseudoBlockGmresSolMgrLinearProblemFailure,
1172  "Belos::PseudoBlockGmresSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
1173 
1174  // Check if we have to create the status tests.
1175  if (!isSTSet_ || (!expResTest_ && !Teuchos::is_null(problem_->getLeftPrec())) ) {
1176  TEUCHOS_TEST_FOR_EXCEPTION( checkStatusTest(),PseudoBlockGmresSolMgrLinearProblemFailure,
1177  "Belos::BlockGmresSolMgr::solve(): Linear problem and requested status tests are incompatible.");
1178  }
1179 
1180  // Create indices for the linear systems to be solved.
1181  int startPtr = 0;
1182  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
1183  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1184 
1185  std::vector<int> currIdx( numCurrRHS );
1186  blockSize_ = numCurrRHS;
1187  for (int i=0; i<numCurrRHS; ++i)
1188  { currIdx[i] = startPtr+i; }
1189 
1190  // Inform the linear problem of the current linear system to solve.
1191  problem_->setLSIndex( currIdx );
1192 
1194  // Parameter list
1195  Teuchos::ParameterList plist;
1196  plist.set("Num Blocks",numBlocks_);
1197 
1198  // Reset the status test.
1199  outputTest_->reset();
1200  loaDetected_ = false;
1201 
1202  // Assume convergence is achieved, then let any failed convergence set this to false.
1203  bool isConverged = true;
1204 
1206  // BlockGmres solver
1207 
1208  Teuchos::RCP<PseudoBlockGmresIter<ScalarType,MV,OP> > block_gmres_iter
1209  = Teuchos::rcp( new PseudoBlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
1210 
1211  // Enter solve() iterations
1212  {
1213 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1214  Teuchos::TimeMonitor slvtimer(*timerSolve_);
1215 #endif
1216 
1217  while ( numRHS2Solve > 0 ) {
1218 
1219  // Reset the active / converged vectors from this block
1220  std::vector<int> convRHSIdx;
1221  std::vector<int> currRHSIdx( currIdx );
1222  currRHSIdx.resize(numCurrRHS);
1223 
1224  // Set the current number of blocks with the pseudo Gmres iteration.
1225  block_gmres_iter->setNumBlocks( numBlocks_ );
1226 
1227  // Reset the number of iterations.
1228  block_gmres_iter->resetNumIters();
1229 
1230  // Reset the number of calls that the status test output knows about.
1231  outputTest_->resetNumCalls();
1232 
1233  // Get a new state struct and initialize the solver.
1235 
1236  // Create the first block in the current Krylov basis for each right-hand side.
1237  std::vector<int> index(1);
1238  Teuchos::RCP<MV> tmpV, R_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
1239  newState.V.resize( blockSize_ );
1240  newState.Z.resize( blockSize_ );
1241  for (int i=0; i<blockSize_; ++i) {
1242  index[0]=i;
1243  tmpV = MVT::CloneViewNonConst( *R_0, index );
1244 
1245  // Get a matrix to hold the orthonormalization coefficients.
1246  Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > tmpZ
1247  = Teuchos::rcp( new Teuchos::SerialDenseVector<int,ScalarType>( 1 ));
1248 
1249  // Orthonormalize the new V_0
1250  int rank = ortho_->normalize( *tmpV, tmpZ );
1251  TEUCHOS_TEST_FOR_EXCEPTION(rank != 1, PseudoBlockGmresSolMgrOrthoFailure,
1252  "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors.");
1253 
1254  newState.V[i] = tmpV;
1255  newState.Z[i] = tmpZ;
1256  }
1257 
1258  newState.curDim = 0;
1259  block_gmres_iter->initialize(newState);
1260  int numRestarts = 0;
1261 
1262  while(1) {
1263 
1264  // tell block_gmres_iter to iterate
1265  try {
1266  block_gmres_iter->iterate();
1267 
1269  //
1270  // check convergence first
1271  //
1273  if ( convTest_->getStatus() == Passed ) {
1274 
1275  if ( expConvTest_->getLOADetected() ) {
1276  // Oops! There was a loss of accuracy (LOA) for one or
1277  // more right-hand sides. That means the implicit
1278  // (a.k.a. "native") residuals claim convergence,
1279  // whereas the explicit (hence expConvTest_, i.e.,
1280  // "explicit convergence test") residuals have not
1281  // converged.
1282  //
1283  // We choose to deal with this situation by deflating
1284  // out the affected right-hand sides and moving on.
1285  loaDetected_ = true;
1286  printer_->stream(Warnings) <<
1287  "Belos::PseudoBlockGmresSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1288  isConverged = false;
1289  }
1290 
1291  // Figure out which linear systems converged.
1292  std::vector<int> convIdx = expConvTest_->convIndices();
1293 
1294  // If the number of converged linear systems is equal to the
1295  // number of current linear systems, then we are done with this block.
1296  if (convIdx.size() == currRHSIdx.size())
1297  break; // break from while(1){block_gmres_iter->iterate()}
1298 
1299  // Get a new state struct and initialize the solver.
1301 
1302  // Inform the linear problem that we are finished with this current linear system.
1303  problem_->setCurrLS();
1304 
1305  // Get the state.
1306  PseudoBlockGmresIterState<ScalarType,MV> oldState = block_gmres_iter->getState();
1307  int curDim = oldState.curDim;
1308 
1309  // Get a new state struct and reset currRHSIdx to have the right-hand sides that
1310  // are left to converge for this block.
1311  int have = 0;
1312  std::vector<int> oldRHSIdx( currRHSIdx );
1313  std::vector<int> defRHSIdx;
1314  for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
1315  bool found = false;
1316  for (unsigned int j=0; j<convIdx.size(); ++j) {
1317  if (currRHSIdx[i] == convIdx[j]) {
1318  found = true;
1319  break;
1320  }
1321  }
1322  if (found) {
1323  defRHSIdx.push_back( i );
1324  }
1325  else {
1326  defState.V.push_back( Teuchos::rcp_const_cast<MV>( oldState.V[i] ) );
1327  defState.Z.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,ScalarType> >( oldState.Z[i] ) );
1328  defState.H.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseMatrix<int,ScalarType> >( oldState.H[i] ) );
1329  defState.sn.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,ScalarType> >( oldState.sn[i] ) );
1330  defState.cs.push_back( Teuchos::rcp_const_cast<Teuchos::SerialDenseVector<int,MagnitudeType> >(oldState.cs[i] ) );
1331  currRHSIdx[have] = currRHSIdx[i];
1332  have++;
1333  }
1334  }
1335  defRHSIdx.resize(currRHSIdx.size()-have);
1336  currRHSIdx.resize(have);
1337 
1338  // Compute the current solution that needs to be deflated if this solver has taken any steps.
1339  if (curDim) {
1340  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1341  Teuchos::RCP<MV> defUpdate = MVT::CloneViewNonConst( *update, defRHSIdx );
1342 
1343  // Set the deflated indices so we can update the solution.
1344  problem_->setLSIndex( convIdx );
1345 
1346  // Update the linear problem.
1347  problem_->updateSolution( defUpdate, true );
1348  }
1349 
1350  // Set the remaining indices after deflation.
1351  problem_->setLSIndex( currRHSIdx );
1352 
1353  // Set the dimension of the subspace, which is the same as the old subspace size.
1354  defState.curDim = curDim;
1355 
1356  // Initialize the solver with the deflated system.
1357  block_gmres_iter->initialize(defState);
1358  }
1360  //
1361  // check for maximum iterations
1362  //
1364  else if ( maxIterTest_->getStatus() == Passed ) {
1365  // we don't have convergence
1366  isConverged = false;
1367  break; // break from while(1){block_gmres_iter->iterate()}
1368  }
1370  //
1371  // check for restarting, i.e. the subspace is full
1372  //
1374  else if ( block_gmres_iter->getCurSubspaceDim() == block_gmres_iter->getMaxSubspaceDim() ) {
1375 
1376  if ( numRestarts >= maxRestarts_ ) {
1377  isConverged = false;
1378  break; // break from while(1){block_gmres_iter->iterate()}
1379  }
1380  numRestarts++;
1381 
1382  printer_->stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
1383 
1384  // Update the linear problem.
1385  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1386  problem_->updateSolution( update, true );
1387 
1388  // Get the state.
1389  PseudoBlockGmresIterState<ScalarType,MV> oldState = block_gmres_iter->getState();
1390 
1391  // Set the new state.
1393  newstate.V.resize(currRHSIdx.size());
1394  newstate.Z.resize(currRHSIdx.size());
1395 
1396  // Compute the restart vectors
1397  // NOTE: Force the linear problem to update the current residual since the solution was updated.
1398  R_0 = MVT::Clone( *(problem_->getInitPrecResVec()), currRHSIdx.size() );
1399  problem_->computeCurrPrecResVec( &*R_0 );
1400  for (unsigned int i=0; i<currRHSIdx.size(); ++i) {
1401  index[0] = i; // index(1) vector declared on line 891
1402 
1403  tmpV = MVT::CloneViewNonConst( *R_0, index );
1404 
1405  // Get a matrix to hold the orthonormalization coefficients.
1406  Teuchos::RCP<Teuchos::SerialDenseVector<int,ScalarType> > tmpZ
1407  = Teuchos::rcp( new Teuchos::SerialDenseVector<int,ScalarType>( 1 ));
1408 
1409  // Orthonormalize the new V_0
1410  int rank = ortho_->normalize( *tmpV, tmpZ );
1411  TEUCHOS_TEST_FOR_EXCEPTION(rank != 1 ,PseudoBlockGmresSolMgrOrthoFailure,
1412  "Belos::PseudoBlockGmresSolMgr::solve(): Failed to compute initial block of orthonormal vectors after the restart.");
1413 
1414  newstate.V[i] = tmpV;
1415  newstate.Z[i] = tmpZ;
1416  }
1417 
1418  // Initialize the solver.
1419  newstate.curDim = 0;
1420  block_gmres_iter->initialize(newstate);
1421 
1422  } // end of restarting
1423 
1425  //
1426  // we returned from iterate(), but none of our status tests Passed.
1427  // something is wrong, and it is probably our fault.
1428  //
1430 
1431  else {
1432  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
1433  "Belos::PseudoBlockGmresSolMgr::solve(): Invalid return from PseudoBlockGmresIter::iterate().");
1434  }
1435  }
1436  catch (const PseudoBlockGmresIterOrthoFailure &e) {
1437 
1438  // Try to recover the most recent least-squares solution
1439  block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
1440 
1441  // Check to see if the most recent least-squares solution yielded convergence.
1442  sTest_->checkStatus( &*block_gmres_iter );
1443  if (convTest_->getStatus() != Passed)
1444  isConverged = false;
1445  break;
1446  }
1447  catch (const std::exception &e) {
1448  printer_->stream(Errors) << "Error! Caught std::exception in PseudoBlockGmresIter::iterate() at iteration "
1449  << block_gmres_iter->getNumIters() << std::endl
1450  << e.what() << std::endl;
1451  throw;
1452  }
1453  }
1454 
1455  // Compute the current solution.
1456  // Update the linear problem.
1457  if (nonnull(userConvStatusTest_)) {
1458  //std::cout << "\nnonnull(userConvStatusTest_)\n";
1459  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1460  problem_->updateSolution( update, true );
1461  }
1462  else if (nonnull(expConvTest_->getSolution())) {
1463  //std::cout << "\nexpConvTest_->getSolution()\n";
1464  Teuchos::RCP<MV> newX = expConvTest_->getSolution();
1465  Teuchos::RCP<MV> curX = problem_->getCurrLHSVec();
1466  MVT::MvAddMv( 0.0, *newX, 1.0, *newX, *curX );
1467  }
1468  else {
1469  //std::cout << "\nblock_gmres_iter->getCurrentUpdate()\n";
1470  Teuchos::RCP<MV> update = block_gmres_iter->getCurrentUpdate();
1471  problem_->updateSolution( update, true );
1472  }
1473 
1474  // Inform the linear problem that we are finished with this block linear system.
1475  problem_->setCurrLS();
1476 
1477  // Update indices for the linear systems to be solved.
1478  startPtr += numCurrRHS;
1479  numRHS2Solve -= numCurrRHS;
1480  if ( numRHS2Solve > 0 ) {
1481  numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1482 
1483  blockSize_ = numCurrRHS;
1484  currIdx.resize( numCurrRHS );
1485  for (int i=0; i<numCurrRHS; ++i)
1486  { currIdx[i] = startPtr+i; }
1487 
1488  // Adapt the status test quorum if we need to.
1489  if (defQuorum_ > blockSize_) {
1490  if (impConvTest_ != Teuchos::null)
1491  impConvTest_->setQuorum( blockSize_ );
1492  if (expConvTest_ != Teuchos::null)
1493  expConvTest_->setQuorum( blockSize_ );
1494  }
1495 
1496  // Set the next indices.
1497  problem_->setLSIndex( currIdx );
1498  }
1499  else {
1500  currIdx.resize( numRHS2Solve );
1501  }
1502 
1503  }// while ( numRHS2Solve > 0 )
1504 
1505  }
1506 
1507  // print final summary
1508  sTest_->print( printer_->stream(FinalSummary) );
1509 
1510  // print timing information
1511 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1512  // Calling summarize() can be expensive, so don't call unless the
1513  // user wants to print out timing details. summarize() will do all
1514  // the work even if it's passed a "black hole" output stream.
1515  if (verbosity_ & TimingDetails)
1516  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1517 #endif
1518 
1519  // get iteration information for this solve
1520  numIters_ = maxIterTest_->getNumIters();
1521 
1522  // Save the convergence test value ("achieved tolerance") for this
1523  // solve. For this solver, convTest_ may either be a single
1524  // residual norm test, or a combination of two residual norm tests.
1525  // In the latter case, the master convergence test convTest_ is a
1526  // SEQ combo of the implicit resp. explicit tests. If the implicit
1527  // test never passes, then the explicit test won't ever be executed.
1528  // This manifests as expConvTest_->getTestValue()->size() < 1. We
1529  // deal with this case by using the values returned by
1530  // impConvTest_->getTestValue().
1531  {
1532  // We'll fetch the vector of residual norms one way or the other.
1533  const std::vector<MagnitudeType>* pTestValues = NULL;
1534  if (expResTest_) {
1535  pTestValues = expConvTest_->getTestValue();
1536  if (pTestValues == NULL || pTestValues->size() < 1) {
1537  pTestValues = impConvTest_->getTestValue();
1538  }
1539  }
1540  else {
1541  // Only the implicit residual norm test is being used.
1542  pTestValues = impConvTest_->getTestValue();
1543  }
1544  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1545  "Belos::PseudoBlockGmresSolMgr::solve(): The implicit convergence test's "
1546  "getTestValue() method returned NULL. Please report this bug to the "
1547  "Belos developers.");
1548  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1549  "Belos::PseudoBlockGmresSolMgr::solve(): The implicit convergence test's "
1550  "getTestValue() method returned a vector of length zero. Please report "
1551  "this bug to the Belos developers.");
1552 
1553  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1554  // achieved tolerances for all vectors in the current solve(), or
1555  // just for the vectors from the last deflation?
1556  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1557  }
1558 
1559  if (!isConverged || loaDetected_) {
1560  return Unconverged; // return from PseudoBlockGmresSolMgr::solve()
1561  }
1562  return Converged; // return from PseudoBlockGmresSolMgr::solve()
1563 }
1564 
1565 
1566 template<class ScalarType, class MV, class OP>
1568 {
1569  std::ostringstream out;
1570 
1571  out << "\"Belos::PseudoBlockGmresSolMgr\": {";
1572  if (this->getObjectLabel () != "") {
1573  out << "Label: " << this->getObjectLabel () << ", ";
1574  }
1575  out << "Num Blocks: " << numBlocks_
1576  << ", Maximum Iterations: " << maxIters_
1577  << ", Maximum Restarts: " << maxRestarts_
1578  << ", Convergence Tolerance: " << convtol_
1579  << "}";
1580  return out.str ();
1581 }
1582 
1583 
1584 template<class ScalarType, class MV, class OP>
1585 void
1587 describe (Teuchos::FancyOStream &out,
1588  const Teuchos::EVerbosityLevel verbLevel) const
1589 {
1590  using Teuchos::TypeNameTraits;
1591  using Teuchos::VERB_DEFAULT;
1592  using Teuchos::VERB_NONE;
1593  using Teuchos::VERB_LOW;
1594  // using Teuchos::VERB_MEDIUM;
1595  // using Teuchos::VERB_HIGH;
1596  // using Teuchos::VERB_EXTREME;
1597  using std::endl;
1598 
1599  // Set default verbosity if applicable.
1600  const Teuchos::EVerbosityLevel vl =
1601  (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
1602 
1603  if (vl != VERB_NONE) {
1604  Teuchos::OSTab tab0 (out);
1605 
1606  out << "\"Belos::PseudoBlockGmresSolMgr\":" << endl;
1607  Teuchos::OSTab tab1 (out);
1608  out << "Template parameters:" << endl;
1609  {
1610  Teuchos::OSTab tab2 (out);
1611  out << "ScalarType: " << TypeNameTraits<ScalarType>::name () << endl
1612  << "MV: " << TypeNameTraits<MV>::name () << endl
1613  << "OP: " << TypeNameTraits<OP>::name () << endl;
1614  }
1615  if (this->getObjectLabel () != "") {
1616  out << "Label: " << this->getObjectLabel () << endl;
1617  }
1618  out << "Num Blocks: " << numBlocks_ << endl
1619  << "Maximum Iterations: " << maxIters_ << endl
1620  << "Maximum Restarts: " << maxRestarts_ << endl
1621  << "Convergence Tolerance: " << convtol_ << endl;
1622  }
1623 }
1624 
1625 } // end Belos namespace
1626 
1627 #endif /* BELOS_PSEUDO_BLOCK_GMRES_SOLMGR_HPP */
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:87
std::vector< Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > > H
The current Hessenberg matrix.
Collection of types and exceptions used within the Belos solvers.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
The current parameters for this solver.
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.
virtual void setUserConvStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &userConvStatusTest)
Set a custom status test.
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params)
Set the parameters the solver manager should use to solve the linear problem.
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:119
Convergence test using the implicit residual norm(s), with an explicit residual norm(s) check for los...
PseudoBlockGmresSolMgrOrthoFailure(const std::string &what_arg)
A factory class for generating StatusTestOutput objects.
std::vector< Teuchos::RCP< const MV > > V
The current Krylov basis.
An implementation of StatusTestResNorm using a family of residual norms.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Return a reference to the linear problem being solved by this solver manager.
Belos::StatusTest class for specifying a maximum number of iterations.
A pure virtual class for defining the status tests for the Belos iterative solvers.
Structure to contain pointers to PseudoBlockGmresIter state variables.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
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 > CloneViewNonConst(MV &mv, const std::vector< int > &index)
Creates a new MV that shares the selected contents of mv (shallow copy).
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. ...
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.
PseudoBlockGmresSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i...
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, MagnitudeType > > > cs
A linear system to solve, and its associated information.
Class which describes the linear problem to be solved by the iterative solver.
Interface to standard and "pseudoblock" GMRES.
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 ...
Iterated Classical Gram-Schmidt (ICGS) implementation of the Belos::OrthoManager class.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
bool isLOADetected() const
Whether a "loss of accuracy" was detected during the last solve().
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
void reset(const ResetType type)
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > Z
The current right-hand side of the least squares system RY = Z.
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 concrete class for performing the pseudo-block GMRES iteration.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Belos::StatusTest for specifying an implicit residual norm stopping criteria that checks for loss of ...
virtual void setDebugStatusTest(const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &debugStatusTest)
Set a debug status test.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
A class for extending the status testing capabilities of Belos via logical combinations.
Orthogonalization manager based on Tall Skinny QR (TSQR)
Class which defines basic traits for the operator type.
int curDim
The current dimension of the reduction.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem to solve.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
std::vector< Teuchos::RCP< const Teuchos::SerialDenseVector< int, ScalarType > > > sn
The current Given&#39;s rotation coefficients.
ReturnType solve()
This method performs possibly repeated calls to the underlying linear solver&#39;s iterate() routine unti...
MatOrthoManager subclass using TSQR or DGKS.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
A list of valid default parameters for this solver.
int getNumIters() const
Iteration count for the most recent call to solve().
Belos header file which uses auto-configuration information to include necessary C++ headers...
std::string description() const
Return a one-line description of this object.
This class implements the pseudo-block GMRES iteration, where a block Krylov subspace is constructed ...
PseudoBlockGmresIterOrthoFailure is thrown when the orthogonalization manager is unable to generate o...
MagnitudeType achievedTol() const
Tolerance achieved by the last solve() invocation.
PseudoBlockGmresSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate...