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

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