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