Belos  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
BelosBlockGCRODRSolMgr.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Belos: Block Linear Solvers Package
4 //
5 // Copyright 2004-2016 NTESS and the Belos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
13 #ifndef BELOS_BLOCK_GCRODR_SOLMGR_HPP
14 #define BELOS_BLOCK_GCRODR_SOLMGR_HPP
15 
16 #include "BelosConfigDefs.hpp"
17 #include "BelosTypes.hpp"
19 #include "BelosLinearProblem.hpp"
20 #include "BelosSolverManager.hpp"
21 #include "BelosGmresIteration.hpp"
22 #include "BelosBlockGCRODRIter.hpp"
23 #include "BelosBlockGmresIter.hpp"
24 #include "BelosBlockFGmresIter.hpp"
27 #include "BelosStatusTestCombo.hpp"
29 #include "BelosOutputManager.hpp"
30 #include "Teuchos_LAPACK.hpp"
31 #include "Teuchos_as.hpp"
32 #ifdef BELOS_TEUCHOS_TIME_MONITOR
33 #include "Teuchos_TimeMonitor.hpp"
34 #endif // BELOS_TEUCHOS_TIME_MONITOR
35 
36 // MLP: Add links to examples here
37 
38 namespace Belos{
39 
41 
42 
48  public:
49  BlockGCRODRSolMgrLinearProblemFailure(const std::string& what_arg) : BelosError(what_arg) {}
50  };
51 
57  public:
58  BlockGCRODRSolMgrOrthoFailure(const std::string& what_arg) : BelosError(what_arg) {}
59  };
60 
66  public:
67  BlockGCRODRSolMgrLAPACKFailure(const std::string& what_arg) : BelosError(what_arg) {}
68  };
69 
75  public:
76  BlockGCRODRSolMgrRecyclingFailure(const std::string& what_arg) : BelosError(what_arg) {}
77  };
78 
80 
92 
93 template<class ScalarType, class MV, class OP>
94 class BlockGCRODRSolMgr : public SolverManager<ScalarType, MV, OP> {
95 private:
96 
100  typedef typename Teuchos::ScalarTraits<ScalarType>::magnitudeType MagnitudeType;
106 
107 public:
109 
110 
117 
145 
147  virtual ~BlockGCRODRSolMgr() {};
149 
152 
154  std::string description() const;
155 
157 
158 
160 
161 
164  return *problem_;
165  }
166 
169 
172 
174  MagnitudeType achievedTol() const {
175  return achievedTol_;
176  }
177 
179  int getNumIters() const { return numIters_; }
180 
182  bool isLOADetected() const { return loaDetected_; }
183 
185 
187 
188 
191  TEUCHOS_TEST_FOR_EXCEPTION(problem.is_null(), std::invalid_argument,
192  "Belos::BlockGCRODRSolMgr::setProblem: The input LinearProblem cannot be null.");
193 
194  // Check state of problem object before proceeding
195  if (! problem->isProblemSet()) {
196  const bool success = problem->setProblem();
197  TEUCHOS_TEST_FOR_EXCEPTION(success, std::runtime_error,
198  "Belos::BlockGCRODRSolMgr::setProblem: Calling the input LinearProblem's setProblem() method failed. This likely means that the "
199  "LinearProblem has a missing (null) matrix A, solution vector X, or right-hand side vector B. Please set these items in the LinearProblem and try again.");
200  }
201 
202  problem_ = problem;
203  }
204 
207 
209 
211 
212 
219  void reset (const ResetType type) {
220  if ((type & Belos::Problem) && ! problem_.is_null ())
221  problem_->setProblem();
222  else if (type & Belos::RecycleSubspace)
223  keff = 0;
224  }
225 
227 
229 
230 
246  // \returns ::ReturnType specifying:
249  ReturnType solve();
250 
252 
253 private:
254 
255  // Called by all constructors; Contains init instructions common to all constructors
256  void init();
257 
258  // Initialize solver state storage
259  void initializeStateStorage();
260 
261  // Recycling Methods
262  // Appending Function name by:
263  // "Kryl" indicates it is specialized for building a recycle space after an
264  // initial run of Block GMRES which generates an initial unaugmented block Krylov subspace
265  // "AugKryl" indicates it is specialized for building a recycle space from the augmented Krylov subspace
266 
267  // Functions which control the building of a recycle space
268  void buildRecycleSpaceKryl(int& keff, Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > block_gmres_iter);
269  void buildRecycleSpaceAugKryl(Teuchos::RCP<BlockGCRODRIter<ScalarType,MV,OP> > gcrodr_iter);
270 
271  // Recycling with Harmonic Ritz Vectors
272  // Computes harmonic eigenpairs of projected matrix created during the priming solve.
273  // The return value is the number of vectors needed to be stored, recycledBlocks or recycledBlocks+1.
274 
275  // HH is the projected problem from the initial cycle of Gmres, it is (at least) of dimension (numBlocks+1)*blockSize x numBlocks.
276  // PP contains the harmonic eigenvectors corresponding to the recycledBlocks eigenvalues of smallest magnitude.
277  int getHarmonicVecsKryl (int m, const SDM& HH, SDM& PP);
278 
279  // HH is the total block projected problem from the GCRO-DR algorithm, it is (at least) of dimension keff+(numBlocks+1)*blockSize x keff+numBlocksm.
280  // VV is the Krylov vectors from the projected GMRES algorithm, which has (at least) (numBlocks+1)*blockSize vectors.
281  // PP contains the harmonic eigenvectors corresponding to the recycledBlocks eigenvalues of smallest magnitude.
282  int getHarmonicVecsAugKryl (int keff,
283  int m,
284  const SDM& HH,
285  const Teuchos::RCP<const MV>& VV,
286  SDM& PP);
287 
288  // Sort list of n floating-point numbers and return permutation vector
289  void sort (std::vector<MagnitudeType>& dlist, int n, std::vector<int>& iperm);
290 
291  // Lapack interface
293 
296 
297  //Output Manager
299  Teuchos::RCP<std::ostream> outputStream_;
300 
301  //Status Test
305  Teuchos::RCP<StatusTestGenResNorm<ScalarType,MV,OP> > expConvTest_, impConvTest_;
307 
309  ortho_factory_type orthoFactory_;
310 
317 
320 
331  mutable Teuchos::RCP<const Teuchos::ParameterList> defaultParams_;
332 
333  //Default Solver Values
334  static const bool adaptiveBlockSize_default_;
335  static const std::string recycleMethod_default_;
336 
337  //Current Solver Values
338  MagnitudeType convTol_, orthoKappa_, achievedTol_;
339  int blockSize_, maxRestarts_, maxIters_, numIters_;
340  int verbosity_, outputStyle_, outputFreq_;
341  bool adaptiveBlockSize_;
342  std::string orthoType_, recycleMethod_;
343  std::string impResScale_, expResScale_;
344  std::string label_;
345 
347  // Solver State Storage
349  //
350  // The number of blocks and recycle blocks (m and k, respectively)
351  int numBlocks_, recycledBlocks_;
352  // Current size of recycled subspace
353  int keff;
354  //
355  // Residual Vector
356  Teuchos::RCP<MV> R_;
357  //
358  // Search Space
359  Teuchos::RCP<MV> V_;
360  //
361  // Recycle subspace and its image under action of the operator
362  Teuchos::RCP<MV> U_, C_;
363  //
364  // Updated recycle Space and its image under action of the operator
365  Teuchos::RCP<MV> U1_, C1_;
366  //
367  // Storage used in constructing recycle space
371  Teuchos::RCP<SDM > PP_;
372  Teuchos::RCP<SDM > HP_;
373  std::vector<ScalarType> tau_;
374  std::vector<ScalarType> work_;
376  std::vector<int> ipiv_;
377 
379  Teuchos::RCP<Teuchos::Time> timerSolve_;
380 
382  bool isSet_;
383 
385  bool loaDetected_;
386 
388  bool builtRecycleSpace_;
389 };
390 
391  //
392  // Set default solver values
393  //
394  template<class ScalarType, class MV, class OP>
395  const bool BlockGCRODRSolMgr<ScalarType,MV,OP>::adaptiveBlockSize_default_ = true;
396 
397  template<class ScalarType, class MV, class OP>
398  const std::string BlockGCRODRSolMgr<ScalarType,MV,OP>::recycleMethod_default_ = "harmvecs";
399 
400  //
401  // Method definitions
402  //
403 
404  template<class ScalarType, class MV, class OP>
406  init();
407  }
408 
409  //Basic Constructor
410  template<class ScalarType, class MV, class OP>
414  // Initialize local pointers to null, and initialize local
415  // variables to default values.
416  init ();
417 
418  TEUCHOS_TEST_FOR_EXCEPTION(problem.is_null(), std::invalid_argument,
419  "Belos::BlockGCRODR constructor: The solver manager's constructor needs "
420  "the linear problem argument 'problem' to be nonnull.");
421 
422  problem_ = problem;
423 
424  // Set the parameters using the list that was passed in. If null, we defer initialization until a non-null list is set (by the
425  // client calling setParameters(), or by calling solve() -- in either case, a null parameter list indicates that default parameters should be used).
426  if (! pl.is_null())
427  setParameters (pl);
428  }
429 
430  template<class ScalarType, class MV, class OP>
432  adaptiveBlockSize_ = adaptiveBlockSize_default_;
433  recycleMethod_ = recycleMethod_default_;
434  isSet_ = false;
435  loaDetected_ = false;
436  builtRecycleSpace_ = false;
437  keff = 0;//Effective Size of Recycle Space
438  //The following are all RCP smart pointers to indicated matrices/vectors.
439  //Some MATLAB notation used in comments.
440  R_ = Teuchos::null;//The Block Residual
441  V_ = Teuchos::null;//Block Arnoldi Vectors
442  U_ = Teuchos::null;//Recycle Space
443  C_ = Teuchos::null;//Image of U Under Action of Operator
444  U1_ = Teuchos::null;//Newly Computed Recycle Space
445  C1_ = Teuchos::null;//Image of Newly Computed Recycle Space
446  PP_ = Teuchos::null;//Coordinates of New Recyc. Vectors in Augmented Space
447  HP_ = Teuchos::null;//H_*PP_ or G_*PP_
448  G_ = Teuchos::null;//G_ such that A*[U V(:,1:numBlocks_*blockSize_)] = [C V_]*G_
449  F_ = Teuchos::null;//Upper Triangular portion of QR factorization for HP_
450  H_ = Teuchos::null;//H_ such that A*V(:,1:numBlocks_*blockSize_) = V_*H_ + C_*C_'*V_
451  B_ = Teuchos::null;//B_ = C_'*V_
452 
453  //THIS BLOCK OF CODE IS COMMENTED OUT AND PLACED ELSEWHERE IN THE CODE
454 /*
455  //WE TEMPORARILY INITIALIZE STATUS TESTS HERE FOR TESTING PURPOSES, BUT THEY SHOULD BE
456  //INITIALIZED SOMEWHERE ELSE, LIKE THE setParameters() FUNCTION
457 
458  //INSTANTIATE AND INITIALIZE TEST OBJECTS AS NEEDED
459  if (maxIterTest_.is_null()){
460  maxIterTest_ = rcp (new StatusTestMaxIters<ScalarType,MV,OP> (maxIters_));
461  }
462  //maxIterTest_->setMaxIters (maxIters_);
463 
464  //INSTANTIATE THE PRINTER
465  if (printer_.is_null()) {
466  printer_ = Teuchos::rcp (new OutputManager<ScalarType> (verbosity_, outputStream_));
467  }
468 
469  if (ortho_.is_null()) // || changedOrthoType) %%%In other codes, this is also triggered if orthogonalization type changed {
470  // Create orthogonalization manager. This requires that the OutputManager (printer_) already be initialized.
471  Teuchos::RCP<const Teuchos::ParameterList> orthoParams = orthoFactory_.getDefaultParameters (orthoType_);
472  ortho_ = orthoFactory_.makeMatOrthoManager (orthoType_, Teuchos::null, printer_, label_, orthoParams);
473  }
474 
475  // Convenience typedefs
476  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
477  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
478 
479  if (impConvTest_.is_null()) {
480  impConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
481  impConvTest_->defineScaleForm (convertStringToScaleType (impResScale_), Belos::TwoNorm);
482  impConvTest_->setShowMaxResNormOnly( true );
483  }
484 
485  if (expConvTest_.is_null()) {
486  expConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
487  expConvTest_->defineResForm (StatusTestResNorm_t::Explicit, Belos::TwoNorm);
488  expConvTest_->defineScaleForm (convertStringToScaleType (expResScale_), Belos::TwoNorm);
489  expConvTest_->setShowMaxResNormOnly( true );
490  }
491 
492  if (convTest_.is_null()) {
493  convTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::SEQ, impConvTest_, expConvTest_));
494  }
495 
496  sTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::OR, maxIterTest_, convTest_));
497 
498  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory (outputStyle_);
499  outputTest_ = stoFactory.create (printer_, sTest_, outputFreq_,
500  Passed+Failed+Undefined);
501 
502  */
503 
504  }
505 
506  // This method requires the solver manager to return a string that describes itself.
507  template<class ScalarType, class MV, class OP>
509  std::ostringstream oss;
510  oss << "Belos::BlockGCRODRSolMgr<" << SCT::name() << ", ...>";
511  oss << "{";
512  oss << "Ortho Type='"<<orthoType_ ;
513  oss << ", Num Blocks=" <<numBlocks_;
514  oss << ", Block Size=" <<blockSize_;
515  oss << ", Num Recycle Blocks=" << recycledBlocks_;
516  oss << ", Max Restarts=" << maxRestarts_;
517  oss << "}";
518  return oss.str();
519  }
520 
521  template<class ScalarType, class MV, class OP>
525  using Teuchos::parameterList;
526  using Teuchos::RCP;
527  using Teuchos::rcpFromRef;
528 
529  if (defaultParams_.is_null()) {
530  RCP<ParameterList> pl = parameterList ();
531 
532  const MagnitudeType convTol = SMT::squareroot (SCT::magnitude (SCT::eps()));
533  const int maxRestarts = 1000;
534  const int maxIters = 5000;
535  const int blockSize = 2;
536  const int numBlocks = 100;
537  const int numRecycledBlocks = 25;
539  const Belos::OutputType outputStyle = Belos::General;
540  const int outputFreq = 1;
541  RCP<std::ostream> outputStream = rcpFromRef (std::cout);
542  const std::string impResScale ("Norm of Preconditioned Initial Residual");
543  const std::string expResScale ("Norm of Initial Residual");
544  const std::string timerLabel ("Belos");
545  const std::string orthoType ("ICGS");
546  RCP<const ParameterList> orthoParams = orthoFactory_.getDefaultParameters (orthoType);
547  //const MagnitudeType orthoKappa = SCT::magnitude (SCT::eps());
548  //
549  // mfh 31 Dec 2011: Negative values of orthoKappa are ignored.
550  // Setting this to a negative value by default ensures that
551  // this parameter is only _not_ ignored if the user
552  // specifically sets a valid value.
553  const MagnitudeType orthoKappa = -SMT::one();
554 
555  // Set all the valid parameters and their default values.
556  pl->set ("Convergence Tolerance", convTol,
557  "The tolerance that the solver needs to achieve in order for "
558  "the linear system(s) to be declared converged. The meaning "
559  "of this tolerance depends on the convergence test details.");
560  pl->set("Maximum Restarts", maxRestarts,
561  "The maximum number of restart cycles (not counting the first) "
562  "allowed for each set of right-hand sides solved.");
563  pl->set("Maximum Iterations", maxIters,
564  "The maximum number of iterations allowed for each set of "
565  "right-hand sides solved.");
566  pl->set("Block Size", blockSize,
567  "How many linear systems to solve at once.");
568  pl->set("Num Blocks", numBlocks,
569  "The maximum number of blocks allowed in the Krylov subspace "
570  "for each set of right-hand sides solved. This includes the "
571  "initial block. Total Krylov basis storage, not counting the "
572  "recycled subspace, is \"Num Blocks\" times \"Block Size\".");
573  pl->set("Num Recycled Blocks", numRecycledBlocks,
574  "The maximum number of vectors in the recycled subspace." );
575  pl->set("Verbosity", verbosity,
576  "What type(s) of solver information should be written "
577  "to the output stream.");
578  pl->set("Output Style", outputStyle,
579  "What style is used for the solver information to write "
580  "to the output stream.");
581  pl->set("Output Frequency", outputFreq,
582  "How often convergence information should be written "
583  "to the output stream.");
584  pl->set("Output Stream", outputStream,
585  "A reference-counted pointer to the output stream where all "
586  "solver output is sent.");
587  pl->set("Implicit Residual Scaling", impResScale,
588  "The type of scaling used in the implicit residual convergence test.");
589  pl->set("Explicit Residual Scaling", expResScale,
590  "The type of scaling used in the explicit residual convergence test.");
591  pl->set("Timer Label", timerLabel,
592  "The string to use as a prefix for the timer labels.");
593  pl->set("Orthogonalization", orthoType,
594  "The orthogonalization method to use. Valid options: " +
595  orthoFactory_.validNamesString());
596  pl->set ("Orthogonalization Parameters", *orthoParams,
597  "Sublist of parameters specific to the orthogonalization method to use.");
598  pl->set("Orthogonalization Constant", orthoKappa,
599  "When using DGKS orthogonalization: the \"depTol\" constant, used "
600  "to determine whether another step of classical Gram-Schmidt is "
601  "necessary. Otherwise ignored. Nonpositive values are ignored.");
602  defaultParams_ = pl;
603  }
604  return defaultParams_;
605  }
606 
607  template<class ScalarType, class MV, class OP>
608  void
611  using Teuchos::isParameterType;
612  using Teuchos::getParameter;
613  using Teuchos::null;
615  using Teuchos::parameterList;
616  using Teuchos::RCP;
617  using Teuchos::rcp;
618  using Teuchos::rcp_dynamic_cast;
619  using Teuchos::rcpFromRef;
623 
624  // The default parameter list contains all parameters that BlockGCRODRSolMgr understands, and none that it doesn't
625  RCP<const ParameterList> defaultParams = getValidParameters();
626 
627  if (params.is_null()) {
628  // Users really should supply a non-null input ParameterList,
629  // so that we can fill it in. However, if they really did give
630  // us a null list, be nice and use default parameters. In this
631  // case, we don't have to do any validation.
632  params_ = parameterList (*defaultParams);
633  }
634  else {
635  // Validate the user's parameter list and fill in any missing parameters with default values.
636  //
637  // Currently, validation is quite strict. The following line
638  // will throw an exception for misspelled or extra parameters.
639  // However, it will _not_ throw an exception if there are
640  // missing parameters: these will be filled in with their
641  // default values. There is additional discussion of other
642  // validation strategies in the comments of this function for
643  // Belos::GCRODRSolMgr.
644  params->validateParametersAndSetDefaults (*defaultParams);
645  // No side effects on the solver until after validation passes.
646  params_ = params;
647  }
648 
649  //
650  // Validate and set parameters relating to the Krylov subspace
651  // dimensions and the number of blocks.
652  //
653  // We try to read and validate all these parameters' values
654  // before setting them in the solver. This is because the
655  // validity of each may depend on the others' values. Our goal
656  // is to avoid incorrect side effects, so we don't set the values
657  // in the solver until we know they are right.
658  //
659  {
660  const int maxRestarts = params_->get<int> ("Maximum Restarts");
661  TEUCHOS_TEST_FOR_EXCEPTION(maxRestarts <= 0, std::invalid_argument,
662  "Belos::BlockGCRODRSolMgr: The \"Maximum Restarts\" parameter "
663  "must be nonnegative, but you specified a negative value of "
664  << maxRestarts << ".");
665 
666  const int maxIters = params_->get<int> ("Maximum Iterations");
667  TEUCHOS_TEST_FOR_EXCEPTION(maxIters <= 0, std::invalid_argument,
668  "Belos::BlockGCRODRSolMgr: The \"Maximum Iterations\" parameter "
669  "must be positive, but you specified a nonpositive value of "
670  << maxIters << ".");
671 
672  const int numBlocks = params_->get<int> ("Num Blocks");
673  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0, std::invalid_argument,
674  "Belos::BlockGCRODRSolMgr: The \"Num Blocks\" parameter must be "
675  "positive, but you specified a nonpositive value of " << numBlocks
676  << ".");
677 
678  const int blockSize = params_->get<int> ("Block Size");
679  TEUCHOS_TEST_FOR_EXCEPTION(blockSize <= 0, std::invalid_argument,
680  "Belos::BlockGCRODRSolMgr: The \"Block Size\" parameter must be "
681  "positive, but you specified a nonpositive value of " << blockSize
682  << ".");
683 
684  const int recycledBlocks = params_->get<int> ("Num Recycled Blocks");
685  TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks <= 0, std::invalid_argument,
686  "Belos::BlockGCRODRSolMgr: The \"Num Recycled Blocks\" parameter must "
687  "be positive, but you specified a nonpositive value of "
688  << recycledBlocks << ".");
689  TEUCHOS_TEST_FOR_EXCEPTION(recycledBlocks >= numBlocks,
690  std::invalid_argument, "Belos::BlockGCRODRSolMgr: The \"Num Recycled "
691  "Blocks\" parameter must be less than the \"Num Blocks\" parameter, "
692  "but you specified \"Num Recycled Blocks\" = " << recycledBlocks
693  << " and \"Num Blocks\" = " << numBlocks << ".");
694 
695  // Now that we've validated the various dimension-related
696  // parameters, we can set them in the solvers.
697  maxRestarts_ = maxRestarts;
698  maxIters_ = maxIters;
699  numBlocks_ = numBlocks;
700  blockSize_ = blockSize;
701  recycledBlocks_ = recycledBlocks;
702  }
703 
704  // Check to see if the timer label changed. If it did, update it in
705  // the parameter list, and create a new timer with that label (if
706  // Belos was compiled with timers enabled).
707  {
708  std::string tempLabel = params_->get<std::string> ("Timer Label");
709  const bool labelChanged = (tempLabel != label_);
710 
711 #ifdef BELOS_TEUCHOS_TIME_MONITOR
712  std::string solveLabel = label_ + ": BlockGCRODRSolMgr total solve time";
713  if (timerSolve_.is_null()) {
714  // We haven't created a timer yet.
715  timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
716  } else if (labelChanged) {
717  // We created a timer before with a different label. In that
718  // case, clear the old timer and create a new timer with the
719  // new label. Clearing old timers prevents them from piling
720  // up, since they are stored statically, possibly without
721  // checking for duplicates.
723  timerSolve_ = Teuchos::TimeMonitor::getNewCounter (solveLabel);
724  }
725 #endif // BELOS_TEUCHOS_TIME_MONITOR
726  }
727 
728  // Check for a change in verbosity level. Just in case, we allow
729  // the type of this parameter to be either int or MsgType (an
730  // enum).
731  if (params_->isParameter ("Verbosity")) {
732  if (isParameterType<int> (*params_, "Verbosity")) {
733  verbosity_ = params_->get<int> ("Verbosity");
734  }
735  else {
736  verbosity_ = (int) getParameter<MsgType> (*params_, "Verbosity");
737  }
738  }
739 
740  // Check for a change in output style. Just in case, we allow
741  // the type of this parameter to be either int or OutputType (an
742  // enum).
743  if (params_->isParameter ("Output Style")) {
744  if (isParameterType<int> (*params_, "Output Style")) {
745  outputStyle_ = params_->get<int> ("Output Style");
746  }
747  else {
748  outputStyle_ = (int) getParameter<OutputType> (*params_, "Output Style");
749  }
750 
751  // Currently, the only way we can update the output style is to
752  // create a new output test. This is because we must first
753  // instantiate a new StatusTestOutputFactory with the new
754  // style, and then use the factory to make a new output test.
755  // Thus, we set outputTest_ to null for now, so we remember
756  // later to create a new one.
757  outputTest_ = null;
758  }
759 
760  // Get the output stream for the output manager.
761  //
762  // It has been pointed out (mfh 28 Feb 2011 in GCRODRSolMgr code)
763  // that including an RCP<std::ostream> parameter makes it
764  // impossible to serialize the parameter list. If you serialize
765  // the list and read it back in from the serialized
766  // representation, you might not even get a valid
767  // RCP<std::ostream>, let alone the same output stream that you
768  // serialized.
769  //
770  // However, existing Belos users depend on setting the output
771  // stream in the parameter list. We retain this behavior for
772  // backwards compatibility.
773  //
774  // In case the output stream can't be read back in, we default to
775  // stdout (std::cout), just to ensure reasonable behavior.
776  if (params_->isParameter ("Output Stream")) {
777  try {
778  outputStream_ = getParameter<RCP<std::ostream> > (*params_, "Output Stream");
779  }
780  catch (InvalidParameter&) {
781  outputStream_ = rcpFromRef (std::cout);
782  }
783  // We assume that a null output stream indicates that the user
784  // doesn't want to print anything, so we replace it with a
785  // "black hole" stream that prints nothing sent to it. (We
786  // can't use outputStream_ = Teuchos::null, since the output
787  // manager assumes that operator<< works on its output stream.)
788  if (outputStream_.is_null()) {
789  outputStream_ = rcp (new Teuchos::oblackholestream);
790  }
791  }
792 
793  outputFreq_ = params_->get<int> ("Output Frequency");
794  if (verbosity_ & Belos::StatusTestDetails) {
795  // Update parameter in our output status test.
796  if (! outputTest_.is_null()) {
797  outputTest_->setOutputFrequency (outputFreq_);
798  }
799  }
800 
801  // Create output manager if we need to, using the verbosity level
802  // and output stream that we fetched above. We do this here because
803  // instantiating an OrthoManager using OrthoManagerFactory requires
804  // a valid OutputManager.
805  if (printer_.is_null()) {
806  printer_ = rcp (new OutputManager<ScalarType> (verbosity_, outputStream_));
807  } else {
808  printer_->setVerbosity (verbosity_);
809  printer_->setOStream (outputStream_);
810  }
811 
812  // Get the orthogonalization manager name ("Orthogonalization").
813  //
814  // Getting default values for the orthogonalization manager
815  // parameters ("Orthogonalization Parameters") requires knowing the
816  // orthogonalization manager name. Save it for later, and also
817  // record whether it's different than before.
818 
819  //bool changedOrthoType = false; // silence warning about not being used
820  if (params_->isParameter ("Orthogonalization")) {
821  const std::string& tempOrthoType =
822  params_->get<std::string> ("Orthogonalization");
823  // Ensure that the specified orthogonalization type is valid.
824  if (! orthoFactory_.isValidName (tempOrthoType)) {
825  std::ostringstream os;
826  os << "Belos::BlockGCRODRSolMgr: Invalid orthogonalization name \""
827  << tempOrthoType << "\". The following are valid options "
828  << "for the \"Orthogonalization\" name parameter: ";
829  orthoFactory_.printValidNames (os);
830  TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, os.str());
831  }
832  if (tempOrthoType != orthoType_) {
833  //changedOrthoType = true;
834  orthoType_ = tempOrthoType;
835  }
836  }
837 
838  // Get any parameters for the orthogonalization
839  // ("Orthogonalization Parameters"). If not supplied, the
840  // orthogonalization manager factory will supply default values.
841  //
842  // NOTE (mfh 12 Jan 2011, summer 2011, 31 Dec 2011) For backwards
843  // compatibility with other Belos GMRES-type solvers, if params
844  // has an "Orthogonalization Constant" parameter and the DGKS
845  // orthogonalization manager is to be used, the value of this
846  // parameter will override DGKS's "depTol" parameter.
847  //
848  // Users must supply the orthogonalization manager parameters as
849  // a sublist (supplying it as an RCP<ParameterList> would make
850  // the resulting parameter list not serializable).
851  RCP<ParameterList> orthoParams = sublist (params, "Orthogonalization Parameters", true);
852  TEUCHOS_TEST_FOR_EXCEPTION(orthoParams.is_null(), std::logic_error,
853  "Failed to get orthogonalization parameters. "
854  "Please report this bug to the Belos developers.");
855 
856  // Check if the desired orthogonalization method changed, or if
857  // the orthogonalization manager has not yet been instantiated.
858  // If either is the case, instantiate a new MatOrthoManager
859  // subclass instance corresponding to the desired
860  // orthogonalization method. We've already fetched the
861  // orthogonalization method name (orthoType_) and its parameters
862  // (orthoParams) above.
863  //
864  // As suggested (by mfh 12 Jan 2011 in Belos::GCRODRSolMgr) In
865  // order to ensure parameter changes get propagated to the
866  // orthomanager we simply reinstantiate the OrthoManager every
867  // time, whether or not the orthogonalization method name or
868  // parameters have changed. This is not efficient for some
869  // orthogonalization managers that keep a lot of internal state.
870  // A more general way to fix this inefficiency would be to
871  //
872  // 1. make each orthogonalization implement
873  // Teuchos::ParameterListAcceptor, and
874  //
875  // 2. make each orthogonalization implement a way to set the
876  // OutputManager instance and timer label.
877  //
878  // Create orthogonalization manager. This requires that the
879  // OutputManager (printer_) already be initialized.
880  ortho_ = orthoFactory_.makeMatOrthoManager (orthoType_, null, printer_,
881  label_, orthoParams);
882 
883  // The DGKS orthogonalization accepts a "Orthogonalization
884  // Constant" parameter (also called kappa in the code, but not in
885  // the parameter list). If its value is provided in the given
886  // parameter list and its value is positive, use it. Ignore
887  // negative values.
888  //
889  // The "Orthogonalization Constant" parameter is retained only
890  // for backwards compatibility.
891 
892  //bool gotValidOrthoKappa = false; // silence warning about not being used
893  if (params_->isParameter ("Orthogonalization Constant")) {
894  const MagnitudeType orthoKappa =
895  params_->get<MagnitudeType> ("Orthogonalization Constant");
896  if (orthoKappa > 0) {
897  orthoKappa_ = orthoKappa;
898  //gotValidOrthoKappa = true;
899  // Only DGKS currently accepts this parameter.
900  if (orthoType_ == "DGKS" && ! ortho_.is_null()) {
901  typedef DGKSOrthoManager<ScalarType, MV, OP> ortho_man_type;
902  // This cast should always succeed; it's a bug otherwise.
903  // (If the cast fails, then orthoType_ doesn't correspond
904  // to the OrthoManager subclass instance that we think we
905  // have, so we initialized the wrong subclass somehow.)
906  rcp_dynamic_cast<ortho_man_type>(ortho_)->setDepTol (orthoKappa_);
907  }
908  }
909  }
910 
911  // Convergence
912  typedef Belos::StatusTestCombo<ScalarType,MV,OP> StatusTestCombo_t;
913  typedef Belos::StatusTestGenResNorm<ScalarType,MV,OP> StatusTestResNorm_t;
914 
915  // Check for convergence tolerance
916  convTol_ = params_->get<MagnitudeType> ("Convergence Tolerance");
917  if (! impConvTest_.is_null()) {
918  impConvTest_->setTolerance (convTol_);
919  }
920  if (! expConvTest_.is_null()) {
921  expConvTest_->setTolerance (convTol_);
922  }
923 
924  // Check for a change in scaling, if so we need to build new residual tests.
925  if (params_->isParameter ("Implicit Residual Scaling")) {
926  std::string tempImpResScale =
927  getParameter<std::string> (*params_, "Implicit Residual Scaling");
928 
929  // Only update the scaling if it's different.
930  if (impResScale_ != tempImpResScale) {
931  ScaleType impResScaleType = convertStringToScaleType (tempImpResScale);
932  impResScale_ = tempImpResScale;
933 
934  if (! impConvTest_.is_null()) {
935  try {
936  impConvTest_->defineScaleForm (impResScaleType, Belos::TwoNorm);
937  }
938  catch (StatusTestError&) {
939  // Delete the convergence test so it gets constructed again.
940  impConvTest_ = null;
941  convTest_ = null;
942  }
943  }
944  }
945  }
946 
947  if (params_->isParameter("Explicit Residual Scaling")) {
948  std::string tempExpResScale =
949  getParameter<std::string> (*params_, "Explicit Residual Scaling");
950 
951  // Only update the scaling if it's different.
952  if (expResScale_ != tempExpResScale) {
953  ScaleType expResScaleType = convertStringToScaleType (tempExpResScale);
954  expResScale_ = tempExpResScale;
955 
956  if (! expConvTest_.is_null()) {
957  try {
958  expConvTest_->defineScaleForm (expResScaleType, Belos::TwoNorm);
959  }
960  catch (StatusTestError&) {
961  // Delete the convergence test so it gets constructed again.
962  expConvTest_ = null;
963  convTest_ = null;
964  }
965  }
966  }
967  }
968 
969  //
970  // Create iteration stopping criteria ("status tests") if we need
971  // to, by combining three different stopping criteria.
972  //
973  // First, construct maximum-number-of-iterations stopping criterion.
974  if (maxIterTest_.is_null()) {
975  maxIterTest_ = rcp (new StatusTestMaxIters<ScalarType,MV,OP> (maxIters_));
976  } else {
977  maxIterTest_->setMaxIters (maxIters_);
978  }
979 
980  // Implicit residual test, using the native residual to determine if
981  // convergence was achieved.
982  if (impConvTest_.is_null()) {
983  impConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
984  impConvTest_->defineScaleForm (convertStringToScaleType (impResScale_),
986  }
987 
988  // Explicit residual test once the native residual is below the tolerance
989  if (expConvTest_.is_null()) {
990  expConvTest_ = rcp (new StatusTestResNorm_t (convTol_));
991  expConvTest_->defineResForm (StatusTestResNorm_t::Explicit, Belos::TwoNorm);
992  expConvTest_->defineScaleForm (convertStringToScaleType (expResScale_),
994  }
995  // Convergence test first tests the implicit residual, then the
996  // explicit residual if the implicit residual test passes.
997  if (convTest_.is_null()) {
998  convTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::SEQ,
999  impConvTest_,
1000  expConvTest_));
1001  }
1002  // Construct the complete iteration stopping criterion:
1003  //
1004  // "Stop iterating if the maximum number of iterations has been
1005  // reached, or if the convergence test passes."
1006  sTest_ = rcp (new StatusTestCombo_t (StatusTestCombo_t::OR,
1007  maxIterTest_, convTest_));
1008  // Create the status test output class.
1009  // This class manages and formats the output from the status test.
1010  StatusTestOutputFactory<ScalarType,MV,OP> stoFactory (outputStyle_);
1011  outputTest_ = stoFactory.create (printer_, sTest_, outputFreq_,
1013 
1014  // Set the solver string for the output test
1015  std::string solverDesc = "Block GCRODR ";
1016  outputTest_->setSolverDesc (solverDesc);
1017 
1018  // Inform the solver manager that the current parameters were set.
1019  isSet_ = true;
1020  }
1021 
1022  // initializeStateStorage.
1023  template<class ScalarType, class MV, class OP>
1024  void
1026  {
1027 
1028  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1029 
1030  // Check if there is any multivector to clone from.
1031  Teuchos::RCP<const MV> rhsMV = problem_->getRHS();
1032 
1033  //The Dimension of the Krylov Subspace
1034  int KrylSpaDim = (numBlocks_ - 1) * blockSize_;
1035 
1036  //Number of columns in [U_ V_(:,1:KrylSpaDim)]
1037  int augSpaDim = KrylSpaDim + recycledBlocks_ + 1;// + 1 is for possible extra recycle vector
1038 
1039  //Number of columns in [C_ V_]
1040  int augSpaImgDim = KrylSpaDim + blockSize_ + recycledBlocks_+1;
1041 
1042  //TEMPORARY SKELETON DEFINITION OF THIS FUNCTION TO GET THINGS WORKING
1043  //NOT EVERYTHING IS INITIALIZE CORRECTLY YET.
1044 
1045  //INITIALIZE RECYCLE SPACE VARIABLES HERE
1046 
1047  //WE DO NOT ALLOCATE V HERE IN THIS SOLVERMANAGER. If no recycle space exists, then block_gmres_iter
1048  //will allocated V for us. If a recycle space already exists, then we will allocate V after updating the
1049  //recycle space for the current problem.
1050  // If the block Krylov subspace has not been initialized before, generate it using the RHS from lp_.
1051  /*if (V_ == Teuchos::null) {
1052  V_ = MVT::Clone( *rhsMV, (numBlocks_+1)*blockSize_ );
1053  }
1054  else{
1055  // Generate V_ by cloning itself ONLY if more space is needed.
1056  if (MVT::GetNumberVecs(*V_) < numBlocks_+1) {
1057  Teuchos::RCP<const MV> tmp = V_;
1058  V_ = MVT::Clone( *tmp, numBlocks_+1 );
1059  }
1060  }*/
1061 
1062  //INTITIALIZE SPACE FOR THE NEWLY COMPUTED RECYCLE SPACE VARIABLES HERE
1063 
1064  if (U_ == Teuchos::null) {
1065  U_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1066  }
1067  else {
1068  // Generate U_ by cloning itself ONLY if more space is needed.
1069  if (MVT::GetNumberVecs(*U_) < recycledBlocks_+1) {
1070  Teuchos::RCP<const MV> tmp = U_;
1071  U_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1072  }
1073  }
1074 
1075  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1076  if (C_ == Teuchos::null) {
1077  C_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1078  }
1079  else {
1080  // Generate C_ by cloning itself ONLY if more space is needed.
1081  if (MVT::GetNumberVecs(*C_) < recycledBlocks_+1) {
1082  Teuchos::RCP<const MV> tmp = C_;
1083  C_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1084  }
1085  }
1086 
1087  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1088  if (U1_ == Teuchos::null) {
1089  U1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1090  }
1091  else {
1092  // Generate U1_ by cloning itself ONLY if more space is needed.
1093  if (MVT::GetNumberVecs(*U1_) < recycledBlocks_+1) {
1094  Teuchos::RCP<const MV> tmp = U1_;
1095  U1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1096  }
1097  }
1098 
1099  // If the subspace has not been initialized before, generate it using the RHS from lp_.
1100  if (C1_ == Teuchos::null) {
1101  C1_ = MVT::Clone( *rhsMV, recycledBlocks_+1 );
1102  }
1103  else {
1104  // Generate C1_ by cloning itself ONLY if more space is needed.
1105  if (MVT::GetNumberVecs(*U1_) < recycledBlocks_+1) {
1106  Teuchos::RCP<const MV> tmp = C1_;
1107  C1_ = MVT::Clone( *tmp, recycledBlocks_+1 );
1108  }
1109  }
1110 
1111  // Generate R_ only if it doesn't exist
1112  if (R_ == Teuchos::null){
1113  R_ = MVT::Clone( *rhsMV, blockSize_ );
1114  }
1115 
1116  //INITIALIZE SOME WORK VARIABLES
1117 
1118  // Generate G_ only if it doesn't exist, otherwise resize it.
1119  if (G_ == Teuchos::null){
1120  G_ = Teuchos::rcp( new SDM( augSpaImgDim, augSpaDim ) );
1121  }
1122  else{
1123  if ( (G_->numRows() != augSpaImgDim) || (G_->numCols() != augSpaDim) )
1124  {
1125  G_->reshape( augSpaImgDim, augSpaDim );
1126  }
1127  G_->putScalar(zero);
1128  }
1129 
1130  // Generate H_ only if it doesn't exist by pointing it to a view of G_.
1131  if (H_ == Teuchos::null){
1132  H_ = Teuchos::rcp (new SDM ( Teuchos::View, *G_, KrylSpaDim + blockSize_, KrylSpaDim, recycledBlocks_+1 ,recycledBlocks_+1 ) );
1133  }
1134 
1135  // Generate F_ only if it doesn't exist, otherwise resize it.
1136  if (F_ == Teuchos::null){
1137  F_ = Teuchos::rcp( new SDM( recycledBlocks_+1, recycledBlocks_+1 ) );
1138  }
1139  else {
1140  if ( (F_->numRows() != recycledBlocks_+1) || (F_->numCols() != recycledBlocks_+1) ){
1141  F_->reshape( recycledBlocks_+1, recycledBlocks_+1 );
1142  }
1143  }
1144  F_->putScalar(zero);
1145 
1146  // Generate PP_ only if it doesn't exist, otherwise resize it.
1147  if (PP_ == Teuchos::null){
1148  PP_ = Teuchos::rcp( new SDM( augSpaImgDim, recycledBlocks_+1 ) );
1149  }
1150  else{
1151  if ( (PP_->numRows() != augSpaImgDim) || (PP_->numCols() != recycledBlocks_+1) ){
1152  PP_->reshape( augSpaImgDim, recycledBlocks_+1 );
1153  }
1154  }
1155 
1156  // Generate HP_ only if it doesn't exist, otherwise resize it.
1157  if (HP_ == Teuchos::null)
1158  HP_ = Teuchos::rcp( new SDM( augSpaImgDim, augSpaDim ) );
1159  else{
1160  if ( (HP_->numRows() != augSpaImgDim) || (HP_->numCols() != augSpaDim) ){
1161  HP_->reshape( augSpaImgDim, augSpaDim );
1162  }
1163  }
1164 
1165  // Size of tau_ will change during computation, so just be sure it starts with appropriate size
1166  tau_.resize(recycledBlocks_+1);
1167 
1168  // Size of work_ will change during computation, so just be sure it starts with appropriate size
1169  work_.resize(recycledBlocks_+1);
1170 
1171  // Size of ipiv_ will change during computation, so just be sure it starts with appropriate size
1172  ipiv_.resize(recycledBlocks_+1);
1173 
1174  }
1175 
1176 template<class ScalarType, class MV, class OP>
1177 void BlockGCRODRSolMgr<ScalarType,MV,OP>::buildRecycleSpaceKryl(int& keff, Teuchos::RCP<BlockGmresIter<ScalarType,MV,OP> > block_gmres_iter){
1178 
1179  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1180  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1181 
1182  int p = block_gmres_iter->getState().curDim; //Dimension of the Krylov space generated
1183  std::vector<int> index(keff);//we use this to index certain columns of U, C, and V to
1184  //get views into pieces of these matrices.
1185 
1186  //GET CORRECT PIECE OF MATRIX H APPROPRIATE TO SIZE OF KRYLOV SUBSPACE
1187  SDM HH(Teuchos::Copy, *H_, p+blockSize_, p);
1188  if(recycledBlocks_ >= p + blockSize_){//keep whole block Krylov subspace
1189  //IF THIS HAS HAPPENED, THIS MEANS WE CONVERGED DURING THIS CYCLE
1190  //THEREFORE, WE DO NOT CONSTRUCT C = A*U;
1191  index.resize(p);
1192  for (int ii=0; ii < p; ++ii) index[ii] = ii;
1193  Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1194  MVT::SetBlock(*V_, index, *Utmp);
1195  keff = p;
1196  }
1197  else{ //use a subspace selection method to get recycle space
1198  int info = 0;
1199  Teuchos::RCP<SDM > PPtmp = Teuchos::rcp (new SDM ( Teuchos::View, *PP_, p, recycledBlocks_+1 ) );
1200  if(recycleMethod_ == "harmvecs"){
1201  keff = getHarmonicVecsKryl(p, HH, *PPtmp);
1202  printer_->stream(Debug) << "keff = " << keff << std::endl;
1203  }
1204 // Hereafter, only keff columns of PP are needed
1205 PPtmp = Teuchos::rcp (new SDM ( Teuchos::View, *PP_, p, keff ) );
1206 // Now get views into C, U, V
1207 index.resize(keff);
1208 for (int ii=0; ii<keff; ++ii) index[ii] = ii;
1209 Teuchos::RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1210 Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1211 Teuchos::RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1212 index.resize(p);
1213 for (int ii=0; ii < p; ++ii) index[ii] = ii;
1214 Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1215 
1216 // Form U (the subspace to recycle)
1217 // U = newstate.V(:,1:p) * PP;
1218 MVT::MvTimesMatAddMv( one, *Vtmp, *PPtmp, zero, *U1tmp );
1219 
1220 // Step #1: Form HP = H*P
1221 SDM HPtmp( Teuchos::View, *HP_, p+blockSize_, keff );
1222 HPtmp.multiply( Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, *H_, *PPtmp, zero );
1223 // Step #1.5: Perform workspace size query for QR factorization of HP
1224 int lwork = -1;
1225 tau_.resize(keff);
1226 lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1227 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure, "Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a workspace size.");
1228 
1229 // Step #2: Compute QR factorization of HP
1230  lwork = static_cast<int> (Teuchos::ScalarTraits<ScalarType>::magnitude (work_[0]));
1231 work_.resize(lwork);
1232 lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1233 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure, "Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a QR factorization.");
1234 
1235 // Step #3: Explicitly construct Q and R factors
1236 // The upper triangular part of HP is copied into R and HP becomes Q.
1237 SDM Rtmp( Teuchos::View, *F_, keff, keff );
1238 for(int ii=0;ii<keff;ii++) { for(int jj=ii;jj<keff;jj++) Rtmp(ii,jj) = HPtmp(ii,jj); }
1239 //lapack.ORGQR(HPtmp.numRows(),HPtmp.numCols(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1240 lapack.UNGQR(HPtmp.numRows(),HPtmp.numCols(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1241 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure, "Belos::BlockGCRODRSolMgr::solve(): LAPACK _UNGQR failed to construct the Q factor.");
1242  // Now we have [Q,R] = qr(H*P)
1243 
1244  // Now compute C = V(:,1:p+blockSize_) * Q
1245  index.resize( p+blockSize_ );
1246  for (int ii=0; ii < (p+blockSize_); ++ii) { index[ii] = ii; }
1247  Vtmp = MVT::CloneView( *V_, index ); // need new view into V (p+blockSize_ vectors now; needed p above)
1248  MVT::MvTimesMatAddMv( one, *Vtmp, HPtmp, zero, *Ctmp );
1249 
1250  // Finally, compute U = U*R^{-1}.
1251  // This unfortuntely requires me to form R^{-1} explicitly and execute U = U * R^{-1}, as
1252  // backsolve capabilities don't exist in the Belos::MultiVec class
1253 
1254 // Step #1: First, compute LU factorization of R
1255 ipiv_.resize(Rtmp.numRows());
1256 lapack.GETRF(Rtmp.numRows(),Rtmp.numCols(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&info);
1257 TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1258 // Step #2: Form inv(R)
1259 lwork = Rtmp.numRows();
1260 work_.resize(lwork);
1261 lapack.GETRI(Rtmp.numRows(),Rtmp.values(),Rtmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1262 //TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure, "Belos::GCRODRSolMgr::solve(): LAPACK _GETRI failed to invert triangular matrix.");
1263 // Step #3: Let U = U * R^{-1}
1264 MVT::MvTimesMatAddMv( one, *U1tmp, Rtmp, zero, *Utmp );
1265 
1266 } //end else from if(recycledBlocks_ >= p + 1)
1267 return;
1268 } // end buildRecycleSpaceKryl defnition
1269 
1270 template<class ScalarType, class MV, class OP>
1271 void BlockGCRODRSolMgr<ScalarType,MV,OP>::buildRecycleSpaceAugKryl(Teuchos::RCP<BlockGCRODRIter<ScalarType,MV,OP> > block_gcrodr_iter){
1272  const MagnitudeType one = Teuchos::ScalarTraits<ScalarType>::one();
1273  const ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1274 
1275  std::vector<MagnitudeType> d(keff);
1276  std::vector<ScalarType> dscalar(keff);
1277  std::vector<int> index(numBlocks_+1);
1278 
1279  // Get the state
1280  BlockGCRODRIterState<ScalarType,MV> oldState = block_gcrodr_iter->getState();
1281  int p = oldState.curDim;
1282 
1283  // insufficient new information to update recycle space
1284  if (p<1) return;
1285 
1286  if(recycledBlocks_ >= keff + p){ //we add new Krylov vectors to existing recycle space
1287  // If this has happened, we have converged during this cycle. Therefore, do not construct C = A*C;
1288  index.resize(p);
1289  for (int ii=0; ii < p; ++ii) { index[ii] = keff+ii; } //get a view after current reycle vectors
1290  Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1291  for (int ii=0; ii < p; ++ii) { index[ii] =ii; }
1292  MVT::SetBlock(*V_, index, *Utmp);
1293  keff += p;
1294  }
1295 
1296  // Take the norm of the recycled vectors
1297  {
1298  index.resize(keff);
1299  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1300  Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1301  d.resize(keff);
1302  dscalar.resize(keff);
1303  MVT::MvNorm( *Utmp, d );
1304  for (int i=0; i<keff; ++i) {
1305  d[i] = one / d[i];
1306  dscalar[i] = (ScalarType)d[i];
1307  }
1308  MVT::MvScale( *Utmp, dscalar );
1309  }
1310 
1311  // Get view into current "full" upper Hessnburg matrix
1312  // note that p describes the dimension of the iter+1 block Krylov space so we have to adjust how we use it to index Gtmp
1313  Teuchos::RCP<SDM> Gtmp = Teuchos::rcp( new SDM( Teuchos::View, *G_, p+keff, p+keff-blockSize_ ) );
1314 
1315  // Insert D into the leading keff x keff block of G
1316  for (int i=0; i<keff; ++i)
1317  (*Gtmp)(i,i) = d[i];
1318 
1319  // Compute the harmoic Ritz pairs for the generalized eigenproblem
1320  // getHarmonicVecsKryl assumes PP has recycledBlocks_+1 columns available
1321  // See previous block of comments for why we subtract p-blockSize_
1322  int keff_new;
1323  {
1324  SDM PPtmp( Teuchos::View, *PP_, p+keff-blockSize_, recycledBlocks_+1 );
1325  keff_new = getHarmonicVecsAugKryl( keff, p-blockSize_, *Gtmp, oldState.V, PPtmp );
1326  }
1327 
1328  // Code to form new U, C
1329  // U = [U V(:,1:p)] * P; (in two steps)
1330 
1331  // U(:,1:keff) = matmul(U(:,1:keff_old),PP(1:keff_old,1:keff)) (step 1)
1332  Teuchos::RCP<MV> U1tmp;
1333  {
1334  index.resize( keff );
1335  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1336  Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1337  index.resize( keff_new );
1338  for (int ii=0; ii<keff_new; ++ii) { index[ii] = ii; }
1339  U1tmp = MVT::CloneViewNonConst( *U1_, index );
1340  SDM PPtmp( Teuchos::View, *PP_, keff, keff_new );
1341  MVT::MvTimesMatAddMv( one, *Utmp, PPtmp, zero, *U1tmp );
1342  }
1343 
1344  // U(:,1:keff) = U(:,1:keff) + matmul(V(:,1:m-k),PP(keff_old+1:m-k+keff_old,1:keff)) (step 2)
1345  {
1346  index.resize(p-blockSize_);
1347  for (int ii=0; ii < p-blockSize_; ii++) { index[ii] = ii; }
1348  Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1349  SDM PPtmp( Teuchos::View, *PP_, p-blockSize_, keff_new, keff );
1350  MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *U1tmp );
1351  }
1352 
1353  // Form GP = G*P
1354  SDM HPtmp( Teuchos::View, *HP_, p+keff, keff_new );
1355  {
1356  SDM PPtmp( Teuchos::View, *PP_, p-blockSize_+keff, keff_new );
1357  HPtmp.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,one,*Gtmp,PPtmp,zero);
1358  }
1359 
1360  // Workspace size query for QR factorization HP= QF (the worksize will be placed in work_[0])
1361  int info = 0, lwork = -1;
1362  tau_.resize(keff_new);
1363  lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1364  TEUCHOS_TEST_FOR_EXCEPTION(info != 0,BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a workspace size.");
1365 
1366  lwork = static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work_[0]));
1367  work_.resize(lwork);
1368  lapack.GEQRF(HPtmp.numRows(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1369  TEUCHOS_TEST_FOR_EXCEPTION(info != 0,BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GEQRF failed to compute a QR factorization.");
1370 
1371  // Explicitly construct Q and F factors
1372  // NOTE: The upper triangular part of HP is copied into F and HP becomes Q.
1373  SDM Ftmp( Teuchos::View, *F_, keff_new, keff_new );
1374  for(int i=0;i<keff_new;i++) { for(int j=i;j<keff_new;j++) Ftmp(i,j) = HPtmp(i,j); }
1375  //lapack.ORGQR(HPtmp.numRows(),HPtmp.numCols(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1376  lapack.UNGQR(HPtmp.numRows(),HPtmp.numCols(),HPtmp.numCols(),HPtmp.values(),HPtmp.stride(),&tau_[0],&work_[0],lwork,&info);
1377  TEUCHOS_TEST_FOR_EXCEPTION(info != 0,BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK _UNGQR failed to construct the Q factor.");
1378 
1379  // Form orthonormalized C and adjust U accordingly so that C = A*U
1380  // C = [C V] * Q;
1381 
1382  // C(:,1:keff) = matmul(C(:,1:keff_old),QQ(1:keff_old,1:keff))
1383  {
1384  Teuchos::RCP<MV> C1tmp;
1385  {
1386  index.resize(keff);
1387  for (int i=0; i < keff; i++) { index[i] = i; }
1388  Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
1389  index.resize(keff_new);
1390  for (int i=0; i < keff_new; i++) { index[i] = i; }
1391  C1tmp = MVT::CloneViewNonConst( *C1_, index );
1392  SDM PPtmp( Teuchos::View, *HP_, keff, keff_new );
1393  MVT::MvTimesMatAddMv( one, *Ctmp, PPtmp, zero, *C1tmp );
1394  }
1395  // Now compute C += V(:,1:p+1) * Q
1396  {
1397  index.resize( p );
1398  for (int i=0; i < p; ++i) { index[i] = i; }
1399  Teuchos::RCP<const MV> Vtmp = MVT::CloneView( *V_, index );
1400  SDM PPtmp( Teuchos::View, *HP_, p, keff_new, keff, 0 );
1401  MVT::MvTimesMatAddMv( one, *Vtmp, PPtmp, one, *C1tmp );
1402  }
1403  }
1404 
1405  // C_ = C1_; (via a swap)
1406  std::swap(C_, C1_);
1407 
1408  // Finally, compute U_ = U_*R^{-1}
1409  // First, compute LU factorization of R
1410  ipiv_.resize(Ftmp.numRows());
1411  lapack.GETRF(Ftmp.numRows(),Ftmp.numCols(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&info);
1412  TEUCHOS_TEST_FOR_EXCEPTION(info != 0,BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1413 
1414  // Now, form inv(R)
1415  lwork = Ftmp.numRows();
1416  work_.resize(lwork);
1417  lapack.GETRI(Ftmp.numRows(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1418  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRI failed to compute an LU factorization.");
1419 
1420  {
1421  index.resize(keff_new);
1422  for (int i=0; i < keff_new; i++) { index[i] = i; }
1423  Teuchos::RCP<MV> Utmp = MVT::CloneViewNonConst( *U_, index );
1424  MVT::MvTimesMatAddMv( one, *U1tmp, Ftmp, zero, *Utmp );
1425  }
1426 
1427  // Set the current number of recycled blocks and subspace
1428  // dimension with the Block GCRO-DR iteration.
1429  if (keff != keff_new) {
1430  keff = keff_new;
1431  block_gcrodr_iter->setSize( keff, numBlocks_ );
1432  // Important to zero this out before next cyle
1433  SDM b1( Teuchos::View, *G_, recycledBlocks_+2, 1, 0, recycledBlocks_ );
1434  b1.putScalar(zero);
1435  }
1436 
1437 } //end buildRecycleSpaceAugKryl definition
1438 
1439 template<class ScalarType, class MV, class OP>
1440 int BlockGCRODRSolMgr<ScalarType,MV,OP>::getHarmonicVecsAugKryl(int keff, int m, const SDM& GG, const Teuchos::RCP<const MV>& VV, SDM& PP){
1441  int i, j;
1442  int m2 = GG.numCols();
1443  bool xtraVec = false;
1444  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1445  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1446  std::vector<int> index;
1447 
1448  // Real and imaginary eigenvalue components
1449  std::vector<MagnitudeType> wr(m2), wi(m2);
1450 
1451  // Magnitude of harmonic Ritz values
1452  std::vector<MagnitudeType> w(m2);
1453 
1454  // Real and imaginary (right) eigenvectors; Don't zero out matrix when constructing
1455  SDM vr(m2,m2,false);
1456 
1457  // Sorted order of harmonic Ritz values
1458  std::vector<int> iperm(m2);
1459 
1460  // Set flag indicating recycle space has been generated this solve
1461  builtRecycleSpace_ = true;
1462 
1463  // Form matrices for generalized eigenproblem
1464 
1465  // B = G' * G; Don't zero out matrix when constructing
1466  SDM B(m2,m2,false);
1467  B.multiply(Teuchos::TRANS,Teuchos::NO_TRANS,one,GG,GG,zero);
1468 
1469  // A_tmp = | C'*U 0 |
1470  // | V_{m+1}'*U I |
1471  SDM A_tmp( keff+m+blockSize_, keff+m );
1472 
1473 
1474  // A_tmp(1:keff,1:keff) = C' * U;
1475  index.resize(keff);
1476  for (int i=0; i<keff; ++i) { index[i] = i; }
1477  Teuchos::RCP<const MV> Ctmp = MVT::CloneView( *C_, index );
1478  Teuchos::RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1479  SDM A11( Teuchos::View, A_tmp, keff, keff );
1480  MVT::MvTransMv( one, *Ctmp, *Utmp, A11 );
1481 
1482  // A_tmp(keff+1:m-k+keff+1,1:keff) = V' * U;
1483  SDM A21( Teuchos::View, A_tmp, m+blockSize_, keff, keff );
1484  index.resize(m+blockSize_);
1485  for (i=0; i < m+blockSize_; i++) { index[i] = i; }
1486  Teuchos::RCP<const MV> Vp = MVT::CloneView( *VV, index );
1487  MVT::MvTransMv( one, *Vp, *Utmp, A21 );
1488 
1489  // A_tmp(keff+1:m-k+keff,keff+1:m-k+keff) = eye(m-k);
1490  for( i=keff; i<keff+m; i++)
1491  A_tmp(i,i) = one;
1492 
1493  // A = G' * A_tmp;
1494  SDM A( m2, A_tmp.numCols() );
1495  A.multiply( Teuchos::TRANS, Teuchos::NO_TRANS, one, GG, A_tmp, zero );
1496 
1497  // Compute k smallest harmonic Ritz pairs
1498  // SUBROUTINE DGGEVX( BALANC, JOBVL, JOBVR, SENSE, N, A, LDA, B, LDB,
1499  // ALPHAR, ALPHAI, BETA, VL, LDVL, VR, LDVR, ILO,
1500  // IHI, LSCALE, RSCALE, ABNRM, BBNRM, RCONDE,
1501  // RCONDV, WORK, LWORK, IWORK, BWORK, INFO )
1502  // MLP: 'SCALING' in DGGEVX generates incorrect eigenvalues. Therefore, only permuting
1503  char balanc='P', jobvl='N', jobvr='V', sense='N';
1504  int ld = A.numRows();
1505  int lwork = 6*ld;
1506  int ldvl = ld, ldvr = ld;
1507  int info = 0,ilo = 0,ihi = 0;
1508  MagnitudeType abnrm = 0.0, bbnrm = 0.0;
1509  ScalarType *vl = 0; // This is never referenced by dggevx if jobvl == 'N'
1510  std::vector<ScalarType> beta(ld);
1511  std::vector<ScalarType> work(lwork);
1512  std::vector<MagnitudeType> rwork(lwork);
1513  std::vector<MagnitudeType> lscale(ld), rscale(ld);
1514  std::vector<MagnitudeType> rconde(ld), rcondv(ld);
1515  std::vector<int> iwork(ld+6);
1516  int *bwork = 0; // If sense == 'N', bwork is never referenced
1517  //lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
1518  // &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
1519  // &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &iwork[0], bwork, &info);
1520  lapack.GGEVX(balanc, jobvl, jobvr, sense, ld, A.values(), ld, B.values(), ld, &wr[0], &wi[0],
1521  &beta[0], vl, ldvl, vr.values(), ldvr, &ilo, &ihi, &lscale[0], &rscale[0],
1522  &abnrm, &bbnrm, &rconde[0], &rcondv[0], &work[0], lwork, &rwork[0],
1523  &iwork[0], bwork, &info);
1524  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure, "Belos::BlockGCRODRSolMgr::solve(): LAPACK GGEVX failed to compute eigensolutions.");
1525 
1526  // Construct magnitude of each harmonic Ritz value
1527  // NOTE : Forming alpha/beta *should* be okay here, given assumptions on construction of matrix pencil above
1528  for( i=0; i<ld; i++ ) // Construct magnitude of each harmonic Ritz value
1529  w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] ) / std::abs( beta[i] );
1530 
1531  this->sort(w,ld,iperm);
1532 
1533  bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
1534 
1535  // Select recycledBlocks_ smallest eigenvectors
1536  for( i=0; i<recycledBlocks_; i++ )
1537  for( j=0; j<ld; j++ )
1538  PP(j,i) = vr(j,iperm[ld-recycledBlocks_+i]);
1539 
1540  if(scalarTypeIsComplex==false) {
1541 
1542  // Determine exact size for PP (i.e., determine if we need to store an additional vector)
1543  if (wi[iperm[ld-recycledBlocks_]] != 0.0) {
1544  int countImag = 0;
1545  for ( i=ld-recycledBlocks_; i<ld; i++ )
1546  if (wi[iperm[i]] != 0.0) countImag++;
1547  // Check to see if this count is even or odd:
1548  if (countImag % 2) xtraVec = true;
1549  }
1550 
1551  if (xtraVec) { // we need to store one more vector
1552  if (wi[iperm[ld-recycledBlocks_]] > 0.0) { // I picked the "real" component
1553  for( j=0; j<ld; j++ ) // so get the "imag" component
1554  PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]+1);
1555  }
1556  else { // I picked the "imag" component
1557  for( j=0; j<ld; j++ ) // so get the "real" component
1558  PP(j,recycledBlocks_) = vr(j,iperm[ld-recycledBlocks_]-1);
1559  }
1560  }
1561 
1562  }
1563 
1564  // Return whether we needed to store an additional vector
1565  if (xtraVec)
1566  return recycledBlocks_+1;
1567  else
1568  return recycledBlocks_;
1569 
1570 } //end getHarmonicVecsAugKryl definition
1571 
1572 template<class ScalarType, class MV, class OP>
1573 int BlockGCRODRSolMgr<ScalarType,MV,OP>::getHarmonicVecsKryl(int m, const SDM& HH, SDM& PP){
1574  bool xtraVec = false;
1575  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1576  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1577 
1578  // Real and imaginary eigenvalue components
1579  std::vector<MagnitudeType> wr(m), wi(m);
1580 
1581  // Real and imaginary (right) eigenvectors; Don't zero out matrix when constructing
1582  SDM vr(m,m,false);
1583 
1584  // Magnitude of harmonic Ritz values
1585  std::vector<MagnitudeType> w(m);
1586 
1587  // Sorted order of harmonic Ritz values, also used for DGEEV
1588  std::vector<int> iperm(m);
1589 
1590  // Output info
1591  int info = 0;
1592 
1593  // Set flag indicating recycle space has been generated this solve
1594  builtRecycleSpace_ = true;
1595 
1596  // Solve linear system: H_m^{-H}*E_m where E_m is the last blockSize_ columns of the identity matrix
1597  SDM HHt( HH, Teuchos::TRANS );
1598  Teuchos::RCP<SDM> harmRitzMatrix = Teuchos::rcp( new SDM( m, blockSize_));
1599 
1600  //Initialize harmRitzMatrix as E_m
1601  for(int i=0; i<=blockSize_-1; i++) (*harmRitzMatrix)[blockSize_-1-i][harmRitzMatrix->numRows()-1-i] = 1;
1602 
1603  //compute harmRitzMatrix <- H_m^{-H}*E_m
1604  lapack.GESV(m, blockSize_, HHt.values(), HHt.stride(), &iperm[0], harmRitzMatrix->values(), harmRitzMatrix->stride(), &info);
1605 
1606  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure, "Belos::BlockGCRODRSolMgr::solve(): LAPACK GESV failed to compute a solution.");
1607  // Compute H_m + H_m^{-H}*E_m*H_lbl^{H}*H_lbl
1608  // H_lbl is bottom-right block of H_, which is a blockSize_ x blockSize_ matrix
1609 
1610  Teuchos::SerialDenseMatrix<int, ScalarType> H_lbl(Teuchos::View, HH, blockSize_, blockSize_, (HH).numRows()-blockSize_, (HH).numCols()-blockSize_ );
1612 
1613  { //So that HTemp will fall out of scope
1614 
1615  // HH_lbl_t <- H_lbl_t*H_lbl
1616  Teuchos::RCP<SDM> Htemp = Teuchos::null;
1617  Htemp = Teuchos::rcp(new SDM(H_lbl_t.numRows(), H_lbl_t.numCols()));
1618  Htemp -> multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, H_lbl_t, H_lbl, zero);
1619  H_lbl_t.assign(*Htemp);
1620  // harmRitzMatrix <- harmRitzMatrix*HH_lbl_t
1621  Htemp = Teuchos::rcp(new SDM(harmRitzMatrix -> numRows(), harmRitzMatrix -> numCols()));
1622  Htemp -> multiply(Teuchos::NO_TRANS, Teuchos::NO_TRANS, one, *harmRitzMatrix, H_lbl_t, zero);
1623  harmRitzMatrix -> assign(*Htemp);
1624 
1625  // We need to add harmRitzMatrix to the last blockSize_ columns of HH and store the result harmRitzMatrix
1626  int harmColIndex, HHColIndex;
1627  Htemp = Teuchos::rcp(new SDM(Teuchos::Copy,HH,HH.numRows()-blockSize_,HH.numCols()));
1628  for(int i = 0; i<blockSize_; i++){
1629  harmColIndex = harmRitzMatrix -> numCols() - i -1;
1630  HHColIndex = m-i-1;
1631  for(int j=0; j<m; j++) (*Htemp)[HHColIndex][j] += (*harmRitzMatrix)[harmColIndex][j];
1632  }
1633  harmRitzMatrix = Htemp;
1634  }
1635 
1636  // Revise to do query for optimal workspace first
1637  // Create simple storage for the left eigenvectors, which we don't care about.
1638 
1639  // Size of workspace and workspace for DGEEV
1640  int lwork = -1;
1641  std::vector<ScalarType> work(1);
1642  std::vector<MagnitudeType> rwork(2*m);
1643 
1644  const int ldvl = 1;
1645  ScalarType* vl = 0;
1646 
1647  // First query GEEV for the optimal workspace size
1648  lapack.GEEV('N', 'V', m, harmRitzMatrix->values(), harmRitzMatrix->stride(), &wr[0], &wi[0],
1649  vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
1650 
1651  lwork = std::abs (static_cast<int> (Teuchos::ScalarTraits<ScalarType>::real (work[0])));
1652  work.resize( lwork );
1653 
1654  lapack.GEEV('N', 'V', m, harmRitzMatrix->values(), harmRitzMatrix->stride(), &wr[0], &wi[0],
1655  vl, ldvl, vr.values(), vr.stride(), &work[0], lwork, &rwork[0], &info);
1656 
1657  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK GEEV failed to compute eigensolutions.");
1658 
1659  // Construct magnitude of each harmonic Ritz value
1660  for( int i=0; i<m; ++i ) w[i] = Teuchos::ScalarTraits<MagnitudeType>::squareroot( wr[i]*wr[i] + wi[i]*wi[i] );
1661 
1662  this->sort(w, m, iperm);
1663 
1664  bool scalarTypeIsComplex = Teuchos::ScalarTraits<ScalarType>::isComplex;
1665 
1666  // Select recycledBlocks_ smallest eigenvectors
1667  for( int i=0; i<recycledBlocks_; ++i )
1668  for(int j=0; j<m; j++ )
1669  PP(j,i) = vr(j,iperm[i]);
1670 
1671  if(scalarTypeIsComplex==false) {
1672 
1673  // Determine exact size for PP (i.e., determine if we need to store an additional vector)
1674  if (wi[iperm[recycledBlocks_-1]] != 0.0) {
1675  int countImag = 0;
1676  for (int i=0; i<recycledBlocks_; ++i )
1677  if (wi[iperm[i]] != 0.0) countImag++;
1678  // Check to see if this count is even or odd:
1679  if (countImag % 2) xtraVec = true;
1680  }
1681 
1682  if (xtraVec) { // we need to store one more vector
1683  if (wi[iperm[recycledBlocks_-1]] > 0.0) { // I picked the "real" component
1684  for(int j=0; j<m; ++j ) // so get the "imag" component
1685  PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]+1);
1686  }
1687  else{ // I picked the "imag" component
1688  for(int j=0; j<m; ++j ) // so get the "real" component
1689  PP(j,recycledBlocks_) = vr(j,iperm[recycledBlocks_-1]-1);
1690  }
1691  }
1692 
1693  }
1694 
1695  // Return whether we needed to store an additional vector
1696  if (xtraVec) {
1697  printer_->stream(Debug) << "Recycled " << recycledBlocks_+1 << " vectors" << std::endl;
1698  return recycledBlocks_+1;
1699  }
1700  else {
1701  printer_->stream(Debug) << "Recycled " << recycledBlocks_ << " vectors" << std::endl;
1702  return recycledBlocks_;
1703  }
1704 } //end getHarmonicVecsKryl
1705 
1706 // This method sorts list of n floating-point numbers and return permutation vector
1707 template<class ScalarType, class MV, class OP>
1708 void BlockGCRODRSolMgr<ScalarType,MV,OP>::sort(std::vector<MagnitudeType>& dlist, int n, std::vector<int>& iperm) {
1709  int l, r, j, i, flag;
1710  int RR2;
1711  MagnitudeType dRR, dK;
1712 
1713  // Initialize the permutation vector.
1714  for(j=0;j<n;j++)
1715  iperm[j] = j;
1716 
1717  if (n <= 1) return;
1718 
1719  l = n / 2 + 1;
1720  r = n - 1;
1721  l = l - 1;
1722  dRR = dlist[l - 1];
1723  dK = dlist[l - 1];
1724 
1725  RR2 = iperm[l - 1];
1726  while (r != 0) {
1727  j = l;
1728  flag = 1;
1729  while (flag == 1) {
1730  i = j;
1731  j = j + j;
1732  if (j > r + 1)
1733  flag = 0;
1734  else {
1735  if (j < r + 1)
1736  if (dlist[j] > dlist[j - 1]) j = j + 1;
1737  if (dlist[j - 1] > dK) {
1738  dlist[i - 1] = dlist[j - 1];
1739  iperm[i - 1] = iperm[j - 1];
1740  }
1741  else {
1742  flag = 0;
1743  }
1744  }
1745  }
1746  dlist[i - 1] = dRR;
1747  iperm[i - 1] = RR2;
1748  if (l == 1) {
1749  dRR = dlist [r];
1750  RR2 = iperm[r];
1751  dK = dlist[r];
1752  dlist[r] = dlist[0];
1753  iperm[r] = iperm[0];
1754  r = r - 1;
1755  }
1756  else {
1757  l = l - 1;
1758  dRR = dlist[l - 1];
1759  RR2 = iperm[l - 1];
1760  dK = dlist[l - 1];
1761  }
1762  }
1763  dlist[0] = dRR;
1764  iperm[0] = RR2;
1765 } //end sort() definition
1766 
1767 template<class ScalarType, class MV, class OP>
1769  using Teuchos::RCP;
1770  using Teuchos::rcp;
1771  using Teuchos::rcp_const_cast;
1772 
1773  // MLP: NEED TO ADD CHECK IF PARAMETERS ARE SET LATER
1774 
1775  ScalarType one = Teuchos::ScalarTraits<ScalarType>::one();
1776  ScalarType zero = Teuchos::ScalarTraits<ScalarType>::zero();
1777  std::vector<int> index(numBlocks_+1);
1778 
1779  // MLP: EXCEPTION TESTING SHOULD GO HERE
1780 
1781  //THE FOLLOWING BLOCK OF CODE IS TO INFORM THE PROBLEM HOW MANY RIGHT HAND
1782  //SIDES WILL BE SOLVED. IF THE CHOSEN BLOCK SIZE IS THE SAME AS THE NUMBER
1783  //OF RIGHT HAND SIDES IN THE PROBLEM THEN WE PASS THE PROBLEM
1784  //INDICES FOR ALL RIGHT HAND SIDES. IF BLOCK SIZE IS GREATER THAN
1785  //THE NUMBER OF RIGHT HAND SIDES AND THE USER HAS ADAPTIVE BLOCK SIZING
1786  //TURNED OFF, THEN THE PROBLEM WILL GENERATE RANDOM RIGHT HAND SIDES SO WE HAVE AS MANY
1787  //RIGHT HAND SIDES AS BLOCK SIZE INDICATES. THIS MAY NEED TO BE CHANGED
1788  //LATER TO ACCOMADATE SMARTER CHOOSING OF FICTITIOUS RIGHT HAND SIDES.
1789  //THIS CODE IS DIRECTLY LIFTED FROM BelosBlockGmresSolMgr.hpp.
1790 
1791  // Create indices for the linear systems to be solved.
1792  int startPtr = 0;
1793  int numRHS2Solve = MVT::GetNumberVecs( *(problem_->getRHS()) );
1794  int numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
1795 
1796  //currIdx holds indices to all right-hand sides to be solved
1797  //or -1's to indicate that random right-hand sides should be generated
1798  std::vector<int> currIdx;
1799 
1800  if ( adaptiveBlockSize_ ) {
1801  blockSize_ = numCurrRHS;
1802  currIdx.resize( numCurrRHS );
1803  for (int i=0; i<numCurrRHS; ++i)
1804  currIdx[i] = startPtr+i;
1805  }
1806  else {
1807  currIdx.resize( blockSize_ );
1808  for (int i=0; i<numCurrRHS; ++i)
1809  currIdx[i] = startPtr+i;
1810  for (int i=numCurrRHS; i<blockSize_; ++i)
1811  currIdx[i] = -1;
1812  }
1813 
1814  // Inform the linear problem of the linear systems to solve/generate.
1815  problem_->setLSIndex( currIdx );
1816 
1817  //ADD ERROR CHECKING TO MAKE SURE SIZE OF BLOCK KRYLOV SUBSPACE NOT LARGER THAN dim
1818  //ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
1819 
1820  // reset loss of orthogonality flag
1821  loaDetected_ = false;
1822 
1823  // Assume convergence is achieved, then let any failed convergence set this to false.
1824  bool isConverged = true;
1825 
1826  // Initialize storage for all state variables
1827  initializeStateStorage();
1828 
1829  // Parameter list MLP: WE WILL NEED TO ASSIGN SOME PARAMETER VALUES LATER
1830  Teuchos::ParameterList plist;
1831 
1832  while (numRHS2Solve > 0){ //This loops through each block of RHS's being solved
1833  // **************************************************************************************************************************
1834  // Begin initial solver preparations. Either update U,C for new operator or generate an initial space using a cycle of GMRES
1835  // **************************************************************************************************************************
1836  //int prime_iterations; // silence warning about not being used
1837  if(keff > 0){//If there is already a subspace to recycle, then use it
1838  // Update U, C for the new operator
1839 
1840 // TEUCHOS_TEST_FOR_EXCEPTION(keff < recycledBlocks_,BlockGCRODRSolMgrRecyclingFailure,
1841 // "Belos::BlockGCRODRSolMgr::solve(): Requested size of recycled subspace is not consistent with the current recycle subspace.");
1842  printer_->stream(Debug) << " Now solving RHS index " << currIdx[0] << " using recycled subspace of dimension " << keff << std::endl << std::endl;
1843 
1844  // Compute image of U_ under the new operator
1845  index.resize(keff);
1846  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1847  RCP<const MV> Utmp = MVT::CloneView( *U_, index );
1848  RCP<MV> Ctmp = MVT::CloneViewNonConst( *C_, index );
1849  problem_->apply( *Utmp, *Ctmp );
1850 
1851  RCP<MV> U1tmp = MVT::CloneViewNonConst( *U1_, index );
1852 
1853  // Orthogonalize this block
1854  // Get a matrix to hold the orthonormalization coefficients.
1855  SDM Ftmp( Teuchos::View, *F_, keff, keff );
1856  int rank = ortho_->normalize(*Ctmp, rcp(&Ftmp,false));
1857  // Throw an error if we could not orthogonalize this block
1858  TEUCHOS_TEST_FOR_EXCEPTION(rank != keff,BlockGCRODRSolMgrOrthoFailure,"Belos::BlockGCRODRSolMgr::solve(): Failed to compute orthonormal basis for initial recycled subspace.");
1859 
1860  // U_ = U_*F^{-1}
1861  // First, compute LU factorization of R
1862  int info = 0;
1863  ipiv_.resize(Ftmp.numRows());
1864  lapack.GETRF(Ftmp.numRows(),Ftmp.numCols(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&info);
1865  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRF failed to compute an LU factorization.");
1866 
1867  // Now, form inv(F)
1868  int lwork = Ftmp.numRows();
1869  work_.resize(lwork);
1870  lapack.GETRI(Ftmp.numRows(),Ftmp.values(),Ftmp.stride(),&ipiv_[0],&work_[0],lwork,&info);
1871  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, BlockGCRODRSolMgrLAPACKFailure,"Belos::BlockGCRODRSolMgr::solve(): LAPACK _GETRI failed to invert triangular matrix.");
1872 
1873  // U_ = U1_; (via a swap)
1874  MVT::MvTimesMatAddMv( one, *Utmp, Ftmp, zero, *U1tmp );
1875  std::swap(U_, U1_);
1876 
1877  // Must reinitialize after swap
1878  index.resize(keff);
1879  for (int ii=0; ii<keff; ++ii) { index[ii] = ii; }
1880  Ctmp = MVT::CloneViewNonConst( *C_, index );
1881  Utmp = MVT::CloneView( *U_, index );
1882 
1883  // Compute C_'*R_
1884  SDM Ctr(keff,blockSize_);
1885  problem_->computeCurrPrecResVec( &*R_ );
1886  MVT::MvTransMv( one, *Ctmp, *R_, Ctr );
1887 
1888  // Update solution ( x += U_*C_'*R_ )
1889  RCP<MV> update = MVT::Clone( *problem_->getCurrLHSVec(), blockSize_ );
1890  MVT::MvInit( *update, 0.0 );
1891  MVT::MvTimesMatAddMv( one, *Utmp, Ctr, one, *update );
1892  problem_->updateSolution( update, true );
1893 
1894  // Update residual norm ( r -= C_*C_'*R_ )
1895  MVT::MvTimesMatAddMv( -one, *Ctmp, Ctr, one, *R_ );
1896 
1897  // We recycled space from previous call
1898  //prime_iterations = 0;
1899 
1900  // Since we have a previous recycle space, we do not need block_gmres_iter
1901  // and we must allocate V ourselves. See comments in initialize() routine for
1902  // further explanation.
1903  if (V_ == Teuchos::null) {
1904  // Check if there is a multivector to clone from.
1905  Teuchos::RCP<const MV> rhsMV = problem_->getRHS();
1906  V_ = MVT::Clone( *rhsMV, (numBlocks_+1)*blockSize_ );
1907  }
1908  else{
1909  // Generate V_ by cloning itself ONLY if more space is needed.
1910  if (MVT::GetNumberVecs(*V_) < (numBlocks_+1)*blockSize_ ) {
1911  Teuchos::RCP<const MV> tmp = V_;
1912  V_ = MVT::Clone( *tmp, (numBlocks_+1)*blockSize_ );
1913  }
1914  }
1915  } //end if(keff > 0)
1916  else{ //if there was no subspace to recycle, then build one
1917  printer_->stream(Debug) << " No recycled subspace available for RHS index " << std::endl << std::endl;
1918 
1919  Teuchos::ParameterList primeList;
1920  //we set this as numBlocks-1 because that is the Krylov dimension we want
1921  //so that the size of the Hessenberg created matches that of what we were expecting.
1922  primeList.set("Num Blocks",numBlocks_-1);
1923  primeList.set("Block Size",blockSize_);
1924  primeList.set("Recycled Blocks",0);
1925  primeList.set("Keep Hessenberg",true);
1926  primeList.set("Initialize Hessenberg",true);
1927 
1928  ptrdiff_t dim = MVT::GetGlobalLength( *(problem_->getRHS()) );
1929  if (blockSize_*static_cast<ptrdiff_t>(numBlocks_) > dim) {//if user has selected a total subspace dimension larger than system dimension
1930  ptrdiff_t tmpNumBlocks = 0;
1931  if (blockSize_ == 1)
1932  tmpNumBlocks = dim / blockSize_; // Allow for a good breakdown.
1933  else{
1934  tmpNumBlocks = ( dim - blockSize_) / blockSize_; // Allow for restarting.
1935  printer_->stream(Warnings) << "Belos::BlockGCRODRSolMgr::solve(): Warning! Requested Krylov subspace dimension is larger than operator dimension!"
1936  << std::endl << "The maximum number of blocks allowed for the Krylov subspace will be adjusted to " << tmpNumBlocks << std::endl;
1937  primeList.set("Num Blocks",Teuchos::as<int>(tmpNumBlocks));
1938  }
1939  }
1940  else{
1941  primeList.set("Num Blocks",numBlocks_-1);
1942  }
1943  //Create Block GMRES iteration object to perform one cycle of GMRES
1945  block_gmres_iter = Teuchos::rcp( new BlockGmresIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,primeList) );
1946 
1947  // MLP: ADD LOGIC TO DEAL WITH USER ASKING TO GENERATE A LARGER SPACE THAN dim AS IN HEIDI'S BlockGmresSolMgr CODE (DIDN'T WE ALREADY DO THIS SOMEWHERE?)
1948  block_gmres_iter->setSize( blockSize_, numBlocks_-1 );
1949 
1950  //Create the first block in the current BLOCK Krylov basis (using residual)
1951  Teuchos::RCP<MV> V_0;
1952  if (currIdx[blockSize_-1] == -1) {//if augmented with random RHS's
1953  V_0 = MVT::Clone( *(problem_->getInitPrecResVec()), blockSize_ );
1954  problem_->computeCurrPrecResVec( &*V_0 );
1955  }
1956  else {
1957  V_0 = MVT::CloneCopy( *(problem_->getInitPrecResVec()), currIdx );
1958  }
1959 
1960  // Get a matrix to hold the orthonormalization coefficients.
1961  Teuchos::RCP<SDM > z_0 = Teuchos::rcp( new SDM( blockSize_, blockSize_ ) );
1962 
1963  // Orthonormalize the new V_0
1964  int rank = ortho_->normalize( *V_0, z_0 );
1965  (void) rank; // silence compiler warning for unused variable
1966  // MLP: ADD EXCEPTION IF INITIAL BLOCK IS RANK DEFFICIENT
1967 
1968  // Set the new state and initialize the iteration.
1970  newstate.V = V_0;
1971  newstate.z = z_0;
1972  newstate.curDim = 0;
1973  block_gmres_iter->initializeGmres(newstate);
1974 
1975  bool primeConverged = false;
1976 
1977  try {
1978  printer_->stream(Debug) << " Preparing to Iterate!!!!" << std::endl << std::endl;
1979  block_gmres_iter->iterate();
1980 
1981  // *********************
1982  // check for convergence
1983  // *********************
1984  if ( convTest_->getStatus() == Passed ) {
1985  printer_->stream(Debug) << "We converged during the prime the pump stage" << std::endl << std::endl;
1986  primeConverged = !(expConvTest_->getLOADetected());
1987  if ( expConvTest_->getLOADetected() ) {
1988  // we don't have convergence
1989  loaDetected_ = true;
1990  printer_->stream(Warnings) << "Belos::BlockGCRODRSolMgr::solve(): Warning! Solver has experienced a loss of accuracy!" << std::endl;
1991  }
1992  }
1993  // *******************************************
1994  // check for maximum iterations of block GMRES
1995  // *******************************************
1996  else if( maxIterTest_->getStatus() == Passed ) {
1997  // we don't have convergence
1998  primeConverged = false;
1999  }
2000  // ****************************************************************
2001  // We need to recycle and continue, print a message indicating this
2002  // ****************************************************************
2003  else{ //debug statements
2004  printer_->stream(Debug) << " We did not converge on priming cycle of Block GMRES. Therefore we recycle and restart. " << std::endl << std::endl;
2005  }
2006  } //end try
2007  catch (const GmresIterationOrthoFailure &e) {
2008  // If the block size is not one, it's not considered a lucky breakdown.
2009  if (blockSize_ != 1) {
2010  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
2011  << block_gmres_iter->getNumIters() << std::endl
2012  << e.what() << std::endl;
2013  if (convTest_->getStatus() != Passed)
2014  primeConverged = false;
2015  }
2016  else {
2017  // If the block size is one, try to recover the most recent least-squares solution
2018  block_gmres_iter->updateLSQR( block_gmres_iter->getCurSubspaceDim() );
2019  // Check to see if the most recent least-squares solution yielded convergence.
2020  sTest_->checkStatus( &*block_gmres_iter );
2021  if (convTest_->getStatus() != Passed)
2022  isConverged = false;
2023  }
2024  } // end catch (const GmresIterationOrthoFailure &e)
2025  catch (const StatusTestNaNError& e) {
2026  // A NaN was detected in the solver. Set the solution to zero and return unconverged.
2027  achievedTol_ = MT::one();
2028  Teuchos::RCP<MV> X = problem_->getLHS();
2029  MVT::MvInit( *X, SCT::zero() );
2030  printer_->stream(Warnings) << "Belos::BlockGCRODRSolMgr::solve(): Warning! NaN has been detected!"
2031  << std::endl;
2032  return Unconverged;
2033  }
2034  catch (const std::exception &e) {
2035  printer_->stream(Errors) << "Error! Caught std::exception in BlockGmresIter::iterate() at iteration "
2036  << block_gmres_iter->getNumIters() << std::endl
2037  << e.what() << std::endl;
2038  throw;
2039  }
2040 
2041  // Record number of iterations in generating initial recycle spacec
2042  //prime_iterations = block_gmres_iter->getNumIters();//instantiated here because it is not needed outside of else{} scope; we'll see if this is true or not
2043 
2044  // Update the linear problem.
2045  RCP<MV> update = block_gmres_iter->getCurrentUpdate();
2046  problem_->updateSolution( update, true );
2047 
2048  // Update the block residual
2049 
2050  problem_->computeCurrPrecResVec( &*R_ );
2051 
2052  // Get the state.
2053  newstate = block_gmres_iter->getState();
2054  int p = newstate.curDim;
2055  (void) p; // silence compiler warning for unused variable
2056  H_->assign(*(newstate.H));//IS THIS A GOOD IDEA? I DID IT SO MEMBERS WOULD HAVE ACCESS TO H.
2057  // Get new Krylov vectors, store in V_
2058 
2059  // Since the block_gmres_iter returns the state
2060  // with const RCP's we need to cast newstate.V as
2061  // a non const RCP. This is okay because block_gmres_iter
2062  // is about to fall out of scope, so we are simply
2063  // rescuing *newstate.V from being destroyed so we can
2064  // use it for future block recycled GMRES cycles
2065  V_ = rcp_const_cast<MV>(newstate.V);
2066  newstate.V.release();
2067  //COMPUTE NEW RECYCLE SPACE SOMEHOW
2068  buildRecycleSpaceKryl(keff, block_gmres_iter);
2069  printer_->stream(Debug) << "Generated recycled subspace using RHS index " << currIdx[0] << " of dimension " << keff << std::endl << std::endl;
2070 
2071  // Return to outer loop if the priming solve
2072  // converged, set the next linear system.
2073  if (primeConverged) {
2074  /* MLP: POSSIBLY INCORRECT CODE WE ARE REPLACING *********************************
2075  // Inform the linear problem that we are
2076  // finished with this block linear system.
2077  problem_->setCurrLS();
2078 
2079  // Update indices for the linear systems to be solved.
2080  // KMS: Fix the numRHS2Solve; Copy from Heidi's BlockGmres code
2081  numRHS2Solve -= 1;
2082  printer_->stream(Debug) << numRHS2Solve << " Right hand sides left to solve" << std::endl;
2083  if ( numRHS2Solve > 0 ) {
2084  currIdx[0]++;
2085  // Set the next indices
2086  problem_->setLSIndex( currIdx );
2087  }
2088  else {
2089  currIdx.resize( numRHS2Solve );
2090  }
2091  *****************************************************************************/
2092 
2093  // Inform the linear problem that we are finished with this block linear system.
2094  problem_->setCurrLS();
2095 
2096  // Update indices for the linear systems to be solved.
2097  startPtr += numCurrRHS;
2098  numRHS2Solve -= numCurrRHS;
2099  if ( numRHS2Solve > 0 ) {
2100  numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
2101  if ( adaptiveBlockSize_ ) {
2102  blockSize_ = numCurrRHS;
2103  currIdx.resize( numCurrRHS );
2104  for (int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2105  }
2106  else {
2107  currIdx.resize( blockSize_ );
2108  for (int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2109  for (int i=numCurrRHS; i<blockSize_; ++i) currIdx[i] = -1;
2110  }
2111  // Set the next indices.
2112  problem_->setLSIndex( currIdx );
2113  }
2114  else {
2115  currIdx.resize( numRHS2Solve );
2116  }
2117  continue;//Begin solving the next block of RHS's
2118  } //end if (primeConverged)
2119  } //end else [if(keff > 0)]
2120 
2121  // *****************************************
2122  // Initial subspace update/construction done
2123  // Start cycles of recycled GMRES
2124  // *****************************************
2125  Teuchos::ParameterList blockgcrodrList;
2126  blockgcrodrList.set("Num Blocks",numBlocks_);
2127  blockgcrodrList.set("Block Size",blockSize_);
2128  blockgcrodrList.set("Recycled Blocks",keff);
2129 
2131  block_gcrodr_iter = Teuchos::rcp( new BlockGCRODRIter<ScalarType,MV,OP>(problem_,printer_,outputTest_,ortho_,blockgcrodrList) );
2133 
2134  index.resize( blockSize_ );
2135  for(int ii = 0; ii < blockSize_; ii++) index[ii] = ii;
2136  Teuchos::RCP<MV> V0 = MVT::CloneViewNonConst( *V_, index );
2137 
2138  // MLP: MVT::SetBlock(*R_,index,*V0);
2139  MVT::Assign(*R_,*V0);
2140 
2141  index.resize(keff);//resize to get appropriate recycle space vectors
2142  for(int i=0; i < keff; i++){ index[i] = i;};
2143  B_ = rcp(new SDM(Teuchos::View, *G_, keff, numBlocks_*blockSize_, 0, keff));
2144  H_ = rcp(new SDM(Teuchos::View, *G_, (numBlocks_-1)*blockSize_ + blockSize_, (numBlocks_-1)*blockSize_, keff ,keff ));
2145 
2146  newstate.V = V_;
2147  newstate.B= B_;
2148  newstate.U = MVT::CloneViewNonConst(*U_, index);
2149  newstate.C = MVT::CloneViewNonConst(*C_, index);
2150  newstate.H = H_;
2151  newstate.curDim = blockSize_;
2152  block_gcrodr_iter -> initialize(newstate);
2153 
2154  int numRestarts = 0;
2155 
2156  while(1){//Each execution of this loop is a cycle of block GCRODR
2157  try{
2158  block_gcrodr_iter -> iterate();
2159 
2160  // ***********************
2161  // Check Convergence First
2162  // ***********************
2163  if( convTest_->getStatus() == Passed ) {
2164  // we have convergence
2165  break;//from while(1)
2166  } //end if converged
2167 
2168  // ***********************************
2169  // Check if maximum iterations reached
2170  // ***********************************
2171  else if(maxIterTest_->getStatus() == Passed ){
2172  // no convergence; hit maxit
2173  isConverged = false;
2174  break; // from while(1)
2175  } //end elseif reached maxit
2176 
2177  // **********************************************
2178  // Check if subspace full; do we need to restart?
2179  // **********************************************
2180  else if (block_gcrodr_iter->getCurSubspaceDim() == block_gcrodr_iter->getMaxSubspaceDim()){
2181 
2182  // Update recycled space even if we have reached max number of restarts
2183 
2184  // Update linear problem
2185  Teuchos::RCP<MV> update = block_gcrodr_iter->getCurrentUpdate();
2186  problem_->updateSolution(update, true);
2187  buildRecycleSpaceAugKryl(block_gcrodr_iter);
2188 
2189  printer_->stream(Debug) << " Generated new recycled subspace using RHS index " << currIdx[0] << " of dimension " << keff << std::endl << std::endl;
2190  // NOTE: If we have hit the maximum number of restarts, then we will quit
2191  if(numRestarts >= maxRestarts_) {
2192  isConverged = false;
2193  break; //from while(1)
2194  } //end if max restarts
2195 
2196  numRestarts++;
2197 
2198  printer_ -> stream(Debug) << " Performing restart number " << numRestarts << " of " << maxRestarts_ << std::endl << std::endl;
2199 
2200  // Create the restart vector (first block in the current Krylov basis)
2201  problem_->computeCurrPrecResVec( &*R_ );
2202  index.resize( blockSize_ );
2203  for (int ii=0; ii<blockSize_; ++ii) index[ii] = ii;
2204  Teuchos::RCP<MV> V0 = MVT::CloneViewNonConst( *V_, index );
2205  MVT::SetBlock(*R_,index,*V0);
2206 
2207  // Set the new state and initialize the solver.
2209  index.resize( numBlocks_*blockSize_ );
2210  for (int ii=0; ii<(numBlocks_*blockSize_); ++ii) index[ii] = ii;
2211  restartState.V = MVT::CloneViewNonConst( *V_, index );
2212  index.resize( keff );
2213  for (int ii=0; ii<keff; ++ii) index[ii] = ii;
2214  restartState.U = MVT::CloneViewNonConst( *U_, index );
2215  restartState.C = MVT::CloneViewNonConst( *C_, index );
2216  B_ = rcp(new SDM(Teuchos::View, *G_, keff, (numBlocks_-1)*blockSize_, 0, keff));
2217  H_ = rcp(new SDM(Teuchos::View, *G_, numBlocks_*blockSize_, (numBlocks_-1)*blockSize_, keff ,keff ));
2218  restartState.B = B_;
2219  restartState.H = H_;
2220  restartState.curDim = blockSize_;
2221  block_gcrodr_iter->initialize(restartState);
2222 
2223  } //end else if need to restart
2224 
2225  // ****************************************************************
2226  // We returned from iterate(), but none of our status tests passed.
2227  // Something is wrong, and it is probably our fault.
2228  // ****************************************************************
2229  else {
2230  TEUCHOS_TEST_FOR_EXCEPTION(true,std::logic_error,"Belos::BlockGCRODRSolMgr::solve(): Invalid return from BlockGCRODRIter::iterate().");
2231  } //end else (no status test passed)
2232 
2233  }//end try
2234  catch(const BlockGCRODRIterOrthoFailure &e){ //there was an exception
2235  // Try to recover the most recent least-squares solution
2236  block_gcrodr_iter->updateLSQR( block_gcrodr_iter->getCurSubspaceDim() );
2237  // Check to see if the most recent least-squares solution yielded convergence.
2238  sTest_->checkStatus( &*block_gcrodr_iter );
2239  if (convTest_->getStatus() != Passed) isConverged = false;
2240  break;
2241  } // end catch orthogonalization failure
2242  catch(const std::exception &e){
2243  printer_->stream(Errors) << "Error! Caught exception in BlockGCRODRIter::iterate() at iteration "
2244  << block_gcrodr_iter->getNumIters() << std::endl
2245  << e.what() << std::endl;
2246  throw;
2247  } //end catch standard exception
2248  } //end while(1)
2249 
2250  // Compute the current solution.
2251  // Update the linear problem.
2252  Teuchos::RCP<MV> update = block_gcrodr_iter->getCurrentUpdate();
2253  problem_->updateSolution( update, true );
2254 
2255  /* MLP: POSSIBLY INCORRECT CODE WE ARE REPLACING *********************************
2256  // Inform the linear problem that we are finished with this block linear system.
2257  problem_->setCurrLS();
2258 
2259  // Update indices for the linear systems to be solved.
2260  numRHS2Solve -= 1;
2261  if ( numRHS2Solve > 0 ) {
2262  currIdx[0]++;
2263  // Set the next indices.
2264  problem_->setLSIndex( currIdx );
2265  }
2266  else {
2267  currIdx.resize( numRHS2Solve );
2268  }
2269  ******************************************************************************/
2270 
2271  // Inform the linear problem that we are finished with this block linear system.
2272  problem_->setCurrLS();
2273 
2274  // Update indices for the linear systems to be solved.
2275  startPtr += numCurrRHS;
2276  numRHS2Solve -= numCurrRHS;
2277  if ( numRHS2Solve > 0 ) {
2278  numCurrRHS = ( numRHS2Solve < blockSize_) ? numRHS2Solve : blockSize_;
2279  if ( adaptiveBlockSize_ ) {
2280  blockSize_ = numCurrRHS;
2281  currIdx.resize( numCurrRHS );
2282  for (int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2283  }
2284  else {
2285  currIdx.resize( blockSize_ );
2286  for (int i=0; i<numCurrRHS; ++i) currIdx[i] = startPtr+i;
2287  for (int i=numCurrRHS; i<blockSize_; ++i) currIdx[i] = -1;
2288  }
2289  // Set the next indices.
2290  problem_->setLSIndex( currIdx );
2291  }
2292  else {
2293  currIdx.resize( numRHS2Solve );
2294  }
2295 
2296  // If we didn't build a recycle space this solve but ran at least k iterations, force build of new recycle space
2297  if (!builtRecycleSpace_) {
2298  buildRecycleSpaceAugKryl(block_gcrodr_iter);
2299  printer_->stream(Debug) << " Generated new recycled subspace using RHS index " << currIdx[0] << " of dimension " << keff << std::endl << std::endl;
2300  }
2301  }//end while (numRHS2Solve > 0)
2302 
2303  // print final summary
2304  sTest_->print( printer_->stream(FinalSummary) );
2305 
2306  // print timing information
2307  #ifdef BELOS_TEUCHOS_TIME_MONITOR
2308  // Calling summarize() can be expensive, so don't call unless the user wants to print out timing details.
2309  // summarize() will do all the work even if it's passed a "black hole" output stream.
2310  if (verbosity_ & TimingDetails) Teuchos::TimeMonitor::summarize( printer_->stream(TimingDetails) );
2311  #endif
2312  // get iteration information for this solve
2313  numIters_ = maxIterTest_->getNumIters();
2314 
2315  // get residual information for this solve
2316  const std::vector<MagnitudeType>* pTestValues = NULL;
2317  pTestValues = impConvTest_->getTestValue();
2318  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues == NULL, std::logic_error,
2319  "Belos::BlockGCRODRSolMgr::solve(): The implicit convergence test's "
2320  "getTestValue() method returned NULL. Please report this bug to the "
2321  "Belos developers.");
2322  TEUCHOS_TEST_FOR_EXCEPTION(pTestValues->size() < 1, std::logic_error,
2323  "Belos::BlockGCRODRSolMgr::solve(): The implicit convergence test's "
2324  "getTestValue() method returned a vector of length zero. Please report "
2325  "this bug to the Belos developers.");
2326  // FIXME (mfh 12 Dec 2011) Does pTestValues really contain the
2327  // achieved tolerances for all vectors in the current solve(), or
2328  // just for the vectors from the last deflation?
2329  achievedTol_ = *std::max_element (pTestValues->begin(), pTestValues->end());
2330 
2331  if (!isConverged) return Unconverged; // return from BlockGCRODRSolMgr::solve()
2332  return Converged; // return from BlockGCRODRSolMgr::solve()
2333 } //end solve()
2334 
2335 } //End Belos Namespace
2336 
2337 #endif /* BELOS_BLOCK_GCRODR_SOLMGR_HPP */
ScaleType convertStringToScaleType(const std::string &scaleType)
Convert the given string to its ScaleType enum value.
Definition: BelosTypes.cpp:74
ScalarType * values() const
Collection of types and exceptions used within the Belos solvers.
const LinearProblem< ScalarType, MV, OP > & getProblem() const
Get current linear problem being solved for in this object.
Belos&#39;s basic output manager for sending information of select verbosity levels to the appropriate ou...
Class which manages the output and verbosity of the Belos solvers.
Structure to contain pointers to BlockGCRODRIter state variables.
virtual ~BlockGCRODRSolMgr()
Destructor.
int getNumIters() const
Get the iteration count for the most recent call to solve().
static void clearCounter(const std::string &name)
ScaleType
The type of scaling to use on the residual norm value.
Definition: BelosTypes.hpp:88
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > B
The projection of the Krylov subspace against the recycled subspace *.
Exception thrown to signal error in a status test during Belos::StatusTest::checkStatus().
Teuchos::RCP< const Teuchos::ParameterList > getCurrentParameters() const
Get a parameter list containing the current parameters for this object.
This class implements the block GMRES iteration, where a block Krylov subspace is constructed...
Teuchos::RCP< const MV > V
The current Krylov basis.
T & get(const std::string &name, T def_value)
BlockGCRODRIterOrthoFailure is thrown when the BlockGCRODRIter object is unable to compute independen...
static RCP< Time > getNewCounter(const std::string &name)
Thrown when the linear problem was not set up correctly.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
A factory class for generating StatusTestOutput objects.
static magnitudeType real(T a)
void setParameters(const Teuchos::RCP< Teuchos::ParameterList > &params)
Set the parameters the solver should use to solve the linear problem.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > z
The current right-hand side of the least squares system RY = Z.
An implementation of StatusTestResNorm using a family of residual norms.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
Structure to contain pointers to GmresIteration state variables.
Belos::StatusTest class for specifying a maximum number of iterations.
BlockGCRODRSolMgrRecyclingFailure(const std::string &what_arg)
void reset(const ResetType type)
Performs a reset of the solver manager specified by the ResetType.
Thrown if any problem occurs in using or creating the recycle subspace.
A factory class for generating StatusTestOutput objects.
Teuchos::RCP< MV > V
The current Krylov basis.
BlockGCRODRSolMgrOrthoFailure(const std::string &what_arg)
int curDim
The current dimension of the reduction.
Traits class which defines basic operations on multivectors.
Belos::StatusTest for logically combining several status tests.
A Belos::StatusTest class for specifying a maximum number of iterations.
ResetType
How to reset the solver.
Definition: BelosTypes.hpp:174
BlockGCRODRSolMgrLAPACKFailure(const std::string &what_arg)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
A solver manager for the Block GCRO-DR (Block Recycling GMRES) linear solver.
Pure virtual base class which describes the basic interface for a solver manager. ...
int curDim
The current dimension of the reduction.
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
MagnitudeType achievedTol() const
Get the residual for the most recent call to solve().
Belos concrete class for performing the block GCRO-DR (block GMRES with recycling) iteration...
BlockGCRODRSolMgr()
Default constructor.
A linear system to solve, and its associated information.
std::string description() const
A description of the Block GCRODR solver manager.
Class which describes the linear problem to be solved by the iterative solver.
void setProblem(const Teuchos::RCP< LinearProblem< ScalarType, MV, OP > > &problem)
Set the linear problem to solve on the next call to solve().
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
ReturnType
Whether the Belos solve converged for all linear systems.
Definition: BelosTypes.hpp:123
The Belos::SolverManager is a templated virtual base class that defines the basic interface that any ...
Thrown when the solution manager was unable to orthogonalize the basis vectors.
Teuchos::RCP< MV > U
The recycled subspace and its projection.
Belos concrete class for performing the block GMRES iteration.
Ptr< T > release()
Implementation of the Block GCRO-DR (Block Recycling GMRES) iteration.
static magnitudeType magnitude(T a)
OrdinalType numCols() const
int setTolerance(MagnitudeType tolerance)
Set the value of the tolerance.
Teuchos::RCP< StatusTestOutput< ScalarType, MV, OP > > create(const Teuchos::RCP< OutputManager< ScalarType > > &printer, Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test, int mod, int printStates)
Create the StatusTestOutput object specified by the outputStyle.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const
Get a parameter list containing the valid parameters for this object.
Belos::StatusTestResNorm for specifying general residual norm stopping criteria.
Teuchos::RCP< Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
An implementation of the Belos::MatOrthoManager that performs orthogonalization using (potentially) m...
A class for extending the status testing capabilities of Belos via logical combinations.
Class which defines basic traits for the operator type.
Parent class to all Belos exceptions.
Definition: BelosTypes.hpp:28
Pure virtual base class which augments the basic interface for a Gmres linear solver iteration...
GmresIterationOrthoFailure is thrown when the GmresIteration object is unable to compute independent ...
Belos header file which uses auto-configuration information to include necessary C++ headers...
BlockGCRODRSolMgrLinearProblemFailure(const std::string &what_arg)
ReturnType solve()
Solve the current linear problem.
Thrown when an LAPACK call fails.
OrdinalType stride() const
OrdinalType numRows() const
Belos concrete class for performing the block, flexible GMRES iteration.
bool isLOADetected() const
Whether a loss of accuracy was detected during the most recent solve.
OutputType
Style of output used to display status test information.
Definition: BelosTypes.hpp:109
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
bool is_null() const

Generated on Fri Dec 20 2024 09:24:48 for Belos by doxygen 1.8.5