Belos  Version of the Day
BelosGCRODRSolMgr.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_GCRODR_SOLMGR_HPP
43 #define BELOS_GCRODR_SOLMGR_HPP
44 
48 
49 #include "BelosConfigDefs.hpp"
51 #include "BelosSolverManager.hpp"
52 #include "BelosLinearProblem.hpp"
53 #include "BelosTypes.hpp"
54 
55 #include "BelosGCRODRIter.hpp"
56 #include "BelosBlockFGmresIter.hpp"
59 #include "BelosStatusTestCombo.hpp"
61 #include "BelosOutputManager.hpp"
62 #include "Teuchos_BLAS.hpp" // includes Teuchos_ConfigDefs.hpp
63 #include "Teuchos_LAPACK.hpp"
64 #include "Teuchos_as.hpp"
65 #ifdef BELOS_TEUCHOS_TIME_MONITOR
66 # include "Teuchos_TimeMonitor.hpp"
67 #endif // BELOS_TEUCHOS_TIME_MONITOR
68 #if defined(HAVE_TEUCHOSCORE_CXX11)
69 # include <type_traits>
70 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
71 
79 namespace Belos {
80 
82 
83 
91  public:
92  GCRODRSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg) {}
93  };
94 
102  public:
103  GCRODRSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg) {}
104  };
105 
113  public:
114  GCRODRSolMgrLAPACKFailure(const std::string& what_arg) : BelosError(what_arg) {}
115  };
116 
124  public:
125  GCRODRSolMgrRecyclingFailure(const std::string& what_arg) : BelosError(what_arg) {}
126  };
127 
129 
154  template<class ScalarType, class MV, class OP,
155  const bool lapackSupportsScalarType =
157  class GCRODRSolMgr :
158  public Details::SolverManagerRequiresLapack<ScalarType,MV,OP>
159  {
160  static const bool requiresLapack =
162  typedef Details::SolverManagerRequiresLapack<ScalarType, MV, OP,
163  requiresLapack> base_type;
164 
165  public:
167  base_type ()
168  {}
169  GCRODRSolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >& problem,
170  const Teuchos::RCP<Teuchos::ParameterList>& pl) :
171  base_type ()
172  {}
173  virtual ~GCRODRSolMgr () {}
174  };
175 
179  template<class ScalarType, class MV, class OP>
180  class GCRODRSolMgr<ScalarType, MV, OP, true> :
181  public Details::SolverManagerRequiresLapack<ScalarType, MV, OP, true>
182  {
183 #if defined(HAVE_TEUCHOSCORE_CXX11)
184 # if defined(HAVE_TEUCHOS_COMPLEX)
185  static_assert (std::is_same<ScalarType, std::complex<float> >::value ||
186  std::is_same<ScalarType, std::complex<double> >::value ||
187  std::is_same<ScalarType, float>::value ||
188  std::is_same<ScalarType, double>::value,
189  "Belos::GCRODRSolMgr: ScalarType must be one of the four "
190  "types (S,D,C,Z) supported by LAPACK.");
191 # else
192  static_assert (std::is_same<ScalarType, float>::value ||
193  std::is_same<ScalarType, double>::value,
194  "Belos::GCRODRSolMgr: ScalarType must be float or double. "
195  "Complex arithmetic support is currently disabled. To "
196  "enable it, set Teuchos_ENABLE_COMPLEX=ON.");
197 # endif // defined(HAVE_TEUCHOS_COMPLEX)
198 #endif // defined(HAVE_TEUCHOSCORE_CXX11)
199 
200  private:
203  typedef Teuchos::ScalarTraits<ScalarType> SCT;
204  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
205  typedef Teuchos::ScalarTraits<MagnitudeType> MT;
207 
208  public:
210 
211 
217  GCRODRSolMgr();
218 
271  GCRODRSolMgr (const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem,
272  const Teuchos::RCP<Teuchos::ParameterList> &pl);
273 
275  virtual ~GCRODRSolMgr() {};
277 
279 
280 
284  return *problem_;
285  }
286 
289  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters() const;
290 
293  Teuchos::RCP<const Teuchos::ParameterList> getCurrentParameters() const {
294  return params_;
295  }
296 
302  Teuchos::Array<Teuchos::RCP<Teuchos::Time> > getTimers() const {
303  return Teuchos::tuple(timerSolve_);
304  }
305 
311  MagnitudeType achievedTol() const {
312  return achievedTol_;
313  }
314 
316  int getNumIters() const {
317  return numIters_;
318  }
319 
322  bool isLOADetected() const { return false; }
323 
325 
327 
328 
330  void setProblem( const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > &problem ) {
331  problem_ = problem;
332  }
333 
335  void setParameters( const Teuchos::RCP<Teuchos::ParameterList> &params );
336 
338 
340 
341 
345  void reset( const ResetType type ) {
346  if ((type & Belos::Problem) && !Teuchos::is_null(problem_)) {
347  bool set = problem_->setProblem();
348  if (!set)
349  throw "Could not set problem.";
350  }
351  else if (type & Belos::RecycleSubspace) {
352  keff = 0;
353  }
354  }
356 
358 
359 
386  ReturnType solve();
387 
389 
391 
393  std::string description() const;
394 
396 
397  private:
398 
399  // Called by all constructors; Contains init instructions common to all constructors
400  void init();
401 
402  // Initialize solver state storage
403  void initializeStateStorage();
404 
405  // Compute updated recycle space given existing recycle space and newly generated Krylov space
406  void buildRecycleSpace2(Teuchos::RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_iter);
407 
408  // Computes harmonic eigenpairs of projected matrix created during the priming solve.
409  // HH is the projected problem from the initial cycle of Gmres, it is (at least) of dimension m+1 x m.
410  // PP contains the harmonic eigenvectors corresponding to the recycledBlocks eigenvalues of smallest magnitude.
411  // The return value is the number of vectors needed to be stored, recycledBlocks or recycledBlocks+1.
412  int getHarmonicVecs1(int m,
413  const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
414  Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
415 
416  // Computes harmonic eigenpairs of projected matrix created during one cycle.
417  // HH is the total block projected problem from the GCRO-DR algorithm, it is (at least) of dimension keff+m+1 x keff+m.
418  // VV is the Krylov vectors from the projected GMRES algorithm, which has (at least) m+1 vectors.
419  // PP contains the harmonic eigenvectors corresponding to the recycledBlocks eigenvalues of smallest magnitude.
420  // The return value is the number of vectors needed to be stored, recycledBlocks or recycledBlocks+1.
421  int getHarmonicVecs2(int keff, int m,
422  const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
423  const Teuchos::RCP<const MV>& VV,
424  Teuchos::SerialDenseMatrix<int,ScalarType>& PP);
425 
426  // Sort list of n floating-point numbers and return permutation vector
427  void sort(std::vector<MagnitudeType>& dlist, int n, std::vector<int>& iperm);
428 
429  // Lapack interface
430  Teuchos::LAPACK<int,ScalarType> lapack;
431 
432  // Linear problem.
433  Teuchos::RCP<LinearProblem<ScalarType,MV,OP> > problem_;
434 
435  // Output manager.
436  Teuchos::RCP<OutputManager<ScalarType> > printer_;
437  Teuchos::RCP<std::ostream> outputStream_;
438 
439  // Status test.
440  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > sTest_;
441  Teuchos::RCP<StatusTestMaxIters<ScalarType,MV,OP> > maxIterTest_;
442  Teuchos::RCP<StatusTest<ScalarType,MV,OP> > convTest_;
443  Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
444  Teuchos::RCP<StatusTestOutput<ScalarType,MV,OP> > outputTest_;
445 
449  Teuchos::RCP<MatOrthoManager<ScalarType,MV,OP> > ortho_;
450 
451  // Current parameter list.
452  Teuchos::RCP<Teuchos::ParameterList> params_;
453 
454  // Default solver values.
455  static const MagnitudeType convTol_default_;
456  static const MagnitudeType orthoKappa_default_;
457  static const int maxRestarts_default_;
458  static const int maxIters_default_;
459  static const int numBlocks_default_;
460  static const int blockSize_default_;
461  static const int recycledBlocks_default_;
462  static const int verbosity_default_;
463  static const int outputStyle_default_;
464  static const int outputFreq_default_;
465  static const std::string impResScale_default_;
466  static const std::string expResScale_default_;
467  static const std::string label_default_;
468  static const std::string orthoType_default_;
469  static const Teuchos::RCP<std::ostream> outputStream_default_;
470 
471  // Current solver values.
472  MagnitudeType convTol_, orthoKappa_, achievedTol_;
473  int maxRestarts_, maxIters_, numIters_;
474  int verbosity_, outputStyle_, outputFreq_;
475  std::string orthoType_;
476  std::string impResScale_, expResScale_;
477 
479  // Solver State Storage
481  //
482  // The number of blocks and recycle blocks (m and k, respectively)
483  int numBlocks_, recycledBlocks_;
484  // Current size of recycled subspace
485  int keff;
486  //
487  // Residual vector
488  Teuchos::RCP<MV> r_;
489  //
490  // Search space
491  Teuchos::RCP<MV> V_;
492  //
493  // Recycled subspace and its image
494  Teuchos::RCP<MV> U_, C_;
495  //
496  // Updated recycle space and its image
497  Teuchos::RCP<MV> U1_, C1_;
498  //
499  // Storage used in constructing harmonic Ritz values/vectors
500  Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H2_;
501  Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H_;
502  Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > B_;
503  Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PP_;
504  Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > HP_;
505  std::vector<ScalarType> tau_;
506  std::vector<ScalarType> work_;
507  Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > R_;
508  std::vector<int> ipiv_;
510 
511  // Timers.
512  std::string label_;
513  Teuchos::RCP<Teuchos::Time> timerSolve_;
514 
515  // Internal state variables.
516  bool isSet_;
517 
518  // Have we generated or regenerated a recycle space yet this solve?
519  bool builtRecycleSpace_;
520  };
521 
522 
523 // Default solver values.
524 template<class ScalarType, class MV, class OP>
527 
528 template<class ScalarType, class MV, class OP>
531 
532 template<class ScalarType, class MV, class OP>
534 
535 template<class ScalarType, class MV, class OP>
537 
538 template<class ScalarType, class MV, class OP>
540 
541 template<class ScalarType, class MV, class OP>
543 
544 template<class ScalarType, class MV, class OP>
546 
547 template<class ScalarType, class MV, class OP>
549 
550 template<class ScalarType, class MV, class OP>
552 
553 template<class ScalarType, class MV, class OP>
555 
556 template<class ScalarType, class MV, class OP>
557 const std::string GCRODRSolMgr<ScalarType,MV,OP,true>::impResScale_default_ = "Norm of Preconditioned Initial Residual";
558 
559 template<class ScalarType, class MV, class OP>
560 const std::string GCRODRSolMgr<ScalarType,MV,OP,true>::expResScale_default_ = "Norm of Initial Residual";
561 
562 template<class ScalarType, class MV, class OP>
563 const std::string GCRODRSolMgr<ScalarType,MV,OP,true>::label_default_ = "Belos";
564 
565 template<class ScalarType, class MV, class OP>
567 
568 template<class ScalarType, class MV, class OP>
569 const Teuchos::RCP<std::ostream>
570 GCRODRSolMgr<ScalarType,MV,OP,true>::outputStream_default_ = Teuchos::rcpFromRef (std::cout);
571 
572 
573 // Empty Constructor
574 template<class ScalarType, class MV, class OP>
576  achievedTol_(0.0),
577  numIters_(0)
578 {
579  init ();
580 }
581 
582 
583 // Basic Constructor
584 template<class ScalarType, class MV, class OP>
586 GCRODRSolMgr(const Teuchos::RCP<LinearProblem<ScalarType,MV,OP> >& problem,
587  const Teuchos::RCP<Teuchos::ParameterList>& pl):
588  achievedTol_(0.0),
589  numIters_(0)
590 {
591  // Initialize local pointers to null, and initialize local variables
592  // to default values.
593  init ();
594 
595  TEUCHOS_TEST_FOR_EXCEPTION(
596  problem == Teuchos::null, std::invalid_argument,
597  "Belos::GCRODRSolMgr constructor: The solver manager's "
598  "constructor needs the linear problem argument 'problem' "
599  "to be non-null.");
600  problem_ = problem;
601 
602  // Set the parameters using the list that was passed in. If null,
603  // we defer initialization until a non-null list is set (by the
604  // client calling setParameters(), or by calling solve() -- in
605  // either case, a null parameter list indicates that default
606  // parameters should be used).
607  if (! pl.is_null ()) {
608  setParameters (pl);
609  }
610 }
611 
612 // Common instructions executed in all constructors
613 template<class ScalarType, class MV, class OP>
615  outputStream_ = outputStream_default_;
616  convTol_ = convTol_default_;
617  orthoKappa_ = orthoKappa_default_;
618  maxRestarts_ = maxRestarts_default_;
619  maxIters_ = maxIters_default_;
620  numBlocks_ = numBlocks_default_;
621  recycledBlocks_ = recycledBlocks_default_;
622  verbosity_ = verbosity_default_;
623  outputStyle_ = outputStyle_default_;
624  outputFreq_ = outputFreq_default_;
625  orthoType_ = orthoType_default_;
626  impResScale_ = impResScale_default_;
627  expResScale_ = expResScale_default_;
628  label_ = label_default_;
629  isSet_ = false;
630  builtRecycleSpace_ = false;
631  keff = 0;
632  r_ = Teuchos::null;
633  V_ = Teuchos::null;
634  U_ = Teuchos::null;
635  C_ = Teuchos::null;
636  U1_ = Teuchos::null;
637  C1_ = Teuchos::null;
638  PP_ = Teuchos::null;
639  HP_ = Teuchos::null;
640  H2_ = Teuchos::null;
641  R_ = Teuchos::null;
642  H_ = Teuchos::null;
643  B_ = Teuchos::null;
644 }
645 
646 template<class ScalarType, class MV, class OP>
647 void
649 setParameters (const Teuchos::RCP<Teuchos::ParameterList> &params)
650 {
651  using Teuchos::isParameterType;
652  using Teuchos::getParameter;
653  using Teuchos::null;
654  using Teuchos::ParameterList;
655  using Teuchos::parameterList;
656  using Teuchos::RCP;
657  using Teuchos::rcp;
658  using Teuchos::rcp_dynamic_cast;
659  using Teuchos::rcpFromRef;
660  using Teuchos::Exceptions::InvalidParameter;
661  using Teuchos::Exceptions::InvalidParameterName;
662  using Teuchos::Exceptions::InvalidParameterType;
663 
664  // The default parameter list contains all parameters that
665  // GCRODRSolMgr understands, and none that it doesn't understand.
666  RCP<const ParameterList> defaultParams = getValidParameters();
667 
668  // Create the internal parameter list if one doesn't already exist.
669  //
670  // (mfh 28 Feb 2011, 10 Mar 2011) At the time this code was written,
671  // ParameterList did not have validators or validateParameters().
672  // This is why the code below carefully validates the parameters one
673  // by one and fills in defaults. This code could be made a lot
674  // shorter by using validators. To do so, we would have to define
675  // appropriate validators for all the parameters. (This would more
676  // or less just move all that validation code out of this routine
677  // into to getValidParameters().)
678  //
679  // For an analogous reason, GCRODRSolMgr defines default parameter
680  // values as class data, as well as in the default ParameterList.
681  // This redundancy could be removed by defining the default
682  // parameter values only in the default ParameterList (which
683  // documents each parameter as well -- handy!).
684  if (params_.is_null()) {
685  params_ = parameterList (*defaultParams);
686  } else {
687  // A common case for setParameters() is for it to be called at the
688  // beginning of the solve() routine. This follows the Belos
689  // pattern of delaying initialization until the last possible
690  // moment (when the user asks Belos to perform the solve). In
691  // this common case, we save ourselves a deep copy of the input
692  // parameter list.
693  if (params_ != params) {
694  // Make a deep copy of the input parameter list. This allows
695  // the caller to modify or change params later, without
696  // affecting the behavior of this solver. This solver will then
697  // only change its internal parameters if setParameters() is
698  // called again.
699  params_ = parameterList (*params);
700  }
701 
702  // Fill in any missing parameters and their default values. Also,
703  // throw an exception if the parameter list has any misspelled or
704  // "extra" parameters. If you don't like this behavior, you'll
705  // want to replace the line of code below with your desired
706  // validation scheme. Note that Teuchos currently only implements
707  // two options:
708  //
709  // 1. validateParameters() requires that params_ has all the
710  // parameters that the default list has, and none that it
711  // doesn't have.
712  //
713  // 2. validateParametersAndSetDefaults() fills in missing
714  // parameters in params_ using the default list, but requires
715  // that any parameter provided in params_ is also in the
716  // default list.
717  //
718  // Here is an easy way to ignore any "extra" or misspelled
719  // parameters: Make a deep copy of the default list, fill in any
720  // "missing" parameters from the _input_ list, and then validate
721  // the input list using the deep copy of the default list. We
722  // show this in code:
723  //
724  // RCP<ParameterList> defaultCopy = parameterList (*getValidParameters ());
725  // defaultCopy->validateParametersAndSetDefaults (params);
726  // params->validateParametersAndSetDefaults (defaultCopy);
727  //
728  // This method is not entirely robust, because the input list may
729  // have incorrect validators set for existing parameters in the
730  // default list. This would then cause "validation" of the
731  // default list to throw an exception. As a result, we've chosen
732  // for now to be intolerant of misspellings and "extra" parameters
733  // in the input list.
734  params_->validateParametersAndSetDefaults (*defaultParams);
735  }
736 
737  // Check for maximum number of restarts.
738  if (params->isParameter ("Maximum Restarts")) {
739  maxRestarts_ = params->get("Maximum Restarts", maxRestarts_default_);
740 
741  // Update parameter in our list.
742  params_->set ("Maximum Restarts", maxRestarts_);
743  }
744 
745  // Check for maximum number of iterations
746  if (params->isParameter ("Maximum Iterations")) {
747  maxIters_ = params->get ("Maximum Iterations", maxIters_default_);
748 
749  // Update parameter in our list and in status test.
750  params_->set ("Maximum Iterations", maxIters_);
751  if (! maxIterTest_.is_null())
752  maxIterTest_->setMaxIters (maxIters_);
753  }
754 
755  // Check for the maximum number of blocks.
756  if (params->isParameter ("Num Blocks")) {
757  numBlocks_ = params->get ("Num Blocks", numBlocks_default_);
758  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks_ <= 0, std::invalid_argument,
759  "Belos::GCRODRSolMgr: The \"Num Blocks\" parameter must "
760  "be strictly positive, but you specified a value of "
761  << numBlocks_ << ".");
762  // Update parameter in our list.
763  params_->set ("Num Blocks", numBlocks_);
764  }
765 
766  // Check for the maximum number of blocks.
767  if (params->isParameter ("Num Recycled Blocks")) {
768  recycledBlocks_ = params->get ("Num Recycled Blocks",
769  recycledBlocks_default_);
770  TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks_ <= 0, std::invalid_argument,
771  "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
772  "parameter must be strictly positive, but you specified "
773  "a value of " << recycledBlocks_ << ".");
774  TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks_ >= numBlocks_, std::invalid_argument,
775  "Belos::GCRODRSolMgr: The \"Num Recycled Blocks\" "
776  "parameter must be less than the \"Num Blocks\" "
777  "parameter, but you specified \"Num Recycled Blocks\" "
778  "= " << recycledBlocks_ << " and \"Num Blocks\" = "
779  << numBlocks_ << ".");
780  // Update parameter in our list.
781  params_->set("Num Recycled Blocks", recycledBlocks_);
782  }
783 
784  // Check to see if the timer label changed. If it did, update it in
785  // the parameter list, and create a new timer with that label (if
786  // Belos was compiled with timers enabled).
787  if (params->isParameter ("Timer Label")) {
788  std::string tempLabel = params->get ("Timer Label", label_default_);
789 
790  // Update parameter in our list and solver timer
791  if (tempLabel != label_) {
792  label_ = tempLabel;
793  params_->set ("Timer Label", label_);
794  std::string solveLabel = label_ + ": GCRODRSolMgr total solve time";
795 #ifdef BELOS_TEUCHOS_TIME_MONITOR
796  timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
797 #endif
798  if (ortho_ != Teuchos::null) {
799  ortho_->setLabel( label_ );
800  }
801  }
802  }
803 
804  // Check for a change in verbosity level
805  if (params->isParameter ("Verbosity")) {
806  if (isParameterType<int> (*params, "Verbosity")) {
807  verbosity_ = params->get ("Verbosity", verbosity_default_);
808  } else {
809  verbosity_ = (int) getParameter<Belos::MsgType> (*params, "Verbosity");
810  }
811  // Update parameter in our list.
812  params_->set ("Verbosity", verbosity_);
813  // If the output manager (printer_) is null, then we will
814  // instantiate it later with the correct verbosity.
815  if (! printer_.is_null())
816  printer_->setVerbosity (verbosity_);
817  }
818 
819  // Check for a change in output style
820  if (params->isParameter ("Output Style")) {
821  if (isParameterType<int> (*params, "Output Style")) {
822  outputStyle_ = params->get ("Output Style", outputStyle_default_);
823  } else {
824  outputStyle_ = (int) getParameter<OutputType> (*params, "Output Style");
825  }
826 
827  // Update parameter in our list.
828  params_->set ("Output Style", outputStyle_);
829  // We will (re)instantiate the output status test afresh below.
830  outputTest_ = null;
831  }
832 
833  // Get the output stream for the output manager.
834  //
835  // While storing the output stream in the parameter list (either as
836  // an RCP or as a nonconst reference) is convenient and safe for
837  // programming, it makes it impossible to serialize the parameter
838  // list, read it back in from the serialized representation, and get
839  // the same output stream as before. This is because output streams
840  // may be arbitrary constructed objects.
841  //
842  // In case the user tried reading in the parameter list from a
843  // serialized representation and the output stream can't be read
844  // back in, we set the output stream to point to std::cout. This
845  // ensures reasonable behavior.
846  if (params->isParameter ("Output Stream")) {
847  try {
848  outputStream_ = getParameter<RCP<std::ostream> > (*params, "Output Stream");
849  } catch (InvalidParameter&) {
850  outputStream_ = rcpFromRef (std::cout);
851  }
852  // We assume that a null output stream indicates that the user
853  // doesn't want to print anything, so we replace it with a "black
854  // hole" stream that prints nothing sent to it. (We can't use a
855  // null output stream, since the output manager always sends
856  // things it wants to print to the output stream.)
857  if (outputStream_.is_null()) {
858  outputStream_ = rcp (new Teuchos::oblackholestream);
859  }
860  // Update parameter in our list.
861  params_->set ("Output Stream", outputStream_);
862  // If the output manager (printer_) is null, then we will
863  // instantiate it later with the correct output stream.
864  if (! printer_.is_null()) {
865  printer_->setOStream (outputStream_);
866  }
867  }
868 
869  // frequency level
870  if (verbosity_ & Belos::StatusTestDetails) {
871  if (params->isParameter ("Output Frequency")) {
872  outputFreq_ = params->get ("Output Frequency", outputFreq_default_);
873  }
874 
875  // Update parameter in out list and output status test.
876  params_->set("Output Frequency", outputFreq_);
877  if (! outputTest_.is_null())
878  outputTest_->setOutputFrequency (outputFreq_);
879  }
880 
881  // Create output manager if we need to, using the verbosity level
882  // and output stream that we fetched above. We do this here because
883  // instantiating an OrthoManager using OrthoManagerFactory requires
884  // a valid OutputManager.
885  if (printer_.is_null()) {
886  printer_ = rcp (new OutputManager<ScalarType> (verbosity_, outputStream_));
887  }
888 
889  // Get the orthogonalization manager name ("Orthogonalization").
890  //
891  // Getting default values for the orthogonalization manager
892  // parameters ("Orthogonalization Parameters") requires knowing the
893  // orthogonalization manager name. Save it for later, and also
894  // record whether it's different than before.
896  bool changedOrthoType = false;
897  if (params->isParameter ("Orthogonalization")) {
898  const std::string& tempOrthoType =
899  params->get ("Orthogonalization", orthoType_default_);
900  // Ensure that the specified orthogonalization type is valid.
901  if (! factory.isValidName (tempOrthoType)) {
902  std::ostringstream os;
903  os << "Belos::GCRODRSolMgr: Invalid orthogonalization name \""
904  << tempOrthoType << "\". The following are valid options "
905  << "for the \"Orthogonalization\" name parameter: ";
906  factory.printValidNames (os);
907  throw std::invalid_argument (os.str());
908  }
909  if (tempOrthoType != orthoType_) {
910  changedOrthoType = true;
911  orthoType_ = tempOrthoType;
912  // Update parameter in our list.
913  params_->set ("Orthogonalization", orthoType_);
914  }
915  }
916 
917  // Get any parameters for the orthogonalization ("Orthogonalization
918  // Parameters"). If not supplied, the orthogonalization manager
919  // factory will supply default values.
920  //
921  // NOTE (mfh 12 Jan 2011) For the sake of backwards compatibility,
922  // if params has an "Orthogonalization Constant" parameter and the
923  // DGKS orthogonalization manager is to be used, the value of this
924  // parameter will override DGKS's "depTol" parameter.
925  //
926  // Users must supply the orthogonalization manager parameters as a
927  // sublist (supplying it as an RCP<ParameterList> would make the
928  // resulting parameter list not serializable).
929  RCP<ParameterList> orthoParams;
930  { // The nonmember function sublist() returns an RCP<ParameterList>,
931  // which is what we want here.
932  using Teuchos::sublist;
933  // Abbreviation to avoid typos.
934  const std::string paramName ("Orthogonalization Parameters");
935 
936  try {
937  orthoParams = sublist (params_, paramName, true);
938  } catch (InvalidParameter&) {
939  // We didn't get the parameter list from params, so get a
940  // default parameter list from the OrthoManagerFactory. Modify
941  // params_ so that it has the default parameter list, and set
942  // orthoParams to ensure it's a sublist of params_ (and not just
943  // a copy of one).
944  params_->set (paramName, factory.getDefaultParameters (orthoType_));
945  orthoParams = sublist (params_, paramName, true);
946  }
947  }
948  TEUCHOS_TEST_FOR_EXCEPTION(orthoParams.is_null(), std::logic_error,
949  "Failed to get orthogonalization parameters. "
950  "Please report this bug to the Belos developers.");
951 
952  // Instantiate a new MatOrthoManager subclass instance if necessary.
953  // If not necessary, then tell the existing instance about the new
954  // parameters.
955  if (ortho_.is_null() || changedOrthoType) {
956  // We definitely need to make a new MatOrthoManager, since either
957  // we haven't made one yet, or we've changed orthogonalization
958  // methods. Creating the orthogonalization manager requires that
959  // the OutputManager (printer_) already be initialized.
960  ortho_ = factory.makeMatOrthoManager (orthoType_, null, printer_,
961  label_, orthoParams);
962  } else {
963  // If the MatOrthoManager implements the ParameterListAcceptor
964  // mix-in interface, we can propagate changes to its parameters
965  // without reinstantiating the MatOrthoManager.
966  //
967  // We recommend that all MatOrthoManager subclasses implement
968  // Teuchos::ParameterListAcceptor, but do not (yet) require this.
969  typedef Teuchos::ParameterListAcceptor PLA;
970  RCP<PLA> pla = rcp_dynamic_cast<PLA> (ortho_);
971  if (pla.is_null()) {
972  // Oops, it's not a ParameterListAcceptor. We have to
973  // reinstantiate the MatOrthoManager in order to pass in the
974  // possibly new parameters.
975  ortho_ = factory.makeMatOrthoManager (orthoType_, null, printer_,
976  label_, orthoParams);
977  } else {
978  pla->setParameterList (orthoParams);
979  }
980  }
981 
982  // The DGKS orthogonalization accepts a "Orthogonalization Constant"
983  // parameter (also called kappa in the code, but not in the
984  // parameter list). If its value is provided in the given parameter
985  // list, and its value is positive, use it. Ignore negative values.
986  //
987  // NOTE (mfh 12 Jan 2011) This overrides the "depTol" parameter that
988  // may have been specified in "Orthogonalization Parameters". We
989  // retain this behavior for backwards compatibility.
990  if (params->isParameter ("Orthogonalization Constant")) {
991  const MagnitudeType orthoKappa =
992  params->get ("Orthogonalization Constant", orthoKappa_default_);
993  if (orthoKappa > 0) {
994  orthoKappa_ = orthoKappa;
995  // Update parameter in our list.
996  params_->set("Orthogonalization Constant", orthoKappa_);
997  // Only DGKS currently accepts this parameter.
998  if (orthoType_ == "DGKS" && ! ortho_.is_null()) {
999  typedef DGKSOrthoManager<ScalarType, MV, OP> ortho_man_type;
1000  // This cast should always succeed; it's a bug
1001  // otherwise. (If the cast fails, then orthoType_
1002  // doesn't correspond to the OrthoManager subclass
1003  // instance that we think we have, so we initialized the
1004  // wrong subclass somehow.)
1005  rcp_dynamic_cast<ortho_man_type>(ortho_)->setDepTol (orthoKappa_);
1006  }
1007  }
1008  }
1009 
1010  // Convergence
1011  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
1012  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
1013 
1014  // Check for convergence tolerance
1015  if (params->isParameter("Convergence Tolerance")) {
1016  convTol_ = params->get ("Convergence Tolerance", convTol_default_);
1017 
1018  // Update parameter in our list and residual tests.
1019  params_->set ("Convergence Tolerance", convTol_);
1020  if (! impConvTest_.is_null())
1021  impConvTest_->setTolerance (convTol_);
1022  if (! expConvTest_.is_null())
1023  expConvTest_->setTolerance (convTol_);
1024  }
1025 
1026  // Check for a change in scaling, if so we need to build new residual tests.
1027  if (params->isParameter ("Implicit Residual Scaling")) {
1028  std::string tempImpResScale =
1029  getParameter<std::string> (*params, "Implicit Residual Scaling");
1030 
1031  // Only update the scaling if it's different.
1032  if (impResScale_ != tempImpResScale) {
1033  ScaleType impResScaleType = convertStringToScaleType (tempImpResScale);
1034  impResScale_ = tempImpResScale;
1035 
1036  // Update parameter in our list and residual tests
1037  params_->set("Implicit Residual Scaling", impResScale_);
1038  // NOTE (mfh 28 Feb 2011) StatusTestImpResNorm only lets you
1039  // call defineScaleForm() once. The code below attempts to call
1040  // defineScaleForm(); if the scale form has already been
1041  // defined, it constructs a new StatusTestImpResNorm instance.
1042  // StatusTestImpResNorm really should not expose the
1043  // defineScaleForm() method, since it's serving an
1044  // initialization purpose; all initialization should happen in
1045  // the constructor whenever possible. In that case, the code
1046  // below could be simplified into a single (re)instantiation.
1047  if (! impConvTest_.is_null()) {
1048  try {
1049  impConvTest_->defineScaleForm (impResScaleType, Belos::TwoNorm);
1050  }
1051  catch (StatusTestError&) {
1052  // Delete the convergence test so it gets constructed again.
1053  impConvTest_ = null;
1054  convTest_ = null;
1055  }
1056  }
1057  }
1058  }
1059 
1060  if (params->isParameter("Explicit Residual Scaling")) {
1061  std::string tempExpResScale =
1062  getParameter<std::string> (*params, "Explicit Residual Scaling");
1063 
1064  // Only update the scaling if it's different.
1065  if (expResScale_ != tempExpResScale) {
1066  ScaleType expResScaleType = convertStringToScaleType (tempExpResScale);
1067  expResScale_ = tempExpResScale;
1068 
1069  // Update parameter in our list and residual tests
1070  params_->set("Explicit Residual Scaling", expResScale_);
1071  // NOTE (mfh 28 Feb 2011) See note above on the (re)construction
1072  // of StatusTestImpResNorm.
1073  if (! expConvTest_.is_null()) {
1074  try {
1075  expConvTest_->defineScaleForm (expResScaleType, Belos::TwoNorm);
1076  }
1077  catch (StatusTestError&) {
1078  // Delete the convergence test so it gets constructed again.
1079  expConvTest_ = null;
1080  convTest_ = null;
1081  }
1082  }
1083  }
1084  }
1085  //
1086  // Create iteration stopping criteria ("status tests") if we need
1087  // to, by combining three different stopping criteria.
1088  //
1089  // First, construct maximum-number-of-iterations stopping criterion.
1090  if (maxIterTest_.is_null())
1091  maxIterTest_ = rcp (new StatusTestMaxIters<ScalarType,MV,OP> (maxIters_));
1092 
1093  // Implicit residual test, using the native residual to determine if
1094  // convergence was achieved.
1095  if (impConvTest_.is_null()) {
1096  impConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
1097  impConvTest_->defineScaleForm (convertStringToScaleType (impResScale_),
1098  Belos::TwoNorm);
1099  }
1100 
1101  // Explicit residual test once the native residual is below the tolerance
1102  if (expConvTest_.is_null()) {
1103  expConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
1104  expConvTest_->defineResForm (StatusTestResNorm_t::Explicit, Belos::TwoNorm);
1105  expConvTest_->defineScaleForm (convertStringToScaleType (expResScale_),
1106  Belos::TwoNorm);
1107  }
1108  // Convergence test first tests the implicit residual, then the
1109  // explicit residual if the implicit residual test passes.
1110  if (convTest_.is_null()) {
1111  convTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::SEQ,
1112  impConvTest_,
1113  expConvTest_));
1114  }
1115  // Construct the complete iteration stopping criterion:
1116  //
1117  // "Stop iterating if the maximum number of iterations has been
1118  // reached, or if the convergence test passes."
1119  sTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::OR,
1120  maxIterTest_,
1121  convTest_));
1122  // Create the status test output class.
1123  // This class manages and formats the output from the status test.
1124  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory (outputStyle_);
1125  outputTest_ = stoFactory.create (printer_, sTest_, outputFreq_,
1127 
1128  // Set the solver string for the output test
1129  std::string solverDesc = " GCRODR ";
1130  outputTest_->setSolverDesc( solverDesc );
1131 
1132  // Create the timer if we need to.
1133  if (timerSolve_.is_null()) {
1134  std::string solveLabel = label_ + ": GCRODRSolMgr total solve time";
1135 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1136  timerSolve_ = Teuchos::TimeMonitor::getNewCounter(solveLabel);
1137 #endif
1138  }
1139 
1140  // Inform the solver manager that the current parameters were set.
1141  isSet_ = true;
1142 }
1143 
1144 
1145 template<class ScalarType, class MV, class OP>
1146 Teuchos::RCP<const Teuchos::ParameterList>
1148 {
1149  using Teuchos::ParameterList;
1150  using Teuchos::parameterList;
1151  using Teuchos::RCP;
1152 
1153  static RCP<const ParameterList> validPL;
1154  if (is_null(validPL)) {
1155  RCP<ParameterList> pl = parameterList ();
1156 
1157  // Set all the valid parameters and their default values.
1158  pl->set("Convergence Tolerance", convTol_default_,
1159  "The relative residual tolerance that needs to be achieved by the\n"
1160  "iterative solver in order for the linear system to be declared converged.");
1161  pl->set("Maximum Restarts", maxRestarts_default_,
1162  "The maximum number of cycles allowed for each\n"
1163  "set of RHS solved.");
1164  pl->set("Maximum Iterations", maxIters_default_,
1165  "The maximum number of iterations allowed for each\n"
1166  "set of RHS solved.");
1167  // mfh 25 Oct 2010: "Block Size" must be 1 because GCRODR is
1168  // currently not a block method: i.e., it does not work on
1169  // multiple right-hand sides at once.
1170  pl->set("Block Size", blockSize_default_,
1171  "Block Size Parameter -- currently must be 1 for GCRODR");
1172  pl->set("Num Blocks", numBlocks_default_,
1173  "The maximum number of vectors allowed in the Krylov subspace\n"
1174  "for each set of RHS solved.");
1175  pl->set("Num Recycled Blocks", recycledBlocks_default_,
1176  "The maximum number of vectors in the recycled subspace." );
1177  pl->set("Verbosity", verbosity_default_,
1178  "What type(s) of solver information should be outputted\n"
1179  "to the output stream.");
1180  pl->set("Output Style", outputStyle_default_,
1181  "What style is used for the solver information outputted\n"
1182  "to the output stream.");
1183  pl->set("Output Frequency", outputFreq_default_,
1184  "How often convergence information should be outputted\n"
1185  "to the output stream.");
1186  pl->set("Output Stream", outputStream_default_,
1187  "A reference-counted pointer to the output stream where all\n"
1188  "solver output is sent.");
1189  pl->set("Implicit Residual Scaling", impResScale_default_,
1190  "The type of scaling used in the implicit residual convergence test.");
1191  pl->set("Explicit Residual Scaling", expResScale_default_,
1192  "The type of scaling used in the explicit residual convergence test.");
1193  pl->set("Timer Label", label_default_,
1194  "The string to use as a prefix for the timer labels.");
1195  {
1197  pl->set("Orthogonalization", orthoType_default_,
1198  "The type of orthogonalization to use. Valid options: " +
1199  factory.validNamesString());
1200  RCP<const ParameterList> orthoParams =
1201  factory.getDefaultParameters (orthoType_default_);
1202  pl->set ("Orthogonalization Parameters", *orthoParams,
1203  "Parameters specific to the type of orthogonalization used.");
1204  }
1205  pl->set("Orthogonalization Constant", orthoKappa_default_,
1206  "When using DGKS orthogonalization: the \"depTol\" constant, used "
1207  "to determine whether another step of classical Gram-Schmidt is "
1208  "necessary. Otherwise ignored.");
1209  validPL = pl;
1210  }
1211  return validPL;
1212 }
1213 
1214 // initializeStateStorage
1215 template<class ScalarType, class MV, class OP>
1217 
1218  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1219 
1220  // Check if there is any multivector to clone from.
1221  Teuchos::RCP<const MV> rhsMV = problem_->getRHS();
1222  if (rhsMV == Teuchos::null) {
1223  // Nothing to do
1224  return;
1225  }
1226  else {
1227 
1228  // Initialize the state storage
1229  TEUCHOS_TEST_FOR_EXCEPTION(static_cast<ptrdiff_t>(numBlocks_) > MVT::GetGlobalLength(*rhsMV),std::invalid_argument,
1230  "Belos::GCRODRSolMgr::initializeStateStorage(): Cannot generate a Krylov basis with dimension larger the operator!");
1231 
1232  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1233  if (U_ == Teuchos::null) {
1234  U_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1235  }
1236  else {
1237  // Generate U_ by cloning itself ONLY if more space is needed.
1238  if (MVT::GetNumberVecs(*U_) < recycledBlocks_+1) {
1239  Teuchos::RCP<const MV> tmp = U_;
1240  U_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1241  }
1242  }
1243 
1244  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1245  if (C_ == Teuchos::null) {
1246  C_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1247  }
1248  else {
1249  // Generate C_ by cloning itself ONLY if more space is needed.
1250  if (MVT::GetNumberVecs(*C_) < recycledBlocks_+1) {
1251  Teuchos::RCP<const MV> tmp = C_;
1252  C_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1253  }
1254  }
1255 
1256  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1257  if (V_ == Teuchos::null) {
1258  V_ = MVT::Clone( *rhsMV, numBlocks_+1 );
1259  }
1260  else {
1261  // Generate V_ by cloning itself ONLY if more space is needed.
1262  if (MVT::GetNumberVecs(*V_) < numBlocks_+1) {
1263  Teuchos::RCP<const MV> tmp = V_;
1264  V_ = MVT::Clone( *tmp, numBlocks_+1 );
1265  }
1266  }
1267 
1268  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1269  if (U1_ == Teuchos::null) {
1270  U1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1271  }
1272  else {
1273  // Generate U1_ by cloning itself ONLY if more space is needed.
1274  if (MVT::GetNumberVecs(*U1_) < recycledBlocks_+1) {
1275  Teuchos::RCP<const MV> tmp = U1_;
1276  U1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1277  }
1278  }
1279 
1280  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1281  if (C1_ == Teuchos::null) {
1282  C1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1283  }
1284  else {
1285  // Generate C1_ by cloning itself ONLY if more space is needed.
1286  if (MVT::GetNumberVecs(*C1_) < recycledBlocks_+1) {
1287  Teuchos::RCP<const MV> tmp = C1_;
1288  C1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1289  }
1290  }
1291 
1292  // Generate r_ only if it doesn't exist
1293  if (r_ == Teuchos::null)
1294  r_ = MVT::Clone( *rhsMV, 1 );
1295 
1296  // Size of tau_ will change during computation, so just be sure it starts with appropriate size
1297  tau_.resize(recycledBlocks_+1);
1298 
1299  // Size of work_ will change during computation, so just be sure it starts with appropriate size
1300  work_.resize(recycledBlocks_+1);
1301 
1302  // Size of ipiv_ will change during computation, so just be sure it starts with appropriate size
1303  ipiv_.resize(recycledBlocks_+1);
1304 
1305  // Generate H2_ only if it doesn't exist, otherwise resize it.
1306  if (H2_ == Teuchos::null)
1307  H2_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 ) );
1308  else {
1309  if ( (H2_->numRows() != numBlocks_+recycledBlocks_+2) || (H2_->numCols() != numBlocks_+recycledBlocks_+1) )
1310  H2_->reshape( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 );
1311  }
1312  H2_->putScalar(zero);
1313 
1314  // Generate R_ only if it doesn't exist, otherwise resize it.
1315  if (R_ == Teuchos::null)
1316  R_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( recycledBlocks_+1, recycledBlocks_+1 ) );
1317  else {
1318  if ( (R_->numRows() != recycledBlocks_+1) || (R_->numCols() != recycledBlocks_+1) )
1319  R_->reshape( recycledBlocks_+1, recycledBlocks_+1 );
1320  }
1321  R_->putScalar(zero);
1322 
1323  // Generate PP_ only if it doesn't exist, otherwise resize it.
1324  if (PP_ == Teuchos::null)
1325  PP_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycledBlocks_+2, recycledBlocks_+1 ) );
1326  else {
1327  if ( (PP_->numRows() != numBlocks_+recycledBlocks_+2) || (PP_->numCols() != recycledBlocks_+1) )
1328  PP_->reshape( numBlocks_+recycledBlocks_+2, recycledBlocks_+1 );
1329  }
1330 
1331  // Generate HP_ only if it doesn't exist, otherwise resize it.
1332  if (HP_ == Teuchos::null)
1333  HP_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 ) );
1334  else {
1335  if ( (HP_->numRows() != numBlocks_+recycledBlocks_+2) || (HP_->numCols() != numBlocks_+recycledBlocks_+1) )
1336  HP_->reshape( numBlocks_+recycledBlocks_+2, numBlocks_+recycledBlocks_+1 );
1337  }
1338 
1339  } // end else
1340 }
1341 
1342 
1343 // solve()
1344 template<class ScalarType, class MV, class OP>
1346  using Teuchos::RCP;
1347  using Teuchos::rcp;
1348 
1349  // Set the current parameters if they were not set before.
1350  // NOTE: This may occur if the user generated the solver manager with the default constructor and
1351  // then didn't set any parameters using setParameters().
1352  if (!isSet_) { setParameters( params_ ); }
1353 
1354  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1355  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1356  std::vector<int> index(numBlocks_+1);
1357 
1358  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,GCRODRSolMgrLinearProblemFailure, "Belos::GCRODRSolMgr::solve(): Linear problem is not a valid object.");
1359 
1360  TEUCHOS_TEST_FOR_EXCEPTION(!problem_->isProblemSet(),GCRODRSolMgrLinearProblemFailure,"Belos::GCRODRSolMgr::solve(): Linear problem is not ready, setProblem() has not been called.");
1361 
1362  // Create indices for the linear systems to be solved.
1363  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
1364  std::vector<int> currIdx(1);
1365  currIdx[0] = 0;
1366 
1367  // Inform the linear problem of the current linear system to solve.
1368  problem_->setLSIndex( currIdx );
1369 
1370  // Check the number of blocks and change them is necessary.
1371  ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
1372  if (static_cast<ptrdiff_t>(numBlocks_) > dim) {
1373  numBlocks_ = Teuchos::as<int>(dim);
1374  printer_->stream(Warnings) <<
1375  "Warning! Requested Krylov subspace dimension is larger than operator dimension!" << std::endl <<
1376  " The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << numBlocks_ << std::endl;
1377  params_->set("Num Blocks", numBlocks_);
1378  }
1379 
1380  // Assume convergence is achieved, then let any failed convergence set this to false.
1381  bool isConverged = true;
1382 
1383  // Initialize storage for all state variables
1384  initializeStateStorage();
1385 
1387  // Parameter list
1388  Teuchos::ParameterList plist;
1389 
1390  plist.set("Num Blocks",numBlocks_);
1391  plist.set("Recycled Blocks",recycledBlocks_);
1392 
1394  // GCRODR solver
1395 
1396  RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_iter;
1397  gcrodr_iter = rcp( new GCRODRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,plist) );
1398  // Number of iterations required to generate initial recycle space (if needed)
1399  int prime_iterations = 0;
1400 
1401  // Enter solve() iterations
1402  {
1403 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1404  Teuchos::TimeMonitor slvtimer(*timerSolve_);
1405 #endif
1406 
1407  while ( numRHS2Solve > 0 ) {
1408 
1409  // Set flag indicating recycle space has not been generated this solve
1410  builtRecycleSpace_ = false;
1411 
1412  // Reset the status test.
1413  outputTest_->reset();
1414 
1416  // Initialize recycled subspace for GCRODR
1417 
1418  // If there is a subspace to recycle, recycle it, otherwise generate the initial recycled subspace.
1419  if (keff > 0) {
1420  TEUCHOS_TEST_FOR_EXCEPTION(keff < recycledBlocks_,GCRODRSolMgrRecyclingFailure,
1421  "Belos::GCRODRSolMgr::solve(): Requested size of recycled subspace is not consistent with the current recycle subspace.");
1422 
1423  printer_->stream(Debug) << " Now solving RHS index " << currIdx[0] << " using recycled subspace of dimension " << keff << std::endl << std::endl;
1424  // Compute image of U_ under the new operator
1425  index.resize(keff);
1426  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1427  RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1428  RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1429  problem_->apply( *Utmp, *Ctmp );
1430 
1431  RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1432 
1433  // Orthogonalize this block
1434  // Get a matrix to hold the orthonormalization coefficients.
1435  Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *R_, keff, keff );
1436  int rank = ortho_->normalize(*Ctmp, rcp(&Rtmp,false));
1437  // Throw an error if we could not orthogonalize this block
1438  TEUCHOS_TEST_FOR_EXCEPTION(rank != keff,GCRODRSolMgrOrthoFailure,"Belos::GCRODRSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
1439 
1440  // U_ = U_*R^{-1}
1441  // First, compute LU factorization of R
1442  int info = 0;
1443  ipiv_.resize(Rtmp.numRows());
1444  lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
1445  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1446 
1447  // Now, form inv(R)
1448  int lwork = Rtmp.numRows();
1449  work_.resize(lwork);
1450  lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1451  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to invert triangular matrix.");
1452 
1453  // U_ = U1_; (via a swap)
1454  MVT::MvTimesMatAddMv( one, *Utmp, Rtmp, zero, *U1tmp );
1455  std::swap(U_, U1_);
1456 
1457  // Must reinitialize after swap
1458  index.resize(keff);
1459  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1460  Ctmp = MVT::CloneViewNonConst( *C_, index );
1461  Utmp = MVT::CloneView( *U_, index );
1462 
1463  // Compute C_'*r_
1464  Teuchos::SerialDenseMatrix<int,ScalarType> Ctr(keff,1);
1465  problem_->computeCurrPrecResVec( &*r_ );
1466  MVT::MvTransMv( one, *Ctmp, *r_, Ctr );
1467 
1468  // Update solution ( x += U_*C_'*r_ )
1469  RCP<MV> update = MVT::Clone( *problem_->getCurrLHSVec(), 1 );
1470  MVT::MvInit( *update, 0.0 );
1471  MVT::MvTimesMatAddMv( one, *Utmp, Ctr, one, *update );
1472  problem_->updateSolution( update, true );
1473 
1474  // Update residual norm ( r -= C_*C_'*r_ )
1475  MVT::MvTimesMatAddMv( -one, *Ctmp, Ctr, one, *r_ );
1476 
1477  // We recycled space from previous call
1478  prime_iterations = 0;
1479 
1480  }
1481  else {
1482 
1483  // Do one cycle of Gmres to "prime the pump" if there is no subspace to recycle
1484  printer_->stream(Debug) << " No recycled subspace available for RHS index " << currIdx[0] << std::endl << std::endl;
1485 
1486  Teuchos::ParameterList primeList;
1487 
1488  // Tell the block solver that the block size is one.
1489  primeList.set("Num Blocks",numBlocks_);
1490  primeList.set("Recycled Blocks",0);
1491 
1492  // Create GCRODR iterator object to perform one cycle of GMRES.
1493  RCP<GCRODRIter<ScalarType,MV,OP> > gcrodr_prime_iter;
1494  gcrodr_prime_iter = rcp( new GCRODRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,primeList) );
1495 
1496  // Create the first block in the current Krylov basis (residual).
1497  problem_->computeCurrPrecResVec( &*r_ );
1498  index.resize( 1 ); index[0] = 0;
1499  RCP<MV> v0 = MVT::CloneViewNonConst( *V_, index );
1500  MVT::SetBlock(*r_,index,*v0); // V(:,0) = r
1501 
1502  // Set the new state and initialize the solver.
1504  index.resize( numBlocks_+1 );
1505  for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1506  newstate.V = MVT::CloneViewNonConst( *V_, index );
1507  newstate.U = Teuchos::null;
1508  newstate.C = Teuchos::null;
1509  newstate.H = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, numBlocks_+1, numBlocks_, recycledBlocks_+1, recycledBlocks_+1 ) );
1510  newstate.B = Teuchos::null;
1511  newstate.curDim = 0;
1512  gcrodr_prime_iter->initialize(newstate);
1513 
1514  // Perform one cycle of GMRES
1515  bool primeConverged = false;
1516  try {
1517  gcrodr_prime_iter->iterate();
1518 
1519  // Check convergence first
1520  if ( convTest_->getStatus() == Passed ) {
1521  // we have convergence
1522  primeConverged = true;
1523  }
1524  }
1525  catch (const GCRODRIterOrthoFailure &e) {
1526  // Try to recover the most recent least-squares solution
1527  gcrodr_prime_iter->updateLSQR( gcrodr_prime_iter->getCurSubspaceDim() );
1528 
1529  // Check to see if the most recent least-squares solution yielded convergence.
1530  sTest_->checkStatus( &*gcrodr_prime_iter );
1531  if (convTest_->getStatus() == Passed)
1532  primeConverged = true;
1533  }
1534  catch (const std::exception &e) {
1535  printer_->stream(Errors) << "Error! Caught exception in GCRODRIter::iterate() at iteration "
1536  << gcrodr_prime_iter->getNumIters() << std::endl
1537  << e.what() << std::endl;
1538  throw;
1539  }
1540  // Record number of iterations in generating initial recycle spacec
1541  prime_iterations = gcrodr_prime_iter->getNumIters();
1542 
1543  // Update the linear problem.
1544  RCP<MV> update = gcrodr_prime_iter->getCurrentUpdate();
1545  problem_->updateSolution( update, true );
1546 
1547  // Get the state.
1548  newstate = gcrodr_prime_iter->getState();
1549  int p = newstate.curDim;
1550 
1551  // Compute harmonic Ritz vectors
1552  // NOTE: The storage for the harmonic Ritz vectors (PP) is made one column larger
1553  // just in case we split a complex conjugate pair.
1554  // NOTE: Generate a recycled subspace only if we have enough vectors. If we converged
1555  // too early, move on to the next linear system and try to generate a subspace again.
1556  if (recycledBlocks_ < p+1) {
1557  int info = 0;
1558  RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > PPtmp = rcp (new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View, *PP_, p, recycledBlocks_+1 ) );
1559  // getHarmonicVecs1 assumes PP has recycledBlocks_+1 columns available
1560  keff = getHarmonicVecs1( p, *newstate.H, *PPtmp );
1561  // Hereafter, only keff columns of PP are needed
1562  PPtmp = rcp (new Teuchos::SerialDenseMatrix<int,ScalarType> ( Teuchos::View, *PP_, p, keff ) );
1563  // Now get views into C, U, V
1564  index.resize(keff);
1565  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1566  RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1567  RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1568  RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1569  index.resize(p);
1570  for (int ii=0; ii < p; ++ii) { index[ii] = ii; }
1571  RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1572 
1573  // Form U (the subspace to recycle)
1574  // U = newstate.V(:,1:p) * PP;
1575  MVT::MvTimesMatAddMv( one, *Vtmp, *PPtmp, zero, *U1tmp );
1576 
1577  // Form orthonormalized C and adjust U so that C = A*U
1578 
1579  // First, compute [Q, R] = qr(H*P);
1580 
1581  // Step #1: Form HP = H*P
1582  Teuchos::SerialDenseMatrix<int,ScalarType> Htmp( Teuchos::View, *H2_, p+1, p, recycledBlocks_+1,recycledBlocks_+1);
1583  Teuchos::SerialDenseMatrix<int,ScalarType> HPtmp( Teuchos::View, *HP_, p+1, keff );
1584  HPtmp.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, Htmp, *PPtmp, zero );
1585 
1586  // Step #1.5: Perform workspace size query for QR
1587  // factorization of HP. On input, lwork must be -1.
1588  // _GEQRF will put the workspace size in work_[0].
1589  int lwork = -1;
1590  tau_.resize (keff);
1591  lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1592  HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
1593  TEUCHOS_TEST_FOR_EXCEPTION(
1594  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve:"
1595  " LAPACK's _GEQRF failed to compute a workspace size.");
1596 
1597  // Step #2: Compute QR factorization of HP
1598  //
1599  // NOTE (mfh 17 Apr 2014) LAPACK promises that the value of
1600  // work_[0] after the workspace query will fit in int. This
1601  // justifies the cast. We call real() first because
1602  // static_cast from std::complex to int doesn't work.
1603  lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work_[0])));
1604  work_.resize (lwork); // Allocate workspace for the QR factorization
1605  lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1606  HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
1607  TEUCHOS_TEST_FOR_EXCEPTION(
1608  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve:"
1609  " LAPACK's _GEQRF failed to compute a QR factorization.");
1610 
1611  // Step #3: Explicitly construct Q and R factors
1612  // NOTE: The upper triangular part of HP is copied into R and HP becomes Q.
1613  Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *R_, keff, keff );
1614  for (int ii = 0; ii < keff; ++ii) {
1615  for (int jj = ii; jj < keff; ++jj) {
1616  Rtmp(ii,jj) = HPtmp(ii,jj);
1617  }
1618  }
1619  // NOTE (mfh 17 Apr 2014): Teuchos::LAPACK's wrapper for
1620  // UNGQR dispatches to the correct Scalar-specific routine.
1621  // It calls {S,D}ORGQR if Scalar is real, and {C,Z}UNGQR if
1622  // Scalar is complex.
1623  lapack.UNGQR (HPtmp.numRows (), HPtmp.numCols (), HPtmp.numCols (),
1624  HPtmp.values (), HPtmp.stride (), &tau_[0], &work_[0],
1625  lwork, &info);
1626  TEUCHOS_TEST_FOR_EXCEPTION(
1627  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1628  "LAPACK's _UNGQR failed to construct the Q factor.");
1629 
1630  // Now we have [Q,R] = qr(H*P)
1631 
1632  // Now compute C = V(:,1:p+1) * Q
1633  index.resize (p + 1);
1634  for (int ii = 0; ii < (p+1); ++ii) {
1635  index[ii] = ii;
1636  }
1637  Vtmp = MVT::CloneView( *V_, index ); // need new view into V (p+1 vectors now; needed p above)
1638  MVT::MvTimesMatAddMv( one, *Vtmp, HPtmp, zero, *Ctmp );
1639 
1640  // Finally, compute U = U*R^{-1}.
1641  // This unfortuntely requires me to form R^{-1} explicitly and execute U = U * R^{-1}, as
1642  // backsolve capabilities don't exist in the Belos::MultiVec class
1643 
1644  // Step #1: First, compute LU factorization of R
1645  ipiv_.resize(Rtmp.numRows());
1646  lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
1647  TEUCHOS_TEST_FOR_EXCEPTION(
1648  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1649  "LAPACK's _GETRF failed to compute an LU factorization.");
1650 
1651  // FIXME (mfh 17 Apr 2014) We have to compute the explicit
1652  // inverse of R here because Belos::MultiVecTraits doesn't
1653  // have a triangular solve (where the triangular matrix is
1654  // globally replicated and the "right-hand side" is the
1655  // distributed MultiVector).
1656 
1657  // Step #2: Form inv(R)
1658  lwork = Rtmp.numRows();
1659  work_.resize(lwork);
1660  lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1661  TEUCHOS_TEST_FOR_EXCEPTION(
1662  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1663  "LAPACK's _GETRI failed to invert triangular matrix.");
1664 
1665  // Step #3: Let U = U * R^{-1}
1666  MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
1667 
1668  printer_->stream(Debug)
1669  << " Generated recycled subspace using RHS index " << currIdx[0]
1670  << " of dimension " << keff << std::endl << std::endl;
1671 
1672  } // if (recycledBlocks_ < p+1)
1673 
1674  // Return to outer loop if the priming solve converged, set the next linear system.
1675  if (primeConverged) {
1676  // Inform the linear problem that we are finished with this block linear system.
1677  problem_->setCurrLS();
1678 
1679  // Update indices for the linear systems to be solved.
1680  numRHS2Solve -= 1;
1681  if (numRHS2Solve > 0) {
1682  currIdx[0]++;
1683  problem_->setLSIndex (currIdx); // Set the next indices
1684  }
1685  else {
1686  currIdx.resize (numRHS2Solve);
1687  }
1688 
1689  continue;
1690  }
1691  } // if (keff > 0) ...
1692 
1693  // Prepare for the Gmres iterations with the recycled subspace.
1694 
1695  // Set the current number of recycled blocks and subspace dimension with the GCRO-DR iteration.
1696  gcrodr_iter->setSize( keff, numBlocks_ );
1697 
1698  // Reset the number of iterations.
1699  gcrodr_iter->resetNumIters(prime_iterations);
1700 
1701  // Reset the number of calls that the status test output knows about.
1702  outputTest_->resetNumCalls();
1703 
1704  // Compute the residual after the priming solve, it will be the first block in the current Krylov basis.
1705  problem_->computeCurrPrecResVec( &*r_ );
1706  index.resize( 1 ); index[0] = 0;
1707  RCP<MV> v0 = MVT::CloneViewNonConst( *V_, index );
1708  MVT::SetBlock(*r_,index,*v0); // V(:,0) = r
1709 
1710  // Set the new state and initialize the solver.
1712  index.resize( numBlocks_+1 );
1713  for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1714  newstate.V = MVT::CloneViewNonConst( *V_, index );
1715  index.resize( keff );
1716  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1717  newstate.C = MVT::CloneViewNonConst( *C_, index );
1718  newstate.U = MVT::CloneViewNonConst( *U_, index );
1719  newstate.B = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, keff, numBlocks_, 0, keff ) );
1720  newstate.H = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, numBlocks_+1, numBlocks_, keff, keff ) );
1721  newstate.curDim = 0;
1722  gcrodr_iter->initialize(newstate);
1723 
1724  // variables needed for inner loop
1725  int numRestarts = 0;
1726  while(1) {
1727 
1728  // tell gcrodr_iter to iterate
1729  try {
1730  gcrodr_iter->iterate();
1731 
1733  //
1734  // check convergence first
1735  //
1737  if ( convTest_->getStatus() == Passed ) {
1738  // we have convergence
1739  break; // break from while(1){gcrodr_iter->iterate()}
1740  }
1742  //
1743  // check for maximum iterations
1744  //
1746  else if ( maxIterTest_->getStatus() == Passed ) {
1747  // we don't have convergence
1748  isConverged = false;
1749  break; // break from while(1){gcrodr_iter->iterate()}
1750  }
1752  //
1753  // check for restarting, i.e. the subspace is full
1754  //
1756  else if ( gcrodr_iter->getCurSubspaceDim() == gcrodr_iter->getMaxSubspaceDim() ) {
1757 
1758  // Update the recycled subspace even if we have hit the maximum number of restarts.
1759 
1760  // Update the linear problem.
1761  RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1762  problem_->updateSolution( update, true );
1763 
1764  buildRecycleSpace2(gcrodr_iter);
1765 
1766  printer_->stream(Debug)
1767  << " Generated new recycled subspace using RHS index "
1768  << currIdx[0] << " of dimension " << keff << std::endl
1769  << std::endl;
1770 
1771  // NOTE: If we have hit the maximum number of restarts then quit
1772  if (numRestarts >= maxRestarts_) {
1773  isConverged = false;
1774  break; // break from while(1){gcrodr_iter->iterate()}
1775  }
1776  numRestarts++;
1777 
1778  printer_->stream(Debug)
1779  << " Performing restart number " << numRestarts << " of "
1780  << maxRestarts_ << std::endl << std::endl;
1781 
1782  // Create the restart vector (first block in the current Krylov basis)
1783  problem_->computeCurrPrecResVec( &*r_ );
1784  index.resize( 1 ); index[0] = 0;
1785  RCP<MV> v00 = MVT::CloneViewNonConst( *V_, index );
1786  MVT::SetBlock(*r_,index,*v00); // V(:,0) = r
1787 
1788  // Set the new state and initialize the solver.
1789  GCRODRIterState<ScalarType,MV> restartState;
1790  index.resize( numBlocks_+1 );
1791  for (int ii=0; ii<(numBlocks_+1); ++ii) { index[ii] = ii; }
1792  restartState.V = MVT::CloneViewNonConst( *V_, index );
1793  index.resize( keff );
1794  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1795  restartState.U = MVT::CloneViewNonConst( *U_, index );
1796  restartState.C = MVT::CloneViewNonConst( *C_, index );
1797  restartState.B = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, keff, numBlocks_, 0, keff ) );
1798  restartState.H = rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, numBlocks_+1, numBlocks_, keff, keff ) );
1799  restartState.curDim = 0;
1800  gcrodr_iter->initialize(restartState);
1801 
1802 
1803  } // end of restarting
1804 
1806  //
1807  // we returned from iterate(), but none of our status tests Passed.
1808  // something is wrong, and it is probably our fault.
1809  //
1811 
1812  else {
1813  TEUCHOS_TEST_FOR_EXCEPTION(
1814  true, std::logic_error, "Belos::GCRODRSolMgr::solve: "
1815  "Invalid return from GCRODRIter::iterate().");
1816  }
1817  }
1818  catch (const GCRODRIterOrthoFailure &e) {
1819  // Try to recover the most recent least-squares solution
1820  gcrodr_iter->updateLSQR( gcrodr_iter->getCurSubspaceDim() );
1821 
1822  // Check to see if the most recent least-squares solution yielded convergence.
1823  sTest_->checkStatus( &*gcrodr_iter );
1824  if (convTest_->getStatus() != Passed)
1825  isConverged = false;
1826  break;
1827  }
1828  catch (const std::exception& e) {
1829  printer_->stream(Errors)
1830  << "Error! Caught exception in GCRODRIter::iterate() at iteration "
1831  << gcrodr_iter->getNumIters() << std::endl << e.what() << std::endl;
1832  throw;
1833  }
1834  }
1835 
1836  // Compute the current solution.
1837  // Update the linear problem.
1838  RCP<MV> update = gcrodr_iter->getCurrentUpdate();
1839  problem_->updateSolution( update, true );
1840 
1841  // Inform the linear problem that we are finished with this block linear system.
1842  problem_->setCurrLS();
1843 
1844  // If we didn't build a recycle space this solve but ran at least k iterations,
1845  // force build of new recycle space
1846 
1847  if (!builtRecycleSpace_) {
1848  buildRecycleSpace2(gcrodr_iter);
1849  printer_->stream(Debug)
1850  << " Generated new recycled subspace using RHS index " << currIdx[0]
1851  << " of dimension " << keff << std::endl << std::endl;
1852  }
1853 
1854  // Update indices for the linear systems to be solved.
1855  numRHS2Solve -= 1;
1856  if (numRHS2Solve > 0) {
1857  currIdx[0]++;
1858  problem_->setLSIndex (currIdx); // Set the next indices
1859  }
1860  else {
1861  currIdx.resize (numRHS2Solve);
1862  }
1863  } // while (numRHS2Solve > 0)
1864  }
1865 
1866  sTest_->print (printer_->stream (FinalSummary)); // print final summary
1867 
1868  // print timing information
1869 #ifdef BELOS_TEUCHOS_TIME_MONITOR
1870  // Calling summarize() can be expensive, so don't call unless the
1871  // user wants to print out timing details. summarize() will do all
1872  // the work even if it's passed a "black hole" output stream.
1873  if (verbosity_ & TimingDetails)
1874  Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
1875 #endif // BELOS_TEUCHOS_TIME_MONITOR
1876 
1877  // get iteration information for this solve
1878  numIters_ = maxIterTest_->getNumIters ();
1879 
1880  // Save the convergence test value ("achieved tolerance") for this
1881  // solve. This solver (unlike BlockGmresSolMgr) always has two
1882  // residual norm status tests: an explicit and an implicit test.
1883  // The master convergence test convTest_ is a SEQ combo of the
1884  // implicit resp. explicit tests. If the implicit test never
1885  // passes, then the explicit test won't ever be executed. This
1886  // manifests as expConvTest_->getTestValue()->size() < 1. We deal
1887  // with this case by using the values returned by
1888  // impConvTest_->getTestValue().
1889  {
1890  const std::vector<MagnitudeType>* pTestValues = expConvTest_->getTestValue();
1891  if (pTestValues == NULL || pTestValues->size() < 1) {
1892  pTestValues = impConvTest_->getTestValue();
1893  }
1894  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
1895  "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1896  "method returned NULL. Please report this bug to the Belos developers.");
1897  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
1898  "Belos::GCRODRSolMgr::solve(): The implicit convergence test's getTestValue() "
1899  "method returned a vector of length zero. Please report this bug to the "
1900  "Belos developers.");
1901 
1902  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
1903  // achieved tolerances for all vectors in the current solve(), or
1904  // just for the vectors from the last deflation?
1905  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
1906  }
1907 
1908  return isConverged ? Converged : Unconverged; // return from solve()
1909 }
1910 
1911 // Given existing recycle space and Krylov space, build new recycle space
1912 template<class ScalarType, class MV, class OP>
1914 
1915  MagnitudeType one = Teuchos::ScalarTraits<MagnitudeType>::one();
1916  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1917 
1918  std::vector<MagnitudeType> d(keff);
1919  std::vector<ScalarType> dscalar(keff);
1920  std::vector<int> index(numBlocks_+1);
1921 
1922  // Get the state
1923  GCRODRIterState<ScalarType,MV> oldState = gcrodr_iter->getState();
1924  int p = oldState.curDim;
1925 
1926  // insufficient new information to update recycle space
1927  if (p<1) return;
1928 
1929  // Take the norm of the recycled vectors
1930  {
1931  index.resize(keff);
1932  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1933  Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1934  d.resize(keff);
1935  dscalar.resize(keff);
1936  MVT::MvNorm( *Utmp, d );
1937  for (int i=0; i<keff; ++i) {
1938  d[i] = one / d[i];
1939  dscalar[i] = (ScalarType)d[i];
1940  }
1941  MVT::MvScale( *Utmp, dscalar );
1942  }
1943 
1944  // Get view into current "full" upper Hessnburg matrix
1945  Teuchos::RCP<Teuchos::SerialDenseMatrix<int,ScalarType> > H2tmp = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::View, *H2_, p+keff+1, p+keff ) );
1946 
1947  // Insert D into the leading keff x keff block of H2
1948  for (int i=0; i<keff; ++i) {
1949  (*H2tmp)(i,i) = d[i];
1950  }
1951 
1952  // Compute the harmoic Ritz pairs for the generalized eigenproblem
1953  // getHarmonicVecs2 assumes PP has recycledBlocks_+1 columns available
1954  int keff_new;
1955  {
1956  Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, p+keff, recycledBlocks_+1 );
1957  keff_new = getHarmonicVecs2( keff, p, *H2tmp, oldState.V, PPtmp );
1958  }
1959 
1960  // Code to form new U, C
1961  // U = [U V(:,1:p)] * P; (in two steps)
1962 
1963  // U(:,1:keff) = matmul(U(:,1:keff_old),PP(1:keff_old,1:keff)) (step 1)
1964  Teuchos::RCP<MV> U1tmp;
1965  {
1966  index.resize( keff );
1967  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1968  Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1969  index.resize( keff_new );
1970  for (int ii=0; ii<keff_new; ++ii) { index[ii] = ii; }
1971  U1tmp = MVT::CloneViewNonConst( *U1_, index );
1972  Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, keff, keff_new );
1973  MVT::MvTimesMatAddMv( one, *Utmp, PPtmp, zero, *U1tmp );
1974  }
1975 
1976  // U(:,1:keff) = U(:,1:keff) + matmul(V(:,1:m-k),PP(keff_old+1:m-k+keff_old,1:keff)) (step 2)
1977  {
1978  index.resize(p);
1979  for (int ii=0; ii < p; ii++) { index[ii] = ii; }
1980  Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1981  Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, p, keff_new, keff );
1982  MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *U1tmp );
1983  }
1984 
1985  // Form HP = H*P
1986  Teuchos::SerialDenseMatrix<int,ScalarType> HPtmp( Teuchos::View, *HP_, p+keff+1, keff_new );
1987  {
1988  Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *PP_, p+keff, keff_new );
1989  HPtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,*H2tmp,PPtmp,zero);
1990  }
1991 
1992  // Workspace size query for QR factorization of HP (the worksize will be placed in work_[0])
1993  int info = 0, lwork = -1;
1994  tau_.resize (keff_new);
1995  lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
1996  HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
1997  TEUCHOS_TEST_FOR_EXCEPTION(
1998  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
1999  "LAPACK's _GEQRF failed to compute a workspace size.");
2000 
2001  // NOTE (mfh 18 Apr 2014) LAPACK promises that the value of work_[0]
2002  // after the workspace query will fit in int. This justifies the
2003  // cast. We call real() first because static_cast from std::complex
2004  // to int doesn't work.
2005  lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work_[0])));
2006  work_.resize (lwork); // Allocate workspace for the QR factorization
2007  lapack.GEQRF (HPtmp.numRows (), HPtmp.numCols (), HPtmp.values (),
2008  HPtmp.stride (), &tau_[0], &work_[0], lwork, &info);
2009  TEUCHOS_TEST_FOR_EXCEPTION(
2010  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
2011  "LAPACK's _GEQRF failed to compute a QR factorization.");
2012 
2013  // Explicitly construct Q and R factors
2014  // NOTE: The upper triangular part of HP is copied into R and HP becomes Q.
2015  Teuchos::SerialDenseMatrix<int,ScalarType> Rtmp( Teuchos::View, *R_, keff_new, keff_new );
2016  for(int i=0;i<keff_new;i++) { for(int j=i;j<keff_new;j++) Rtmp(i,j) = HPtmp(i,j); }
2017 
2018  // NOTE (mfh 18 Apr 2014): Teuchos::LAPACK's wrapper for UNGQR
2019  // dispatches to the correct Scalar-specific routine. It calls
2020  // {S,D}ORGQR if Scalar is real, and {C,Z}UNGQR if Scalar is
2021  // complex.
2022  lapack.UNGQR (HPtmp.numRows (), HPtmp.numCols (), HPtmp.numCols (),
2023  HPtmp.values (), HPtmp.stride (), &tau_[0], &work_[0],
2024  lwork, &info);
2025  TEUCHOS_TEST_FOR_EXCEPTION(
2026  info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve: "
2027  "LAPACK's _UNGQR failed to construct the Q factor.");
2028 
2029  // Form orthonormalized C and adjust U accordingly so that C = A*U
2030  // C = [C V] * Q;
2031 
2032  // C(:,1:keff) = matmul(C(:,1:keff_old),QQ(1:keff_old,1:keff))
2033  {
2034  Teuchos::RCP<MV> C1tmp;
2035  {
2036  index.resize(keff);
2037  for (int i=0; i < keff; i++) { index[i] = i; }
2038  Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
2039  index.resize(keff_new);
2040  for (int i=0; i < keff_new; i++) { index[i] = i; }
2041  C1tmp = MVT::CloneViewNonConst( *C1_, index );
2042  Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *HP_, keff, keff_new );
2043  MVT::MvTimesMatAddMv( one, *Ctmp, PPtmp, zero, *C1tmp );
2044  }
2045  // Now compute C += V(:,1:p+1) * Q
2046  {
2047  index.resize( p+1 );
2048  for (int i=0; i < p+1; ++i) { index[i] = i; }
2049  Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
2050  Teuchos::SerialDenseMatrix<int,ScalarType> PPtmp( Teuchos::View, *HP_, p+1, keff_new, keff, 0 );
2051  MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *C1tmp );
2052  }
2053  }
2054 
2055  // C_ = C1_; (via a swap)
2056  std::swap(C_, C1_);
2057 
2058  // Finally, compute U_ = U_*R^{-1}
2059  // First, compute LU factorization of R
2060  ipiv_.resize(Rtmp.numRows());
2061  lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
2062  TEUCHOS_TEST_FOR_EXCEPTION(info != 0,GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
2063 
2064  // Now, form inv(R)
2065  lwork = Rtmp.numRows();
2066  work_.resize(lwork);
2067  lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
2068  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
2069 
2070  {
2071  index.resize(keff_new);
2072  for (int i=0; i < keff_new; i++) { index[i] = i; }
2073  Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
2074  MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
2075  }
2076 
2077  // Set the current number of recycled blocks and subspace dimension with the GCRO-DR iteration.
2078  if (keff != keff_new) {
2079  keff = keff_new;
2080  gcrodr_iter->setSize( keff, numBlocks_ );
2081  // Important to zero this out before next cyle
2082  Teuchos::SerialDenseMatrix<int,ScalarType> b1( Teuchos::View, *H2_, recycledBlocks_+2, 1, 0, recycledBlocks_ );
2083  b1.putScalar(zero);
2084  }
2085 
2086 }
2087 
2088 
2089 // Compute the harmonic eigenpairs of the projected, dense system.
2090 template<class ScalarType, class MV, class OP>
2092  const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
2093  Teuchos::SerialDenseMatrix<int,ScalarType>& PP) {
2094  int i, j;
2095  bool xtraVec = false;
2096  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
2097  //ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
2098 
2099  // Real and imaginary eigenvalue components
2100  std::vector<MagnitudeType> wr(m), wi(m);
2101 
2102  // Real and imaginary (right) eigenvectors; Don't zero out matrix when constructing
2103  Teuchos::SerialDenseMatrix<int,ScalarType> vr(m,m,false);
2104 
2105  // Magnitude of harmonic Ritz values
2106  std::vector<MagnitudeType> w(m);
2107 
2108  // Sorted order of harmonic Ritz values, also used for DGEEV
2109  std::vector<int> iperm(m);
2110 
2111  // Size of workspace and workspace for DGEEV
2112  int lwork = 4*m;
2113  std::vector<ScalarType> work(lwork);
2114  std::vector<MagnitudeType> rwork(lwork);
2115 
2116  // Output info
2117  int info = 0;
2118 
2119  // Set flag indicating recycle space has been generated this solve
2120  builtRecycleSpace_ = true;
2121 
2122  // Solve linear system: H_m^{-H}*e_m
2123  Teuchos::SerialDenseMatrix<int, ScalarType> HHt( HH, Teuchos::TRANS );
2124  Teuchos::SerialDenseVector<int, ScalarType> e_m( m );
2125  e_m[m-1] = one;
2126  lapack.GESV(m, 1, HHt.values(), HHt.stride(), &iperm[0], e_m.values(), e_m.stride(), &info);
2127  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve(): LAPACK GESV failed to compute a solution.");
2128 
2129  // Compute H_m + d*H_m^{-H}*e_m*e_m^H
2130  ScalarType d = HH(m, m-1) * HH(m, m-1);
2131  Teuchos::SerialDenseMatrix<int, ScalarType> harmHH( Teuchos::Copy, HH, HH.numRows(), HH.numCols() );
2132  for( i=0; i<m; ++i )
2133  harmHH(i, m-1) += d * e_m[i];
2134 
2135  // Revise to do query for optimal workspace first
2136  // Create simple storage for the left eigenvectors, which we don't care about.
2137  const int ldvl = m;
2138  ScalarType* vl = 0;
2139  //lapack.GEEV('N', 'V', m, harmHH.values(), harmHH.stride(), &wr[0], &wi[0],
2140  // vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &info);
2141  lapack.GEEV('N', 'V', m, harmHH.values(), harmHH.stride(), &wr[0], &wi[0],
2142  vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
2143  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure,"Belos::GCRODRSolMgr::solve(): LAPACK GEEV failed to compute eigensolutions.");
2144 
2145  // Construct magnitude of each harmonic Ritz value
2146  for( i=0; i<m; ++i )
2147  w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] );
2148 
2149  // Construct magnitude of each harmonic Ritz value
2150  this->sort(w, m, iperm);
2151 
2152  const bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
2153 
2154  // Select recycledBlocks_ smallest eigenvectors
2155  for( i=0; i<recycledBlocks_; ++i ) {
2156  for( j=0; j<m; j++ ) {
2157  PP(j,i) = vr(j,iperm[i]);
2158  }
2159  }
2160 
2161  if(!scalarTypeIsComplex) {
2162 
2163  // Determine exact size for PP (i.e., determine if we need to store an additional vector)
2164  if (wi[iperm[recycledBlocks_-1]] != 0.0) {
2165  int countImag = 0;
2166  for ( i=0; i<recycledBlocks_; ++i ) {
2167  if (wi[iperm[i]] != 0.0)
2168  countImag++;
2169  }
2170  // Check to see if this count is even or odd:
2171  if (countImag % 2)
2172  xtraVec = true;
2173  }
2174 
2175  if (xtraVec) { // we need to store one more vector
2176  if (wi[iperm[recycledBlocks_-1]] > 0.0) { // I picked the "real" component
2177  for( j=0; j<m; ++j ) { // so get the "imag" component
2178  PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]+1);
2179  }
2180  }
2181  else { // I picked the "imag" component
2182  for( j=0; j<m; ++j ) { // so get the "real" component
2183  PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]-1);
2184  }
2185  }
2186  }
2187 
2188  }
2189 
2190  // Return whether we needed to store an additional vector
2191  if (xtraVec) {
2192  return recycledBlocks_+1;
2193  }
2194  else {
2195  return recycledBlocks_;
2196  }
2197 
2198 }
2199 
2200 // Compute the harmonic eigenpairs of the projected, dense system.
2201 template<class ScalarType, class MV, class OP>
2203  const Teuchos::SerialDenseMatrix<int,ScalarType>& HH,
2204  const Teuchos::RCP<const MV>& VV,
2205  Teuchos::SerialDenseMatrix<int,ScalarType>& PP) {
2206  int i, j;
2207  int m2 = HH.numCols();
2208  bool xtraVec = false;
2209  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
2210  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
2211  std::vector<int> index;
2212 
2213  // Real and imaginary eigenvalue components
2214  std::vector<MagnitudeType> wr(m2), wi(m2);
2215 
2216  // Magnitude of harmonic Ritz values
2217  std::vector<MagnitudeType> w(m2);
2218 
2219  // Real and imaginary (right) eigenvectors; Don't zero out matrix when constructing
2220  Teuchos::SerialDenseMatrix<int,ScalarType> vr(m2,m2,false);
2221 
2222  // Sorted order of harmonic Ritz values
2223  std::vector<int> iperm(m2);
2224 
2225  // Set flag indicating recycle space has been generated this solve
2226  builtRecycleSpace_ = true;
2227 
2228  // Form matrices for generalized eigenproblem
2229 
2230  // B = H2' * H2; Don't zero out matrix when constructing
2231  Teuchos::SerialDenseMatrix<int,ScalarType> B(m2,m2,false);
2232  B.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,HH,HH,zero);
2233 
2234  // A_tmp = | C'*U 0 |
2235  // | V_{m+1}'*U I |
2236  Teuchos::SerialDenseMatrix<int,ScalarType> A_tmp( keffloc+m+1, keffloc+m );
2237 
2238  // A_tmp(1:keffloc,1:keffloc) = C' * U;
2239  index.resize(keffloc);
2240  for (i=0; i<keffloc; ++i) { index[i] = i; }
2241  Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
2242  Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
2243  Teuchos::SerialDenseMatrix<int,ScalarType> A11( Teuchos::View, A_tmp, keffloc, keffloc );
2244  MVT::MvTransMv( one, *Ctmp, *Utmp, A11 );
2245 
2246  // A_tmp(keffloc+1:m-k+keffloc+1,1:keffloc) = V' * U;
2247  Teuchos::SerialDenseMatrix<int,ScalarType> A21( Teuchos::View, A_tmp, m+1, keffloc, keffloc );
2248  index.resize(m+1);
2249  for (i=0; i < m+1; i++) { index[i] = i; }
2250  Teuchos::RCP<const MV> Vp = MVT::CloneView( *VV, index );
2251  MVT::MvTransMv( one, *Vp, *Utmp, A21 );
2252 
2253  // A_tmp(keffloc+1:m-k+keffloc,keffloc+1:m-k+keffloc) = eye(m-k);
2254  for( i=keffloc; i<keffloc+m; i++ ) {
2255  A_tmp(i,i) = one;
2256  }
2257 
2258  // A = H2' * A_tmp;
2259  Teuchos::SerialDenseMatrix<int,ScalarType> A( m2, A_tmp.numCols() );
2260  A.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, one, HH, A_tmp, zero );
2261 
2262  // Compute k smallest harmonic Ritz pairs
2263  // SUBROUTINE DGGEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
2264  // ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO,
2265  // IHI, LSCALE, RSCALE, ABNRM, BBNRM, RCONDE,
2266  // RCONDV, WORK, LWORK, IWORK, BWORK, INFO )
2267  // MLP: 'SCALING' in DGGEVX generates incorrect eigenvalues. Therefore, only permuting
2268  char balanc='P', jobvl='N', jobvr='V', sense='N';
2269  int ld = A.numRows();
2270  int lwork = 6*ld;
2271  int ldvl = ld, ldvr = ld;
2272  int info = 0,ilo = 0,ihi = 0;
2273  MagnitudeType abnrm = 0.0, bbnrm = 0.0;
2274  ScalarType *vl = 0; // This is never referenced by dggevx if jobvl == 'N'
2275  std::vector<ScalarType> beta(ld);
2276  std::vector<ScalarType> work(lwork);
2277  std::vector<MagnitudeType> rwork(lwork);
2278  std::vector<MagnitudeType> lscale(ld), rscale(ld);
2279  std::vector<MagnitudeType> rconde(ld), rcondv(ld);
2280  std::vector<int> iwork(ld+6);
2281  int *bwork = 0; // If sense == 'N', bwork is never referenced
2282  //lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
2283  // &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
2284  // &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &iwork[0], bwork, &info);
2285  lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
2286  &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
2287  &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &rwork[0],
2288  &iwork[0], bwork, &info);
2289  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, GCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve(): LAPACK GGEVX failed to compute eigensolutions.");
2290 
2291  // Construct magnitude of each harmonic Ritz value
2292  // NOTE : Forming alpha/beta *should* be okay here, given assumptions on construction of matrix pencil above
2293  for( i=0; i<ld; i++ ) {
2294  w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot (wr[i]*wr[i] + wi[i]*wi[i]) /
2295  Teuchos::ScalarTraits<ScalarType>::magnitude (beta[i]);
2296  }
2297 
2298  // Construct magnitude of each harmonic Ritz value
2299  this->sort(w,ld,iperm);
2300 
2301  const bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
2302 
2303  // Select recycledBlocks_ smallest eigenvectors
2304  for( i=0; i<recycledBlocks_; i++ ) {
2305  for( j=0; j<ld; j++ ) {
2306  PP(j,i) = vr(j,iperm[ld-recycledBlocks_+i]);
2307  }
2308  }
2309 
2310  if(!scalarTypeIsComplex) {
2311 
2312  // Determine exact size for PP (i.e., determine if we need to store an additional vector)
2313  if (wi[iperm[ld-recycledBlocks_]] != 0.0) {
2314  int countImag = 0;
2315  for ( i=ld-recycledBlocks_; i<ld; i++ ) {
2316  if (wi[iperm[i]] != 0.0)
2317  countImag++;
2318  }
2319  // Check to see if this count is even or odd:
2320  if (countImag % 2)
2321  xtraVec = true;
2322  }
2323 
2324  if (xtraVec) { // we need to store one more vector
2325  if (wi[iperm[ld-recycledBlocks_]] > 0.0) { // I picked the "real" component
2326  for( j=0; j<ld; j++ ) { // so get the "imag" component
2327  PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]+1);
2328  }
2329  }
2330  else { // I picked the "imag" component
2331  for( j=0; j<ld; j++ ) { // so get the "real" component
2332  PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]-1);
2333  }
2334  }
2335  }
2336 
2337  }
2338 
2339  // Return whether we needed to store an additional vector
2340  if (xtraVec) {
2341  return recycledBlocks_+1;
2342  }
2343  else {
2344  return recycledBlocks_;
2345  }
2346 
2347 }
2348 
2349 
2350 // This method sorts list of n floating-point numbers and return permutation vector
2351 template<class ScalarType, class MV, class OP>
2352 void GCRODRSolMgr<ScalarType,MV,OP,true>::sort(std::vector<MagnitudeType>& dlist, int n, std::vector<int>& iperm) {
2353  int l, r, j, i, flag;
2354  int RR2;
2355  MagnitudeType dRR, dK;
2356 
2357  // Initialize the permutation vector.
2358  for(j=0;j<n;j++)
2359  iperm[j] = j;
2360 
2361  if (n <= 1) return;
2362 
2363  l = n / 2 + 1;
2364  r = n - 1;
2365  l = l - 1;
2366  dRR = dlist[l - 1];
2367  dK = dlist[l - 1];
2368 
2369  RR2 = iperm[l - 1];
2370  while (r != 0) {
2371  j = l;
2372  flag = 1;
2373 
2374  while (flag == 1) {
2375  i = j;
2376  j = j + j;
2377 
2378  if (j > r + 1)
2379  flag = 0;
2380  else {
2381  if (j < r + 1)
2382  if (dlist[j] > dlist[j - 1]) j = j + 1;
2383 
2384  if (dlist[j - 1] > dK) {
2385  dlist[i - 1] = dlist[j - 1];
2386  iperm[i - 1] = iperm[j - 1];
2387  }
2388  else {
2389  flag = 0;
2390  }
2391  }
2392  }
2393  dlist[i - 1] = dRR;
2394  iperm[i - 1] = RR2;
2395 
2396  if (l == 1) {
2397  dRR = dlist [r];
2398  RR2 = iperm[r];
2399  dK = dlist[r];
2400  dlist[r] = dlist[0];
2401  iperm[r] = iperm[0];
2402  r = r - 1;
2403  }
2404  else {
2405  l = l - 1;
2406  dRR = dlist[l - 1];
2407  RR2 = iperm[l - 1];
2408  dK = dlist[l - 1];
2409  }
2410  }
2411  dlist[0] = dRR;
2412  iperm[0] = RR2;
2413 }
2414 
2415 
2416 template<class ScalarType, class MV, class OP>
2418  std::ostringstream out;
2419  out << "Belos::GCRODRSolMgr<...,"<<Teuchos::ScalarTraits<ScalarType>::name()<<">";
2420  out << "{";
2421  out << "Ortho Type: \"" << orthoType_ << "\"";
2422  out << ", Num Blocks: " <<numBlocks_;
2423  out << ", Num Recycle Blocks: " << recycledBlocks_;
2424  out << ", Max Restarts: " << maxRestarts_;
2425  out << "}";
2426  return out.str ();
2427 }
2428 
2429 } // namespace Belos
2430 
2431 #endif /* BELOS_GCRODR_SOLMGR_HPP */
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:87
Collection of types and exceptions used within the Belos solvers.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
GCRODRSolMgrLAPACKFailure(const std::string &what_arg)
Class which manages the output and verbosity of the Belos solvers.
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:119
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
GCRODRSolMgrLAPACKFailure is thrown when a nonzero value is retuned from an LAPACK call...
MagnitudeType achievedTol() const
Tolerance achieved by the last solve() invocation.
A factory class for generating StatusTestOutput objects.
void reset(const ResetType type)
Performs a reset of the solver manager specified by the ResetType. This informs the solver manager th...
An implementation of StatusTestResNorm using a family of residual norms.
GCRODRSolMgrLinearProblemFailure is thrown when the linear problem is not setup (i.e.
GCRODRSolMgrRecyclingFailure(const std::string &what_arg)
This class implements the GCRODR iteration, where a single-std::vector Krylov subspace is constructed...
std::ostream & printValidNames(std::ostream &out) const
Print all recognized MatOrthoManager names to the given ostream.
Belos::StatusTest class for specifying a maximum number of iterations.
int curDim
The current dimension of the reduction.
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< Belos::MatOrthoManager< Scalar, MV, OP > > makeMatOrthoManager(const std::string &ortho, const Teuchos::RCP< const OP > &M, const Teuchos::RCP< OutputManager< Scalar > > &outMan, const std::string &label, const Teuchos::RCP< Teuchos::ParameterList > &params)
Return an instance of the specified MatOrthoManager subclass.
bool isLOADetected() const
Return whether a loss of accuracy was detected by this solver during the most current solve...
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
std::string validNamesString() const
List (as a string) of recognized MatOrthoManager names.
A Belos::StatusTest class for specifying a maximum number of iterations.
GCRODRSolMgrOrthoFailure(const std::string &what_arg)
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:205
Pure virtual base class which describes the basic interface for a solver manager. ...
Teuchos::RCP< MV > V
The current Krylov basis.
Teuchos::RCP< MV > U
The recycled subspace and its projection.
Teuchos::RCP< const Teuchos::ParameterList > getDefaultParameters(const std::string &name) const
Default parameters for the given MatOrthoManager subclass.
A linear system to solve, and its associated information.
GCRODRSolMgrOrthoFailure is thrown when the orthogonalization manager is unable to generate orthonorm...
Class which describes the linear problem to be solved by the iterative solver.
GCRODRIterOrthoFailure is thrown when the GCRODRIter object is unable to compute independent directio...
Implementation of the GCRODR (Recycling GMRES) iterative linear solver.
Type traits class that says whether Teuchos::LAPACK has a valid implementation for the given ScalarTy...
bool isValidName(const std::string &name) const
Whether this factory recognizes the MatOrthoManager with the given name.
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:154
Belos concrete class for performing the GCRO-DR iteration.
Structure to contain pointers to GCRODRIter state variables.
GCRODRSolMgrRecyclingFailure is thrown when any problem occurs in using/creating the recycling subspa...
int getNumIters() const
Get the iteration count for the most recent call to solve().
int setTolerance(MagnitudeType tolerance)
Set the value of the tolerance.
GCRODRSolMgrLinearProblemFailure(const std::string &what_arg)
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
GCRODRSolMgr(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< Teuchos::ParameterList > &pl)
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > create(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, int mod, int printStates)
Create the StatusTestOutput object specified by the outputStyle.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem that needs to be solved.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace.
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.
Teuchos::RCP< MV > C
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:60
Base class for Belos::SolverManager subclasses which normally can only compile with ScalarType types ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
Teuchos::Array< Teuchos::RCP< Teuchos::Time > > getTimers() const
Return the timers for this object.
Belos concrete class for performing the block, flexible GMRES iteration.

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