Belos  Version of the Day
BelosPCPGSolMgr.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_PCPG_SOLMGR_HPP
43 #define BELOS_PCPG_SOLMGR_HPP
44 
48 
49 #include "BelosConfigDefs.hpp"
50 #include "BelosTypes.hpp"
51 #include "BelosLinearProblem.hpp"
52 #include "BelosSolverManager.hpp"
53 
54 #include "BelosPCPGIter.hpp"
55 
61 #include "BelosStatusTestCombo.hpp"
63 #include "BelosOutputManager.hpp"
64 #include "Teuchos_BLAS.hpp"
65 #include "Teuchos_LAPACK.hpp"
66 #ifdef BELOS_TEUCHOS_TIME_MONITOR
67 # include "Teuchos_TimeMonitor.hpp"
68 #endif
69 #if defined(HAVE_TEUCHOSCORE_CXX11)
70 # include <type_traits>
71 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
72 #include "Teuchos_TypeTraits.hpp"
73 
74 namespace Belos {
75 
77 
78 
86  PCPGSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg)
87  {}};
88 
94  class PCPGSolMgrOrthoFailure : public BelosError {public:
95  PCPGSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg)
96  {}};
97 
104  class PCPGSolMgrLAPACKFailure : public BelosError {public:
105  PCPGSolMgrLAPACKFailure(const std::string& what_arg) : BelosError(what_arg)
106  {}};
107 
114  class PCPGSolMgrRecyclingFailure : public BelosError {public:
115  PCPGSolMgrRecyclingFailure(const std::string& what_arg) : BelosError(what_arg)
116  {}};
117 
119 
120 
144 
145  // Partial specialization for complex ScalarType.
146  // This contains a trivial implementation.
147  // See discussion in the class documentation above.
148  //
149  // FIXME (mfh 09 Sep 2015) This also is a stub for types other than
150  // float or double.
151  template<class ScalarType, class MV, class OP,
152  const bool supportsScalarType =
154  ! Teuchos::ScalarTraits<ScalarType>::isComplex>
155  class PCPGSolMgr :
156  public Details::SolverManagerRequiresRealLapack<ScalarType, MV, OP,
157  Belos::Details::LapackSupportsScalar<ScalarType>::value &&
158  ! Teuchos::ScalarTraits<ScalarType>::isComplex>
159  {
160  static const bool scalarTypeIsSupported =
162  ! Teuchos::ScalarTraits<ScalarType>::isComplex;
163  typedef Details::SolverManagerRequiresRealLapack<ScalarType, MV, OP,
164  scalarTypeIsSupported> base_type;
165 
166  public:
168  base_type ()
169  {}
170  PCPGSolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
171  const Teuchos::RCP<Teuchos::ParameterList> &pl) :
172  base_type ()
173  {}
174  virtual ~PCPGSolMgr () {}
175  };
176 
177  template<class ScalarType, class MV, class OP>
178  class PCPGSolMgr<ScalarType, MV, OP, true> :
179  public Details::SolverManagerRequiresRealLapack<ScalarType, MV, OP, true> {
180  private:
183  typedef Teuchos::ScalarTraits<ScalarType> SCT;
184  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
185  typedef Teuchos::ScalarTraits<MagnitudeType> MT;
186 
187  public:
189 
190 
197  PCPGSolMgr();
198 
234  PCPGSolMgr( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
235  const Teuchos::RCP<Teuchos::ParameterList> &pl );
236 
238  virtual ~PCPGSolMgr() {};
240 
242 
243 
247  return *problem_;
248  }
249 
252  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
253 
256  Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const { return params_; }
257 
263  Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
264  return Teuchos::tuple(timerSolve_);
265  }
266 
272  MagnitudeType achievedTol() const {
273  return achievedTol_;
274  }
275 
277  int getNumIters() const {
278  return numIters_;
279  }
280 
283  bool isLOADetected() const { return false; }
284 
286 
288 
289 
291  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) { problem_ = problem; }
292 
294  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params );
295 
297 
299 
300 
304  void reset( const ResetType type ) { if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) problem_->setProblem(); }
306 
308 
309 
327  ReturnType solve();
328 
330 
333 
335  std::string description() const;
336 
338 
339  private:
340 
341  // In the A-inner product, perform an RRQR decomposition without using A unless absolutely necessary. Given
342  // the seed space U and C = A U, find U1 and C1 with span(U1)=span(U) such that C1'U1 = I maintaining C=AU.
343  int ARRQR(int numVecs, int numOrthVecs, const Teuchos::SerialDenseMatrix<int,ScalarType>& D);
344 
345  // Linear problem.
346  Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
347 
348  // Output manager.
349  Teuchos::RCP<OutputManager<ScalarType> > printer_;
350  Teuchos::RCP<std::ostream> outputStream_;
351 
352  // Status test.
353  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
354  Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
355  Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > convTest_;
356  Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
357 
358  // Orthogonalization manager.
359  Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
360 
361  // Current parameter list.
362  Teuchos::RCP<Teuchos::ParameterList> params_;
363 
364  // Default solver values.
365  static const MagnitudeType convtol_default_;
366  static const MagnitudeType orthoKappa_default_;
367  static const int maxIters_default_;
368  static const int deflatedBlocks_default_;
369  static const int savedBlocks_default_;
370  static const int verbosity_default_;
371  static const int outputStyle_default_;
372  static const int outputFreq_default_;
373  static const std::string label_default_;
374  static const std::string orthoType_default_;
375  static const Teuchos::RCP<std::ostream> outputStream_default_;
376 
377  //
378  // Current solver values.
379  //
380 
382  MagnitudeType convtol_;
383 
385  MagnitudeType orthoKappa_;
386 
388  MagnitudeType achievedTol_;
389 
391  int numIters_;
392 
394  int maxIters_;
395 
396  int deflatedBlocks_, savedBlocks_, verbosity_, outputStyle_, outputFreq_;
397  std::string orthoType_;
398 
399  // Recycled subspace, its image and the residual
400  Teuchos::RCP<MV> U_, C_, R_;
401 
402  // Actual dimension of current recycling subspace (<= savedBlocks_ )
403  int dimU_;
404 
405  // Timers.
406  std::string label_;
407  Teuchos::RCP<Teuchos::Time> timerSolve_;
408 
409  // Internal state variables.
410  bool isSet_;
411  };
412 
413 
414 // Default solver values.
415 template<class ScalarType, class MV, class OP>
418 
419 template<class ScalarType, class MV, class OP>
422 
423 template<class ScalarType, class MV, class OP>
425 
426 template<class ScalarType, class MV, class OP>
428 
429 template<class ScalarType, class MV, class OP>
431 
432 template<class ScalarType, class MV, class OP>
434 
435 template<class ScalarType, class MV, class OP>
437 
438 template<class ScalarType, class MV, class OP>
440 
441 template<class ScalarType, class MV, class OP>
442 const std::string PCPGSolMgr<ScalarType,MV,OP,true>::label_default_ = "Belos";
443 
444 template<class ScalarType, class MV, class OP>
446 
447 template<class ScalarType, class MV, class OP>
448 const Teuchos::RCP<std::ostream> PCPGSolMgr<ScalarType,MV,OP,true>::outputStream_default_ = Teuchos::rcp(&std::cout,false);
449 
450 
451 // Empty Constructor
452 template<class ScalarType, class MV, class OP>
454  outputStream_(outputStream_default_),
455  convtol_(convtol_default_),
456  orthoKappa_(orthoKappa_default_),
457  achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
458  numIters_(0),
459  maxIters_(maxIters_default_),
460  deflatedBlocks_(deflatedBlocks_default_),
461  savedBlocks_(savedBlocks_default_),
462  verbosity_(verbosity_default_),
463  outputStyle_(outputStyle_default_),
464  outputFreq_(outputFreq_default_),
465  orthoType_(orthoType_default_),
466  dimU_(0),
467  label_(label_default_),
468  isSet_(false)
469 {}
470 
471 
472 // Basic Constructor
473 template<class ScalarType, class MV, class OP>
475  const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
476  const Teuchos::RCP<Teuchos::ParameterList> &pl ) :
477  problem_(problem),
478  outputStream_(outputStream_default_),
479  convtol_(convtol_default_),
480  orthoKappa_(orthoKappa_default_),
481  achievedTol_(Teuchos::ScalarTraits<MagnitudeType>::zero()),
482  numIters_(0),
483  maxIters_(maxIters_default_),
484  deflatedBlocks_(deflatedBlocks_default_),
485  savedBlocks_(savedBlocks_default_),
486  verbosity_(verbosity_default_),
487  outputStyle_(outputStyle_default_),
488  outputFreq_(outputFreq_default_),
489  orthoType_(orthoType_default_),
490  dimU_(0),
491  label_(label_default_),
492  isSet_(false)
493 {
494  TEUCHOS_TEST_FOR_EXCEPTION(
495  problem_.is_null (), std::invalid_argument,
496  "Belos::PCPGSolMgr two-argument constructor: "
497  "'problem' is null. You must supply a non-null Belos::LinearProblem "
498  "instance when calling this constructor.");
499 
500  if (! pl.is_null ()) {
501  // Set the parameters using the list that was passed in.
502  setParameters (pl);
503  }
504 }
505 
506 
507 template<class ScalarType, class MV, class OP>
508 void PCPGSolMgr<ScalarType,MV,OP,true>::setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params )
509 {
510  // Create the internal parameter list if ones doesn't already exist.
511  if (params_ == Teuchos::null) {
512  params_ = Teuchos::rcp( new Teuchos::ParameterList(*getValidParameters()) );
513  }
514  else {
515  params->validateParameters(*getValidParameters());
516  }
517 
518  // Check for maximum number of iterations
519  if (params->isParameter("Maximum Iterations")) {
520  maxIters_ = params->get("Maximum Iterations",maxIters_default_);
521 
522  // Update parameter in our list and in status test.
523  params_->set("Maximum Iterations", maxIters_);
524  if (maxIterTest_!=Teuchos::null)
525  maxIterTest_->setMaxIters( maxIters_ );
526  }
527 
528  // Check for the maximum numbers of saved and deflated blocks.
529  if (params->isParameter("Num Saved Blocks")) {
530  savedBlocks_ = params->get("Num Saved Blocks",savedBlocks_default_);
531  TEUCHOS_TEST_FOR_EXCEPTION(savedBlocks_ <= 0, std::invalid_argument,
532  "Belos::PCPGSolMgr: \"Num Saved Blocks\" must be strictly positive.");
533 
534  // savedBlocks > number of matrix rows and columns, not known in parameters.
535  //TEUCHOS_TEST_FOR_EXCEPTION(savedBlocks_ >= maxIters_, std::invalid_argument,
536  //"Belos::PCPGSolMgr: \"Num Saved Blocks\" must be less than \"Maximum Iterations\".");
537 
538  // Update parameter in our list.
539  params_->set("Num Saved Blocks", savedBlocks_);
540  }
541  if (params->isParameter("Num Deflated Blocks")) {
542  deflatedBlocks_ = params->get("Num Deflated Blocks",deflatedBlocks_default_);
543  TEUCHOS_TEST_FOR_EXCEPTION(deflatedBlocks_ < 0, std::invalid_argument,
544  "Belos::PCPGSolMgr: \"Num Deflated Blocks\" must be positive.");
545 
546  TEUCHOS_TEST_FOR_EXCEPTION(deflatedBlocks_ > savedBlocks_, std::invalid_argument,
547  "Belos::PCPGSolMgr: \"Num Deflated Blocks\" must be <= \"Num Saved Blocks\".");
548 
549  // Update parameter in our list.
550  params_->set("Num Deflated Blocks", deflatedBlocks_);
551  }
552 
553  // Check to see if the timer label changed.
554  if (params->isParameter("Timer Label")) {
555  std::string tempLabel = params->get("Timer Label", label_default_);
556 
557  // Update parameter in our list and solver timer
558  if (tempLabel != label_) {
559  label_ = tempLabel;
560  params_->set("Timer Label", label_);
561  std::string solveLabel = label_ + ": PCPGSolMgr total solve time";
562 #ifdef BELOS_TEUCHOS_TIME_MONITOR
563  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
564 #endif
565  if (ortho_ != Teuchos::null) {
566  ortho_->setLabel( label_ );
567  }
568  }
569  }
570 
571  // Check if the orthogonalization changed.
572  if (params->isParameter("Orthogonalization")) {
573  std::string tempOrthoType = params->get("Orthogonalization",orthoType_default_);
574  TEUCHOS_TEST_FOR_EXCEPTION( tempOrthoType != "DGKS" && tempOrthoType != "ICGS" && tempOrthoType != "IMGS",
575  std::invalid_argument,
576  "Belos::PCPGSolMgr: \"Orthogonalization\" must be either \"DGKS\", \"ICGS\", or \"IMGS\".");
577  if (tempOrthoType != orthoType_) {
578  orthoType_ = tempOrthoType;
579  params_->set("Orthogonalization", orthoType_);
580  // Create orthogonalization manager
581  if (orthoType_=="DGKS") {
582  if (orthoKappa_ <= 0) {
583  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
584  }
585  else {
586  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
587  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
588  }
589  }
590  else if (orthoType_=="ICGS") {
591  ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
592  }
593  else if (orthoType_=="IMGS") {
594  ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
595  }
596  }
597  }
598 
599  // Check which orthogonalization constant to use.
600  if (params->isParameter("Orthogonalization Constant")) {
601  orthoKappa_ = params->get("Orthogonalization Constant",orthoKappa_default_);
602 
603  // Update parameter in our list.
604  params_->set("Orthogonalization Constant",orthoKappa_);
605  if (orthoType_=="DGKS") {
606  if (orthoKappa_ > 0 && ortho_ != Teuchos::null) {
607  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
608  }
609  }
610  }
611 
612  // Check for a change in verbosity level
613  if (params->isParameter("Verbosity")) {
614  if (Teuchos::isParameterType<int>(*params,"Verbosity")) {
615  verbosity_ = params->get("Verbosity", verbosity_default_);
616  } else {
617  verbosity_ = (int)Teuchos::getParameter<Belos::MsgType>(*params,"Verbosity");
618  }
619 
620  // Update parameter in our list.
621  params_->set("Verbosity", verbosity_);
622  if (printer_ != Teuchos::null)
623  printer_->setVerbosity(verbosity_);
624  }
625 
626  // Check for a change in output style
627  if (params->isParameter("Output Style")) {
628  if (Teuchos::isParameterType<int>(*params,"Output Style")) {
629  outputStyle_ = params->get("Output Style", outputStyle_default_);
630  } else {
631  outputStyle_ = (int)Teuchos::getParameter<Belos::OutputType>(*params,"Output Style");
632  }
633 
634  // Reconstruct the convergence test if the explicit residual test is not being used.
635  params_->set("Output Style", outputStyle_);
636  outputTest_ = Teuchos::null;
637  }
638 
639  // output stream
640  if (params->isParameter("Output Stream")) {
641  outputStream_ = Teuchos::getParameter<Teuchos::RCP<std::ostream> >(*params,"Output Stream");
642 
643  // Update parameter in our list.
644  params_->set("Output Stream", outputStream_);
645  if (printer_ != Teuchos::null)
646  printer_->setOStream( outputStream_ );
647  }
648 
649  // frequency level
650  if (verbosity_ & Belos::StatusTestDetails) {
651  if (params->isParameter("Output Frequency")) {
652  outputFreq_ = params->get("Output Frequency", outputFreq_default_);
653  }
654 
655  // Update parameter in out list and output status test.
656  params_->set("Output Frequency", outputFreq_);
657  if (outputTest_ != Teuchos::null)
658  outputTest_->setOutputFrequency( outputFreq_ );
659  }
660 
661  // Create output manager if we need to.
662  if (printer_ == Teuchos::null) {
663  printer_ = Teuchos::rcp( new OutputManager<ScalarType>(verbosity_, outputStream_) );
664  }
665 
666  // Convergence
667  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
668  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
669 
670  // Check for convergence tolerance
671  if (params->isParameter("Convergence Tolerance")) {
672  convtol_ = params->get("Convergence Tolerance",convtol_default_);
673 
674  // Update parameter in our list and residual tests.
675  params_->set("Convergence Tolerance", convtol_);
676  if (convTest_ != Teuchos::null)
677  convTest_->setTolerance( convtol_ );
678  }
679 
680  // Create status tests if we need to.
681 
682  // Basic test checks maximum iterations and native residual.
683  if (maxIterTest_ == Teuchos::null)
684  maxIterTest_ = Teuchos::rcp( new StatusTestMaxIters<ScalarType,MV,OP>( maxIters_ ) );
685 
686  if (convTest_ == Teuchos::null)
687  convTest_ = Teuchos::rcp( new StatusTestResNorm_t( convtol_, 1 ) );
688 
689  sTest_ = Teuchos::rcp( new StatusTestCombo_t( StatusTestCombo_t::OR, maxIterTest_, convTest_ ) );
690 
691  // Create the status test output class.
692  // This class manages and formats the output from the status test.
693  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory( outputStyle_ );
694  outputTest_ = stoFactory.create( printer_, sTest_, outputFreq_, Passed+Failed+Undefined );
695 
696  // Set the solver string for the output test
697  std::string solverDesc = " PCPG ";
698  outputTest_->setSolverDesc( solverDesc );
699 
700 
701  // Create orthogonalization manager if we need to.
702  if (ortho_ == Teuchos::null) {
703  params_->set("Orthogonalization", orthoType_);
704  if (orthoType_=="DGKS") {
705  if (orthoKappa_ <= 0) {
706  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
707  }
708  else {
709  ortho_ = Teuchos::rcp( new DGKSOrthoManager<ScalarType,MV,OP>( label_ ) );
710  Teuchos::rcp_dynamic_cast<DGKSOrthoManager<ScalarType,MV,OP> >(ortho_)->setDepTol( orthoKappa_ );
711  }
712  }
713  else if (orthoType_=="ICGS") {
714  ortho_ = Teuchos::rcp( new ICGSOrthoManager<ScalarType,MV,OP>( label_ ) );
715  }
716  else if (orthoType_=="IMGS") {
717  ortho_ = Teuchos::rcp( new IMGSOrthoManager<ScalarType,MV,OP>( label_ ) );
718  }
719  else {
720  TEUCHOS_TEST_FOR_EXCEPTION(orthoType_!="ICGS"&&orthoType_!="DGKS"&&orthoType_!="IMGS",std::logic_error,
721  "Belos::PCPGSolMgr(): Invalid orthogonalization type.");
722  }
723  }
724 
725  // Create the timer if we need to.
726  if (timerSolve_ == Teuchos::null) {
727  std::string solveLabel = label_ + ": PCPGSolMgr total solve time";
728 #ifdef BELOS_TEUCHOS_TIME_MONITOR
729  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
730 #endif
731  }
732 
733  // Inform the solver manager that the current parameters were set.
734  isSet_ = true;
735 }
736 
737 
738 template<class ScalarType, class MV, class OP>
739 Teuchos::RCP<const Teuchos::ParameterList>
741 {
742  static Teuchos::RCP<const Teuchos::ParameterList> validPL;
743  if (is_null(validPL)) {
744  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
745  // Set all the valid parameters and their default values.
746  pl->set("Convergence Tolerance", convtol_default_,
747  "The relative residual tolerance that needs to be achieved by the\n"
748  "iterative solver in order for the linear system to be declared converged.");
749  pl->set("Maximum Iterations", maxIters_default_,
750  "The maximum number of iterations allowed for each\n"
751  "set of RHS solved.");
752  pl->set("Num Deflated Blocks", deflatedBlocks_default_,
753  "The maximum number of vectors in the seed subspace." );
754  pl->set("Num Saved Blocks", savedBlocks_default_,
755  "The maximum number of vectors saved from old Krylov subspaces." );
756  pl->set("Verbosity", verbosity_default_,
757  "What type(s) of solver information should be outputted\n"
758  "to the output stream.");
759  pl->set("Output Style", outputStyle_default_,
760  "What style is used for the solver information outputted\n"
761  "to the output stream.");
762  pl->set("Output Frequency", outputFreq_default_,
763  "How often convergence information should be outputted\n"
764  "to the output stream.");
765  pl->set("Output Stream", outputStream_default_,
766  "A reference-counted pointer to the output stream where all\n"
767  "solver output is sent.");
768  pl->set("Timer Label", label_default_,
769  "The string to use as a prefix for the timer labels.");
770  pl->set("Orthogonalization", orthoType_default_,
771  "The type of orthogonalization to use: DGKS, ICGS, IMGS");
772  pl->set("Orthogonalization Constant",orthoKappa_default_,
773  "The constant used by DGKS orthogonalization to determine\n"
774  "whether another step of classical Gram-Schmidt is necessary.");
775  validPL = pl;
776  }
777  return validPL;
778 }
779 
780 
781 // solve()
782 template<class ScalarType, class MV, class OP>
784 
785  // Set the current parameters if are not set already.
786  if (!isSet_) { setParameters( params_ ); }
787 
788  Teuchos::BLAS<int,ScalarType> blas;
789  Teuchos::LAPACK<int,ScalarType> lapack;
790  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
791  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
792 
793  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,PCPGSolMgrLinearProblemFailure,
794  "Belos::PCPGSolMgr::solve(): Linear problem is not a valid object.");
795 
796  TEUCHOS_TEST_FOR_EXCEPTION(!problem_->isProblemSet(),PCPGSolMgrLinearProblemFailure,
797  "Belos::PCPGSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
798 
799  // Create indices for the linear systems to be solved.
800  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
801  std::vector<int> currIdx(1);
802  currIdx[0] = 0;
803 
804  bool debug = false;
805 
806  // Inform the linear problem of the current linear system to solve.
807  problem_->setLSIndex( currIdx ); // block size == 1
808 
809  // Assume convergence is achieved, then let any failed convergence set this to false.
810  bool isConverged = true;
811 
813  // PCPG iteration parameter list
814  Teuchos::ParameterList plist;
815  plist.set("Saved Blocks", savedBlocks_);
816  plist.set("Block Size", 1);
817  plist.set("Keep Diagonal", true);
818  plist.set("Initialize Diagonal", true);
819 
821  // PCPG solver
822 
823  Teuchos::RCP<PCPGIter<ScalarType,MV,OP> > pcpg_iter;
824  pcpg_iter = Teuchos::rcp( new PCPGIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
825  // Number of iterations required to generate initial recycle space (if needed)
826 
827  // Enter solve() iterations
828  {
829 #ifdef BELOS_TEUCHOS_TIME_MONITOR
830  Teuchos::TimeMonitor slvtimer(*timerSolve_);
831 #endif
832  while ( numRHS2Solve > 0 ) { // test for quick return
833 
834  // Reset the status test.
835  outputTest_->reset();
836 
837  // Create the first block in the current Krylov basis (residual).
838  if (R_ == Teuchos::null)
839  R_ = MVT::Clone( *(problem_->getRHS()), 1 );
840 
841  problem_->computeCurrResVec( &*R_ );
842 
843 
844  // Hypothesis: if U_ is not null, then neither is C_ and furthermore U'C= I.
845  // TODO: ensure hypothesis right here ... I have to think about use cases.
846 
847  if( U_ != Teuchos::null ){
848  // Hypothesis: if U_ is not null, then neither is C_ and furthermore U'C= I.
849 
850  // possibly over solved equation ... I want residual norms
851  // relative to the initial residual, not what I am about to compute.
852  Teuchos::RCP<MV> cur_soln_vec = problem_->getCurrLHSVec();
853  std::vector<MagnitudeType> rnorm0(1);
854  MVT::MvNorm( *R_, rnorm0 ); // rnorm0 = norm(R_);
855 
856  // Z := U_'*R_; xo += U_*Z ;R_ -= C_*Z
857  std::cout << "Solver Manager: dimU_ = " << dimU_ << std::endl;
858  Teuchos::SerialDenseMatrix<int,ScalarType> Z( dimU_, 1 );
859 
860  Teuchos::RCP<const MV> Uactive, Cactive;
861  std::vector<int> active_columns( dimU_ );
862  for (int i=0; i < dimU_; ++i) active_columns[i] = i;
863  Uactive = MVT::CloneView(*U_, active_columns);
864  Cactive = MVT::CloneView(*C_, active_columns);
865 
866  if( debug ){
867  std::cout << " Solver Manager : check duality of seed basis " << std::endl;
868  Teuchos::SerialDenseMatrix<int,ScalarType> H( dimU_, dimU_ );
869  MVT::MvTransMv( one, *Uactive, *Cactive, H );
870  H.print( std::cout );
871  }
872 
873  MVT::MvTransMv( one, *Uactive, *R_, Z );
874  Teuchos::RCP<MV> tempU = MVT::Clone( *R_, 1 );
875  MVT::MvTimesMatAddMv( one, *Uactive, Z, zero, *tempU ); // UZ
876  MVT::MvAddMv( one, *tempU, one, *cur_soln_vec, *cur_soln_vec ); // xo += tmp;
877  MVT::MvTimesMatAddMv( one, *Cactive, Z, zero, *tempU ); // CZ
878  MVT::MvAddMv( -one, *tempU, one, *R_, *R_ ); // R_ -= tmp;
879  std::vector<MagnitudeType> rnorm(1);
880  MVT::MvNorm( *R_, rnorm );
881  if( rnorm[0] < rnorm0[0] * .001 ){ //reorthogonalize
882  MVT::MvTransMv( one, *Uactive, *R_, Z );
883  MVT::MvTimesMatAddMv( one, *Uactive, Z, zero, *tempU );
884  MVT::MvAddMv( one, *tempU, one, *cur_soln_vec, *cur_soln_vec ); // xo += UZ;
885  MVT::MvTimesMatAddMv( one, *Cactive, Z, zero, *tempU );
886  MVT::MvAddMv( -one, *tempU, one, *R_, *R_ ); // R_ -= CZ;
887  }
888  Uactive = Teuchos::null;
889  Cactive = Teuchos::null;
890  tempU = Teuchos::null;
891  }
892  else {
893  dimU_ = 0;
894  }
895 
896 
897  // Set the new state and initialize the solver.
898  PCPGIterState<ScalarType,MV> pcpgState; // fails if R == null.
899 
900  pcpgState.R = R_;
901  if( U_ != Teuchos::null ) pcpgState.U = U_;
902  if( C_ != Teuchos::null ) pcpgState.C = C_;
903  if( dimU_ > 0 ) pcpgState.curDim = dimU_;
904  pcpg_iter->initialize(pcpgState);
905 
906  // treat initialize() exceptions here? how to use try-catch-throw? DMD
907 
908  // Get the current number of deflated blocks with the PCPG iteration
909  dimU_ = pcpgState.curDim;
910  if( !dimU_ ) printer_->stream(Debug) << " No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
911  pcpg_iter->resetNumIters();
912 
913  if( dimU_ > savedBlocks_ )
914  std::cout << "Error: dimU_ = " << dimU_ << " > savedBlocks_ = " << savedBlocks_ << std::endl;
915 
916  while(1) { // dummy loop for break
917 
918  // tell pcpg_iter to iterate
919  try {
920  if( debug ) printf("********** Calling iterate...\n");
921  pcpg_iter->iterate();
922 
924  //
925  // check convergence first
926  //
928  if ( convTest_->getStatus() == Passed ) {
929  // we have convergence
930  break; // break from while(1){pcpg_iter->iterate()}
931  }
933  //
934  // check for maximum iterations
935  //
937  else if ( maxIterTest_->getStatus() == Passed ) {
938  // we don't have convergence
939  isConverged = false;
940  break; // break from while(1){pcpg_iter->iterate()}
941  }
942  else {
943 
945  //
946  // we returned from iterate(), but none of our status tests Passed.
947  // Something is wrong, and it is probably the developers fault.
948  //
950 
951  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,
952  "Belos::PCPGSolMgr::solve(): Invalid return from PCPGIter::iterate().");
953  } // end if
954  } // end try
955  catch (const PCPGIterOrthoFailure &e) {
956 
957  // Check to see if the most recent solution yielded convergence.
958  sTest_->checkStatus( &*pcpg_iter );
959  if (convTest_->getStatus() != Passed)
960  isConverged = false;
961  break;
962  }
963  catch (const std::exception &e) {
964  printer_->stream(Errors) << "Error! Caught exception in PCPGIter::iterate() at iteration "
965  << pcpg_iter->getNumIters() << std::endl
966  << e.what() << std::endl;
967  throw;
968  }
969  } // end of while(1)
970 
971  // Update the linear problem.
972  Teuchos::RCP<MV> update = pcpg_iter->getCurrentUpdate();
973  problem_->updateSolution( update, true );
974 
975  // Inform the linear problem that we are finished with this block linear system.
976  problem_->setCurrLS();
977 
978  // Get the state. How did pcpgState die?
979  PCPGIterState<ScalarType,MV> oldState = pcpg_iter->getState();
980 
981  dimU_ = oldState.curDim;
982  int q = oldState.prevUdim;
983 
984  std::cout << "SolverManager: dimU_ " << dimU_ << " prevUdim= " << q << std::endl;
985 
986  if( q > deflatedBlocks_ )
987  std::cout << "SolverManager: Error deflatedBlocks = " << deflatedBlocks_ << std::endl;
988 
989  int rank;
990  if( dimU_ > q ){ // Orthogonalize [U;C](:,prevUdim:dimU_)
991  //Given the seed space U and C = A U for some symmetric positive definite A,
992  //find U1 and C1 with span(U1)=span(U) such that C1'U1 = I maintaining C=AU
993 
994  //oldState.D->print( std::cout ); D = diag( C'*U )
995 
996  U_ = oldState.U; //MVT::MvPrint( *U_, std::cout );
997  C_ = oldState.C; //MVT::MvPrint( *C_, std::cout );
998  rank = ARRQR(dimU_,q, *oldState.D );
999  if( rank < dimU_ ) {
1000  std::cout << " rank decreased in ARRQR, something to do? " << std::endl;
1001  }
1002  dimU_ = rank;
1003 
1004  } // Now U_ and C_ = AU are dual bases.
1005 
1006  if( dimU_ > deflatedBlocks_ ){
1007 
1008  if( !deflatedBlocks_ ){
1009  U_ = Teuchos::null;
1010  C_ = Teuchos::null;
1011  dimU_ = deflatedBlocks_;
1012  break;
1013  }
1014 
1015  bool Harmonic = false; // (Harmonic) Ritz vectors
1016 
1017  Teuchos::RCP<MV> Uorth;
1018 
1019  std::vector<int> active_cols( dimU_ );
1020  for (int i=0; i < dimU_; ++i) active_cols[i] = i;
1021 
1022  if( Harmonic ){
1023  Uorth = MVT::CloneCopy(*C_, active_cols);
1024  }
1025  else{
1026  Uorth = MVT::CloneCopy(*U_, active_cols);
1027  }
1028 
1029  // Explicitly construct Q and R factors
1030  Teuchos::SerialDenseMatrix<int,ScalarType> R(dimU_,dimU_);
1031  rank = ortho_->normalize(*Uorth, Teuchos::rcp(&R,false));
1032  Uorth = Teuchos::null;
1033  // TODO: During the previous solve, the matrix that normalizes U(1:q) was computed and discarded.
1034  // One might save it, reuse it here, and just normalize columns U(q+1:dimU_) here.
1035 
1036  // throw an error if U is both A-orthonormal and rank deficient
1037  TEUCHOS_TEST_FOR_EXCEPTION(rank < dimU_,PCPGSolMgrOrthoFailure,
1038  "Belos::PCPGSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
1039 
1040 
1041  // R VT' = Ur S,
1042  Teuchos::SerialDenseMatrix<int,ScalarType> VT; // Not referenced
1043  Teuchos::SerialDenseMatrix<int,ScalarType> Ur; // Not referenced
1044  int lwork = 5*dimU_; // minimal, extra computation < 67*dimU_
1045  int info = 0; // Hermite
1046  int lrwork = 1;
1047  if( problem_->isHermitian() ) lrwork = dimU_;
1048  std::vector<ScalarType> work(lwork); //
1049  std::vector<ScalarType> Svec(dimU_); //
1050  std::vector<ScalarType> rwork(lrwork);
1051  lapack.GESVD('N', 'O',
1052  R.numRows(),R.numCols(),R.values(), R.numRows(),
1053  &Svec[0],
1054  Ur.values(),1,
1055  VT.values(),1, // Output: VT stored in R
1056  &work[0], lwork,
1057  &rwork[0], &info);
1058 
1059  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, PCPGSolMgrLAPACKFailure,
1060  "Belos::PCPGSolMgr::solve(): LAPACK _GESVD failed to compute singular values.");
1061 
1062  if( work[0] != 67. * dimU_ )
1063  std::cout << " SVD " << dimU_ << " lwork " << work[0] << std::endl;
1064  for( int i=0; i< dimU_; i++)
1065  std::cout << i << " " << Svec[i] << std::endl;
1066 
1067  Teuchos::SerialDenseMatrix<int,ScalarType> wholeV( R, Teuchos::TRANS);
1068 
1069  int startRow = 0, startCol = 0;
1070  if( Harmonic )
1071  startCol = dimU_ - deflatedBlocks_;
1072 
1073  Teuchos::SerialDenseMatrix<int,ScalarType> V(Teuchos::Copy,
1074  wholeV,
1075  wholeV.numRows(),
1076  deflatedBlocks_,
1077  startRow,
1078  startCol);
1079  std::vector<int> active_columns( dimU_ );
1080  std::vector<int> def_cols( deflatedBlocks_ );
1081  for (int i=0; i < dimU_; ++i) active_columns[i] = i;
1082  for (int i=0; i < deflatedBlocks_; ++i) def_cols[i] = i;
1083 
1084  Teuchos::RCP<MV> Uactive = MVT::CloneViewNonConst(*U_, def_cols);
1085  Teuchos::RCP<MV> Ucopy = MVT::CloneCopy( *U_, active_columns );
1086  MVT::MvTimesMatAddMv( one, *Ucopy, V, zero, *Uactive ); // U:= U*V
1087  Ucopy = Teuchos::null;
1088  Uactive = Teuchos::null;
1089  Teuchos::RCP<MV> Cactive = MVT::CloneViewNonConst(*C_, def_cols);
1090  Teuchos::RCP<MV> Ccopy = MVT::CloneCopy( *C_, active_columns );
1091  MVT::MvTimesMatAddMv( one, *Ccopy, V, zero, *Cactive ); // C:= C*V
1092  Ccopy = Teuchos::null;
1093  Cactive = Teuchos::null;
1094  dimU_ = deflatedBlocks_;
1095  }
1096  printer_->stream(Debug) << " Generated recycled subspace using RHS index " << currIdx[0] << " of dimension " << dimU_ << std::endl << std::endl;
1097 
1098  // Inform the linear problem that we are finished with this block linear system.
1099  problem_->setCurrLS();
1100 
1101  // Update indices for the linear systems to be solved.
1102  numRHS2Solve -= 1;
1103  if ( numRHS2Solve > 0 ) {
1104  currIdx[0]++;
1105 
1106  // Set the next indices.
1107  problem_->setLSIndex( currIdx );
1108  }
1109  else {
1110  currIdx.resize( numRHS2Solve );
1111  }
1112  }// while ( numRHS2Solve > 0 )
1113  }
1114 
1115  // print final summary
1116  sTest_->print( printer_->stream(FinalSummary) );
1117 
1118  // print timing information
1119 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1120  // Calling summarize() can be expensive, so don't call unless the
1121  // user wants to print out timing details. summarize() will do all
1122  // the work even if it's passed a "black hole" output stream.
1123  if (verbosity_ & TimingDetails)
1124  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1125 #endif
1126 
1127  // Save the convergence test value ("achieved tolerance") for this solve.
1128  {
1129  using Teuchos::rcp_dynamic_cast;
1130  typedef StatusTestGenResNorm<ScalarType,MV,OP> conv_test_type;
1131  // testValues is nonnull and not persistent.
1132  const std::vector<MagnitudeType>* pTestValues =
1133  rcp_dynamic_cast<conv_test_type>(convTest_)->getTestValue();
1134 
1135  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1136  "Belos::PCPGSolMgr::solve(): The convergence test's getTestValue() "
1137  "method returned NULL. Please report this bug to the Belos developers.");
1138 
1139  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1140  "Belos::PCPGSolMgr::solve(): The convergence test's getTestValue() "
1141  "method returned a vector of length zero. Please report this bug to the "
1142  "Belos developers.");
1143 
1144  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1145  // achieved tolerances for all vectors in the current solve(), or
1146  // just for the vectors from the last deflation?
1147  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1148  }
1149 
1150  // get iteration information for this solve
1151  numIters_ = maxIterTest_->getNumIters();
1152 
1153  if (!isConverged) {
1154  return Unconverged; // return from PCPGSolMgr::solve()
1155  }
1156  return Converged; // return from PCPGSolMgr::solve()
1157 }
1158 
1159 // A-orthogonalize the Seed Space
1160 // Note that Anasazi::GenOrthoManager provides simplified versions of the algorithm,
1161 // that are not rank revealing, and are not designed for PCPG in other ways too.
1162 template<class ScalarType, class MV, class OP>
1163 int PCPGSolMgr<ScalarType,MV,OP,true>::ARRQR(int p, int q, const Teuchos::SerialDenseMatrix<int,ScalarType>& D)
1164 {
1165  using Teuchos::RCP;
1166  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1167  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1168 
1169  // Allocate memory for scalars.
1170  Teuchos::SerialDenseMatrix<int,ScalarType> alpha( 1, 1 );
1171  Teuchos::SerialDenseMatrix<int,ScalarType> gamma( 1, 1 );
1172  Teuchos::SerialDenseMatrix<int,ScalarType> anorm( 1, 1 );
1173  std::vector<int> curind(1);
1174  std::vector<int> ipiv(p - q); // RRQR Pivot indices
1175  std::vector<ScalarType> Pivots(p); // RRQR Pivots
1176  int i, imax, j, l;
1177  ScalarType rteps = 1.5e-8;
1178 
1179  // Scale such that diag( U'C) = I
1180  for( i = q ; i < p ; i++ ){
1181  ipiv[i-q] = i;
1182  curind[0] = i;
1183  RCP<MV> P = MVT::CloneViewNonConst(*U_,curind);
1184  RCP<MV> AP = MVT::CloneViewNonConst(*C_,curind);
1185  anorm(0,0) = one / Teuchos::ScalarTraits<ScalarType>::squareroot( D(i-q,i-q) ) ;
1186  MVT::MvAddMv( anorm(0,0), *P, zero, *AP, *P );
1187  MVT::MvAddMv( zero, *P, anorm(0,0), *AP, *AP );
1188  Pivots[i] = one;
1189  }
1190 
1191  for( i = q ; i < p ; i++ ){
1192  if( q < i && i < p-1 ){ // Find the largest pivot
1193  imax = i;
1194  l = ipiv[imax-q];
1195  for( j = i+1 ; j < p ; j++ ){
1196  const int k = ipiv[j-q];
1197  if( Pivots[k] > Pivots[l] ){
1198  imax = j;
1199  l = k;
1200  }
1201  } // end for
1202  if( imax > i ){
1203  l = ipiv[imax-q]; // swap ipiv( imax ) and ipiv(i+1)
1204  ipiv[imax-q] = ipiv[i-q];
1205  ipiv[i-q] = l;
1206  }
1207  } // largest pivot found
1208  int k = ipiv[i-q];
1209 
1210  if( Pivots[k] > 1.5625e-2 ){
1211  anorm(0,0) = Pivots[k]; // A-norm of u
1212  }
1213  else{ // anorm(0,0) = sqrt( U(:,k)'*C(:,k) );
1214  curind[0] = k;
1215  RCP<const MV> P = MVT::CloneView(*U_,curind);
1216  RCP<const MV> AP = MVT::CloneView(*C_,curind);
1217  MVT::MvTransMv( one, *P, *AP, anorm );
1218  anorm(0,0) = Teuchos::ScalarTraits<ScalarType>::squareroot( anorm(0,0) ) ;
1219  }
1220  if( rteps <= anorm(0,0) && anorm(0,0) < 9.765625e-4){
1221  /*
1222  C(:,k) = A*U(:,k); % Change C
1223  fixC = U(:, ipiv(1:i-1) )'*C(:,k);
1224  U(:,k) = U(:,k) - U(:, ipiv(1:i-1) )*fixC;
1225  C(:,k) = C(:,k) - C(:, ipiv(1:i-1) )*fixC;
1226  anorm = sqrt( U(:,k)'*C(:,k) );
1227  */
1228  std::cout << "ARRQR: Bad case not implemented" << std::endl;
1229  }
1230  if( anorm(0,0) < rteps ){ // rank [U;C] = i-1
1231  std::cout << "ARRQR : deficient case not implemented " << std::endl;
1232  //U = U(:, ipiv(1:i-1) );
1233  //C = C(:, ipiv(1:i-1) );
1234  p = q + i;
1235  // update curDim_ in State
1236  break;
1237  }
1238  curind[0] = k;
1239  RCP<MV> P = MVT::CloneViewNonConst(*U_,curind);
1240  RCP<MV> AP = MVT::CloneViewNonConst(*C_,curind);
1241  MVT::MvAddMv( anorm(0,0), *P, zero, *AP, *P ); // U(:,k) = U(:,k)/anorm;
1242  MVT::MvAddMv( zero, *P, anorm(0,0), *AP, *AP ); // C(:,k) = C(:,k)/anorm;
1243  P = Teuchos::null;
1244  AP = Teuchos::null;
1245  Pivots[k] = one; // delete, for diagonostic purposes
1246  P = MVT::CloneViewNonConst(*U_,curind); // U(:,k)
1247  AP = MVT::CloneViewNonConst(*C_,curind); // C(:,k)
1248  for( j = i+1 ; j < p ; j++ ){
1249  l = ipiv[j-q]; // ahhh
1250  curind[0] = l;
1251  RCP<MV> Q = MVT::CloneViewNonConst(*U_,curind); // segmentation fault, j=i+1=5
1252  MVT::MvTransMv( one, *Q, *AP, alpha); // alpha(0,0) = U(:,l)'*C(:,k);
1253  MVT::MvAddMv( -alpha(0,0), *P, one, *Q, *Q ); // U(:,l) -= U(:,k) * alpha(0,0);
1254  Q = Teuchos::null;
1255  RCP<MV> AQ = MVT::CloneViewNonConst(*C_,curind);
1256  MVT::MvAddMv( -alpha(0,0), *AP, one, *AQ, *AQ ); // C(:,l) -= C(:,l) - C(:,k) * alpha(0,0);
1257  AQ = Teuchos::null;
1258  gamma(0,0) = ( Pivots[l] - alpha(0,0))*( Pivots[l] + alpha(0,0));
1259  if( gamma(0,0) > 0){
1260  Pivots[l] = Teuchos::ScalarTraits<ScalarType>::squareroot( gamma(0,0) );
1261  }
1262  else {
1263  Pivots[l] = zero; //rank deficiency revealed
1264  }
1265  }
1266  }
1267  return p;
1268 }
1269 
1270 // The method returns a string describing the solver manager.
1271 template<class ScalarType, class MV, class OP>
1273 {
1274  std::ostringstream oss;
1275  oss << "Belos::PCPGSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
1276  oss << "{";
1277  oss << "Ortho Type='"<<orthoType_;
1278  oss << "}";
1279  return oss.str();
1280 }
1281 
1282 } // end Belos namespace
1283 
1284 #endif /* BELOS_PCPG_SOLMGR_HPP */
bool isLOADetected() const
Return whether a loss of accuracy was detected by this solver during the most current solve...
Collection of types and exceptions used within the Belos solvers.
int prevUdim
Number of block columns in matrices C and U before current iteration.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
PCPGSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
Class which manages the output and verbosity of the Belos solvers.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Belos concrete class to iterate Preconditioned Conjugate Projected Gradients.
PCPGSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonormal...
Teuchos::RCP< MV > R
The current residual.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
Teuchos::RCP< MV > C
C = AU, U spans recycled subspace.
Base class for Belos::SolverManager subclasses which normally can only compile with real ScalarType t...
PCPGSolMgrRecyclingFailure is thrown when any problem occurs in using/creating the recycling subspace...
A factory class for generating StatusTestOutput objects.
int curDim
The current dimension of the reduction.
An implementation of StatusTestResNorm using a family of residual norms.
PCPGSolMgrOrthoFailure(const std::string &what_arg)
Structure to contain pointers to PCPGIter state variables.
PCPGSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
Belos::StatusTest class for specifying a maximum number of iterations.
A factory class for generating StatusTestOutput objects.
int getNumIters() const
Get the iteration count for the most recent call to solve().
Iterated Modified Gram-Schmidt (IMGS) implementation of the Belos::OrthoManager class.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > D
The current Hessenberg matrix.
Classical Gram-Schmidt (with DGKS correction) implementation of the Belos::OrthoManager class...
A Belos::StatusTest class for specifying a maximum number of iterations.
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:205
Pure virtual base class which describes the basic interface for a solver manager. ...
void reset(const ResetType type)
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
A linear system to solve, and its associated information.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
Class which describes the linear problem to be solved by the iterative solver.
PCPGSolMgrLAPACKFailure is thrown when a nonzero value is retuned from an LAPACK call.
PCPGSolMgrLAPACKFailure(const std::string &what_arg)
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:154
MagnitudeType achievedTol() const
Tolerance achieved by the last solve() invocation.
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< MV > U
The recycled subspace.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using multiple steps ...
int setTolerance(MagnitudeType tolerance)
Set the value of the tolerance.
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.
PCPGIterOrthoFailure is thrown when the PCPGIter object is unable to compute independent direction ve...
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
PCPG iterative linear solver.
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
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.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
PCPGSolMgrRecyclingFailure(const std::string &what_arg)
This class implements the PCPG iteration, where a single-std::vector Krylov subspace is constructed...
Belos header file which uses auto-configuration information to include necessary C++ headers...
PCPGSolMgrLinearProblemFailure(const std::string &what_arg)

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