Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnasaziBlockKrylovSchur.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Anasazi: Block Eigensolvers Package
4 //
5 // Copyright 2004 NTESS and the Anasazi contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
14 #ifndef ANASAZI_BLOCK_KRYLOV_SCHUR_HPP
15 #define ANASAZI_BLOCK_KRYLOV_SCHUR_HPP
16 
17 #include "AnasaziTypes.hpp"
18 
19 #include "AnasaziEigensolver.hpp"
22 #include "Teuchos_ScalarTraits.hpp"
23 #include "AnasaziHelperTraits.hpp"
24 
25 #include "AnasaziOrthoManager.hpp"
26 
27 #include "Teuchos_LAPACK.hpp"
28 #include "Teuchos_BLAS.hpp"
31 #include "Teuchos_TimeMonitor.hpp"
32 
47 namespace Anasazi {
48 
50 
51 
56  template <class ScalarType, class MulVec>
62  int curDim;
75  BlockKrylovSchurState() : curDim(0), V(Teuchos::null),
76  H(Teuchos::null), S(Teuchos::null),
77  Q(Teuchos::null) {}
78  };
79 
81 
83 
84 
97  class BlockKrylovSchurInitFailure : public AnasaziError {public:
98  BlockKrylovSchurInitFailure(const std::string& what_arg) : AnasaziError(what_arg)
99  {}};
100 
108  BlockKrylovSchurOrthoFailure(const std::string& what_arg) : AnasaziError(what_arg)
109  {}};
110 
112 
113 
114  template <class ScalarType, class MV, class OP>
115  class BlockKrylovSchur : public Eigensolver<ScalarType,MV,OP> {
116  public:
118 
119 
133  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
136  Teuchos::ParameterList &params
137  );
138 
140  virtual ~BlockKrylovSchur() {};
142 
143 
145 
146 
168  void iterate();
169 
192 
196  void initialize();
197 
206  bool isInitialized() const { return initialized_; }
207 
217  state.curDim = curDim_;
218  state.V = V_;
219  state.H = H_;
220  state.Q = Q_;
221  state.S = schurH_;
222  return state;
223  }
224 
226 
227 
229 
230 
232  int getNumIters() const { return(iter_); }
233 
235  void resetNumIters() { iter_=0; }
236 
244  Teuchos::RCP<const MV> getRitzVectors() { return ritzVectors_; }
245 
253  std::vector<Value<ScalarType> > getRitzValues() {
254  std::vector<Value<ScalarType> > ret = ritzValues_;
255  ret.resize(ritzIndex_.size());
256  return ret;
257  }
258 
264  std::vector<int> getRitzIndex() { return ritzIndex_; }
265 
270  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getResNorms() {
271  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> ret(0);
272  return ret;
273  }
274 
279  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRes2Norms() {
280  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> ret(0);
281  return ret;
282  }
283 
288  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRitzRes2Norms() {
289  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> ret = ritzResiduals_;
290  ret.resize(ritzIndex_.size());
291  return ret;
292  }
293 
295 
297 
298 
301 
304 
306  const Eigenproblem<ScalarType,MV,OP>& getProblem() const { return(*problem_); };
307 
314  void setSize(int blockSize, int numBlocks);
315 
317  void setBlockSize(int blockSize);
318 
320  void setStepSize(int stepSize);
321 
323  void setNumRitzVectors(int numRitzVecs);
324 
326  int getStepSize() const { return(stepSize_); }
327 
329  int getBlockSize() const { return(blockSize_); }
330 
332  int getNumRitzVectors() const { return(numRitzVecs_); }
333 
339  int getCurSubspaceDim() const {
340  if (!initialized_) return 0;
341  return curDim_;
342  }
343 
345  int getMaxSubspaceDim() const { return (problem_->isHermitian()?blockSize_*numBlocks_:blockSize_*numBlocks_+1); }
346 
347 
360  void setAuxVecs(const Teuchos::Array<Teuchos::RCP<const MV> > &auxvecs);
361 
364 
366 
368 
369 
371  void currentStatus(std::ostream &os);
372 
374 
376 
377 
379  bool isRitzVecsCurrent() const { return ritzVecsCurrent_; }
380 
382  bool isRitzValsCurrent() const { return ritzValsCurrent_; }
383 
385  bool isSchurCurrent() const { return schurCurrent_; }
386 
388 
390 
391 
393  void computeRitzVectors();
394 
396  void computeRitzValues();
397 
399  void computeSchurForm( const bool sort = true );
400 
402 
403  private:
404  //
405  // Convenience typedefs
406  //
407  typedef MultiVecTraits<ScalarType,MV> MVT;
410  typedef typename SCT::magnitudeType MagnitudeType;
411  typedef typename std::vector<ScalarType>::iterator STiter;
412  typedef typename std::vector<MagnitudeType>::iterator MTiter;
413  const MagnitudeType MT_ONE;
414  const MagnitudeType MT_ZERO;
415  const MagnitudeType NANVAL;
416  const ScalarType ST_ONE;
417  const ScalarType ST_ZERO;
418  //
419  // Internal structs
420  //
421  struct CheckList {
422  bool checkV;
423  bool checkArn;
424  bool checkAux;
425  CheckList() : checkV(false), checkArn(false), checkAux(false) {};
426  };
427  //
428  // Internal methods
429  //
430  std::string accuracyCheck(const CheckList &chk, const std::string &where) const;
431  void sortSchurForm( Teuchos::SerialDenseMatrix<int,ScalarType>& H,
433  std::vector<int>& order );
434  //
435  // Classes inputed through constructor that define the eigenproblem to be solved.
436  //
442  //
443  // Information obtained from the eigenproblem
444  //
446  //
447  // Internal timers
448  //
449  Teuchos::RCP<Teuchos::Time> timerOp_, timerSortRitzVal_,
450  timerCompSF_, timerSortSF_,
451  timerCompRitzVec_, timerOrtho_;
452  //
453  // Counters
454  //
455  int count_ApplyOp_;
456 
457  //
458  // Algorithmic parameters.
459  //
460  // blockSize_ is the solver block size; it controls the number of eigenvectors that
461  // we compute, the number of residual vectors that we compute, and therefore the number
462  // of vectors added to the basis on each iteration.
463  int blockSize_;
464  // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
465  int numBlocks_;
466  // stepSize_ dictates how many iterations are performed before eigenvectors and eigenvalues
467  // are computed again
468  int stepSize_;
469 
470  //
471  // Current solver state
472  //
473  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
474  // is capable of running; _initialize is controlled by the initialize() member method
475  // For the implications of the state of initialized_, please see documentation for initialize()
476  bool initialized_;
477  //
478  // curDim_ reflects how much of the current basis is valid
479  // NOTE: for Hermitian, 0 <= curDim_ <= blockSize_*numBlocks_
480  // for non-Hermitian, 0 <= curDim_ <= blockSize_*numBlocks_ + 1
481  // this also tells us how many of the values in _theta are valid Ritz values
482  int curDim_;
483  //
484  // State Multivecs
485  Teuchos::RCP<MV> ritzVectors_, V_;
486  int numRitzVecs_;
487  //
488  // Projected matrices
489  // H_ : Projected matrix from the Krylov-Schur factorization AV = VH + FB^T
490  //
492  //
493  // Schur form of Projected matrices (these are only updated when the Ritz values/vectors are updated).
494  // schurH_: Schur form reduction of H
495  // Q_: Schur vectors of H
498  //
499  // Auxiliary vectors
501  int numAuxVecs_;
502  //
503  // Number of iterations that have been performed.
504  int iter_;
505  //
506  // State flags
507  bool ritzVecsCurrent_, ritzValsCurrent_, schurCurrent_;
508  //
509  // Current eigenvalues, residual norms
510  std::vector<Value<ScalarType> > ritzValues_;
511  std::vector<MagnitudeType> ritzResiduals_;
512  //
513  // Current index vector for Ritz values and vectors
514  std::vector<int> ritzIndex_; // computed by BKS
515  std::vector<int> ritzOrder_; // returned from sort manager
516  //
517  // Number of Ritz pairs to be printed upon output, if possible
518  int numRitzPrint_;
519  };
520 
521 
523  // Constructor
524  template <class ScalarType, class MV, class OP>
528  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
531  Teuchos::ParameterList &params
532  ) :
533  MT_ONE(Teuchos::ScalarTraits<MagnitudeType>::one()),
534  MT_ZERO(Teuchos::ScalarTraits<MagnitudeType>::zero()),
535  NANVAL(Teuchos::ScalarTraits<MagnitudeType>::nan()),
536  ST_ONE(Teuchos::ScalarTraits<ScalarType>::one()),
537  ST_ZERO(Teuchos::ScalarTraits<ScalarType>::zero()),
538  // problem, tools
539  problem_(problem),
540  sm_(sorter),
541  om_(printer),
542  tester_(tester),
543  orthman_(ortho),
544  // timers, counters
545 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
546  timerOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchur::Operation Op*x")),
547  timerSortRitzVal_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchur::Sorting Ritz values")),
548  timerCompSF_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchur::Computing Schur form")),
549  timerSortSF_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchur::Sorting Schur form")),
550  timerCompRitzVec_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchur::Computing Ritz vectors")),
551  timerOrtho_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockKrylovSchur::Orthogonalization")),
552 #endif
553  count_ApplyOp_(0),
554  // internal data
555  blockSize_(0),
556  numBlocks_(0),
557  stepSize_(0),
558  initialized_(false),
559  curDim_(0),
560  numRitzVecs_(0),
561  auxVecs_( Teuchos::Array<Teuchos::RCP<const MV> >(0) ),
562  numAuxVecs_(0),
563  iter_(0),
564  ritzVecsCurrent_(false),
565  ritzValsCurrent_(false),
566  schurCurrent_(false),
567  numRitzPrint_(0)
568  {
569  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,std::invalid_argument,
570  "Anasazi::BlockKrylovSchur::constructor: user specified null problem pointer.");
571  TEUCHOS_TEST_FOR_EXCEPTION(sm_ == Teuchos::null,std::invalid_argument,
572  "Anasazi::BlockKrylovSchur::constructor: user passed null sort manager pointer.");
573  TEUCHOS_TEST_FOR_EXCEPTION(om_ == Teuchos::null,std::invalid_argument,
574  "Anasazi::BlockKrylovSchur::constructor: user passed null output manager pointer.");
575  TEUCHOS_TEST_FOR_EXCEPTION(tester_ == Teuchos::null,std::invalid_argument,
576  "Anasazi::BlockKrylovSchur::constructor: user passed null status test pointer.");
577  TEUCHOS_TEST_FOR_EXCEPTION(orthman_ == Teuchos::null,std::invalid_argument,
578  "Anasazi::BlockKrylovSchur::constructor: user passed null orthogonalization manager pointer.");
579  TEUCHOS_TEST_FOR_EXCEPTION(problem_->isProblemSet() == false, std::invalid_argument,
580  "Anasazi::BlockKrylovSchur::constructor: user specified problem is not set.");
581  TEUCHOS_TEST_FOR_EXCEPTION(sorter == Teuchos::null,std::invalid_argument,
582  "Anasazi::BlockKrylovSchur::constructor: user specified null sort manager pointer.");
583  TEUCHOS_TEST_FOR_EXCEPTION(printer == Teuchos::null,std::invalid_argument,
584  "Anasazi::BlockKrylovSchur::constructor: user specified null output manager pointer.");
585  TEUCHOS_TEST_FOR_EXCEPTION(tester == Teuchos::null,std::invalid_argument,
586  "Anasazi::BlockKrylovSchur::constructor: user specified null status test pointer.");
587  TEUCHOS_TEST_FOR_EXCEPTION(ortho == Teuchos::null,std::invalid_argument,
588  "Anasazi::BlockKrylovSchur::constructor: user specified null ortho manager pointer.");
589 
590  // Get problem operator
591  Op_ = problem_->getOperator();
592 
593  // get the step size
594  TEUCHOS_TEST_FOR_EXCEPTION(!params.isParameter("Step Size"), std::invalid_argument,
595  "Anasazi::BlockKrylovSchur::constructor: mandatory parameter 'Step Size' is not specified.");
596  int ss = params.get("Step Size",numBlocks_);
597  setStepSize(ss);
598 
599  // set the block size and allocate data
600  int bs = params.get("Block Size", 1);
601  int nb = params.get("Num Blocks", 3*problem_->getNEV());
602  setSize(bs,nb);
603 
604  // get the number of Ritz vectors to compute and allocate data.
605  // --> if this parameter is not specified in the parameter list, then it's assumed that no Ritz vectors will be computed.
606  int numRitzVecs = params.get("Number of Ritz Vectors", 0);
607  setNumRitzVectors( numRitzVecs );
608 
609  // get the number of Ritz values to print out when currentStatus is called.
610  numRitzPrint_ = params.get("Print Number of Ritz Values", bs);
611  }
612 
613 
615  // Set the block size
616  // This simply calls setSize(), modifying the block size while retaining the number of blocks.
617  template <class ScalarType, class MV, class OP>
619  {
620  setSize(blockSize,numBlocks_);
621  }
622 
623 
625  // Set the step size.
626  template <class ScalarType, class MV, class OP>
628  {
629  TEUCHOS_TEST_FOR_EXCEPTION(stepSize <= 0, std::invalid_argument, "Anasazi::BlockKrylovSchur::setStepSize(): new step size must be positive and non-zero.");
630  stepSize_ = stepSize;
631  }
632 
634  // Set the number of Ritz vectors to compute.
635  template <class ScalarType, class MV, class OP>
637  {
638  // This routine only allocates space; it doesn't not perform any computation
639  // any change in size will invalidate the state of the solver.
640 
641  TEUCHOS_TEST_FOR_EXCEPTION(numRitzVecs < 0, std::invalid_argument, "Anasazi::BlockKrylovSchur::setNumRitzVectors(): number of Ritz vectors to compute must be positive.");
642 
643  // Check to see if the number of requested Ritz vectors has changed.
644  if (numRitzVecs != numRitzVecs_) {
645  if (numRitzVecs) {
646  ritzVectors_ = Teuchos::null;
647  ritzVectors_ = MVT::Clone(*V_, numRitzVecs);
648  } else {
649  ritzVectors_ = Teuchos::null;
650  }
651  numRitzVecs_ = numRitzVecs;
652  ritzVecsCurrent_ = false;
653  }
654  }
655 
657  // Set the block size and make necessary adjustments.
658  template <class ScalarType, class MV, class OP>
659  void BlockKrylovSchur<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks)
660  {
661  // This routine only allocates space; it doesn't not perform any computation
662  // any change in size will invalidate the state of the solver.
663 
664  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks <= 0 || blockSize <= 0, std::invalid_argument, "Anasazi::BlockKrylovSchur::setSize was passed a non-positive argument.");
665  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks < 3, std::invalid_argument, "Anasazi::BlockKrylovSchur::setSize(): numBlocks must be at least three.");
666  if (blockSize == blockSize_ && numBlocks == numBlocks_) {
667  // do nothing
668  return;
669  }
670 
671  blockSize_ = blockSize;
672  numBlocks_ = numBlocks;
673 
675  // grab some Multivector to Clone
676  // in practice, getInitVec() should always provide this, but it is possible to use a
677  // Eigenproblem with nothing in getInitVec() by manually initializing with initialize();
678  // in case of that strange scenario, we will try to Clone from V_; first resort to getInitVec(),
679  // because we would like to clear the storage associated with V_ so we have room for the new V_
680  if (problem_->getInitVec() != Teuchos::null) {
681  tmp = problem_->getInitVec();
682  }
683  else {
684  tmp = V_;
685  TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
686  "Anasazi::BlockKrylovSchur::setSize(): eigenproblem did not specify initial vectors to clone from.");
687  }
688 
689 
691  // blockSize*numBlocks dependent
692  //
693  int newsd;
694  if (problem_->isHermitian()) {
695  newsd = blockSize_*numBlocks_;
696  } else {
697  newsd = blockSize_*numBlocks_+1;
698  }
699  // check that new size is valid
700  TEUCHOS_TEST_FOR_EXCEPTION(static_cast<ptrdiff_t>(newsd) > MVT::GetGlobalLength(*tmp),std::invalid_argument,
701  "Anasazi::BlockKrylovSchur::setSize(): maximum basis size is larger than problem dimension.");
702 
703  ritzValues_.resize(newsd);
704  ritzResiduals_.resize(newsd,MT_ONE);
705  ritzOrder_.resize(newsd);
706  V_ = Teuchos::null;
707  V_ = MVT::Clone(*tmp,newsd+blockSize_);
708  H_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(newsd+blockSize_,newsd) );
710 
711  initialized_ = false;
712  curDim_ = 0;
713  }
714 
715 
717  // Set the auxiliary vectors
718  template <class ScalarType, class MV, class OP>
720  typedef typename Teuchos::Array<Teuchos::RCP<const MV> >::iterator tarcpmv;
721 
722  // set new auxiliary vectors
723  auxVecs_ = auxvecs;
724 
725  if (om_->isVerbosity( Debug ) ) {
726  // Check almost everything here
727  CheckList chk;
728  chk.checkAux = true;
729  om_->print( Debug, accuracyCheck(chk, ": in setAuxVecs()") );
730  }
731 
732  numAuxVecs_ = 0;
733  for (tarcpmv i=auxVecs_.begin(); i != auxVecs_.end(); i++) {
734  numAuxVecs_ += MVT::GetNumberVecs(**i);
735  }
736 
737  // If the solver has been initialized, X and P are not necessarily orthogonal to new auxiliary vectors
738  if (numAuxVecs_ > 0 && initialized_) {
739  initialized_ = false;
740  }
741  }
742 
744  /* Initialize the state of the solver
745  *
746  * POST-CONDITIONS:
747  *
748  * V_ is orthonormal, orthogonal to auxVecs_, for first curDim_ vectors
749  *
750  */
751 
752  template <class ScalarType, class MV, class OP>
754  {
755  // NOTE: memory has been allocated by setBlockSize(). Use SetBlock below; do not Clone
756 
757  std::vector<int> bsind(blockSize_);
758  for (int i=0; i<blockSize_; i++) bsind[i] = i;
759 
760  // in BlockKrylovSchur, V and H are required.
761  // if either doesn't exist, then we will start with the initial vector.
762  //
763  // inconsistent multivectors widths and lengths will not be tolerated, and
764  // will be treated with exceptions.
765  //
766  std::string errstr("Anasazi::BlockKrylovSchur::initialize(): specified multivectors must have a consistent length and width.");
767 
768  // set up V,H: if the user doesn't specify both of these these,
769  // we will start over with the initial vector.
770  if (newstate.V != Teuchos::null && newstate.H != Teuchos::null) {
771 
772  // initialize V_,H_, and curDim_
773 
774  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_),
775  std::invalid_argument, errstr );
776  if (newstate.V != V_) {
777  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < blockSize_,
778  std::invalid_argument, errstr );
779  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) > getMaxSubspaceDim(),
780  std::invalid_argument, errstr );
781  }
782  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim > getMaxSubspaceDim(),
783  std::invalid_argument, errstr );
784 
785  curDim_ = newstate.curDim;
786  int lclDim = MVT::GetNumberVecs(*newstate.V);
787 
788  // check size of H
789  TEUCHOS_TEST_FOR_EXCEPTION(newstate.H->numRows() < curDim_ || newstate.H->numCols() < curDim_, std::invalid_argument, errstr);
790 
791  if (curDim_ == 0 && lclDim > blockSize_) {
792  om_->stream(Warnings) << "Anasazi::BlockKrylovSchur::initialize(): the solver was initialized with a kernel of " << lclDim << std::endl
793  << "The block size however is only " << blockSize_ << std::endl
794  << "The last " << lclDim - blockSize_ << " vectors of the kernel will be overwritten on the first call to iterate()." << std::endl;
795  }
796 
797 
798  // copy basis vectors from newstate into V
799  if (newstate.V != V_) {
800  std::vector<int> nevind(lclDim);
801  for (int i=0; i<lclDim; i++) nevind[i] = i;
802  MVT::SetBlock(*newstate.V,nevind,*V_);
803  }
804 
805  // put data into H_, make sure old information is not still hanging around.
806  if (newstate.H != H_) {
807  H_->putScalar( ST_ZERO );
808  Teuchos::SerialDenseMatrix<int,ScalarType> newH(Teuchos::View,*newstate.H,curDim_+blockSize_,curDim_);
810  lclH = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*H_,curDim_+blockSize_,curDim_) );
811  lclH->assign(newH);
812 
813  // done with local pointers
814  lclH = Teuchos::null;
815  }
816 
817  }
818  else {
819  // user did not specify a basis V
820  // get vectors from problem or generate something, projectAndNormalize, call initialize() recursively
821  Teuchos::RCP<const MV> ivec = problem_->getInitVec();
822  TEUCHOS_TEST_FOR_EXCEPTION(ivec == Teuchos::null,std::invalid_argument,
823  "Anasazi::BlockKrylovSchur::initialize(): eigenproblem did not specify initial vectors to clone from.");
824 
825  int lclDim = MVT::GetNumberVecs(*ivec);
826  bool userand = false;
827  if (lclDim < blockSize_) {
828  // we need at least blockSize_ vectors
829  // use a random multivec
830  userand = true;
831  }
832 
833  if (userand) {
834  // make an index
835  std::vector<int> dimind2(lclDim);
836  for (int i=0; i<lclDim; i++) { dimind2[i] = i; }
837 
838  // alloc newV as a view of the first block of V
839  Teuchos::RCP<MV> newV1 = MVT::CloneViewNonConst(*V_,dimind2);
840 
841  // copy the initial vectors into the first lclDim vectors of V
842  MVT::SetBlock(*ivec,dimind2,*newV1);
843 
844  // resize / reinitialize the index vector
845  dimind2.resize(blockSize_-lclDim);
846  for (int i=0; i<blockSize_-lclDim; i++) { dimind2[i] = lclDim + i; }
847 
848  // initialize the rest of the vectors with random vectors
849  Teuchos::RCP<MV> newV2 = MVT::CloneViewNonConst(*V_,dimind2);
850  MVT::MvRandom(*newV2);
851  }
852  else {
853  // alloc newV as a view of the first block of V
854  Teuchos::RCP<MV> newV1 = MVT::CloneViewNonConst(*V_,bsind);
855 
856  // get a view of the first block of initial vectors
857  Teuchos::RCP<const MV> ivecV = MVT::CloneView(*ivec,bsind);
858 
859  // assign ivec to first part of newV
860  MVT::SetBlock(*ivecV,bsind,*newV1);
861  }
862 
863  // get pointer into first block of V
864  Teuchos::RCP<MV> newV = MVT::CloneViewNonConst(*V_,bsind);
865 
866  // remove auxVecs from newV and normalize newV
867  if (auxVecs_.size() > 0) {
868 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
869  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
870 #endif
871 
873  int rank = orthman_->projectAndNormalize(*newV,auxVecs_);
875  "Anasazi::BlockKrylovSchur::initialize(): couldn't generate initial basis of full rank." );
876  }
877  else {
878 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
879  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
880 #endif
881 
882  int rank = orthman_->normalize(*newV);
884  "Anasazi::BlockKrylovSchur::initialize(): couldn't generate initial basis of full rank." );
885  }
886 
887  // set curDim
888  curDim_ = 0;
889 
890  // clear pointer
891  newV = Teuchos::null;
892  }
893 
894  // The Ritz vectors/values and Schur form are no longer current.
895  ritzVecsCurrent_ = false;
896  ritzValsCurrent_ = false;
897  schurCurrent_ = false;
898 
899  // the solver is initialized
900  initialized_ = true;
901 
902  if (om_->isVerbosity( Debug ) ) {
903  // Check almost everything here
904  CheckList chk;
905  chk.checkV = true;
906  chk.checkArn = true;
907  chk.checkAux = true;
908  om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
909  }
910 
911  // Print information on current status
912  if (om_->isVerbosity(Debug)) {
913  currentStatus( om_->stream(Debug) );
914  }
915  else if (om_->isVerbosity(IterationDetails)) {
916  currentStatus( om_->stream(IterationDetails) );
917  }
918  }
919 
920 
922  // initialize the solver with default state
923  template <class ScalarType, class MV, class OP>
925  {
927  initialize(empty);
928  }
929 
930 
932  // Perform BlockKrylovSchur iterations until the StatusTest tells us to stop.
933  template <class ScalarType, class MV, class OP>
935  {
936  //
937  // Allocate/initialize data structures
938  //
939  if (initialized_ == false) {
940  initialize();
941  }
942 
943  // Compute the current search dimension.
944  // If the problem is non-Hermitian and the blocksize is one, let the solver use the extra vector.
945  int searchDim = blockSize_*numBlocks_;
946  if (problem_->isHermitian() == false) {
947  searchDim++;
948  }
949 
951  // iterate until the status test tells us to stop.
952  //
953  // also break if our basis is full
954  //
955  while (tester_->checkStatus(this) != Passed && curDim_+blockSize_ <= searchDim) {
956 
957  iter_++;
958 
959  // F can be found at the curDim_ block, but the next block is at curDim_ + blockSize_.
960  int lclDim = curDim_ + blockSize_;
961 
962  // Get the current part of the basis.
963  std::vector<int> curind(blockSize_);
964  for (int i=0; i<blockSize_; i++) { curind[i] = lclDim + i; }
965  Teuchos::RCP<MV> Vnext = MVT::CloneViewNonConst(*V_,curind);
966 
967  // Get a view of the previous vectors
968  // this is used for orthogonalization and for computing V^H K H
969  for (int i=0; i<blockSize_; i++) { curind[i] = curDim_ + i; }
970  Teuchos::RCP<const MV> Vprev = MVT::CloneView(*V_,curind);
971  // om_->stream(Debug) << "Vprev: " << std::endl;
972  // MVT::MvPrint(*Vprev,om_->stream(Debug));
973 
974  // Compute the next vector in the Krylov basis: Vnext = Op*Vprev
975  {
976 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
977  Teuchos::TimeMonitor lcltimer( *timerOp_ );
978 #endif
979  OPT::Apply(*Op_,*Vprev,*Vnext);
980  count_ApplyOp_ += blockSize_;
981  }
982  // om_->stream(Debug) << "Vnext: " << std::endl;
983  // MVT::MvPrint(*Vnext,om_->stream(Debug));
984  Vprev = Teuchos::null;
985 
986  // Remove all previous Krylov-Schur basis vectors and auxVecs from Vnext
987  {
988 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
989  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
990 #endif
991 
992  // Get a view of all the previous vectors
993  std::vector<int> prevind(lclDim);
994  for (int i=0; i<lclDim; i++) { prevind[i] = i; }
995  Vprev = MVT::CloneView(*V_,prevind);
996  Teuchos::Array<Teuchos::RCP<const MV> > AVprev(1, Vprev);
997 
998  // Get a view of the part of the Hessenberg matrix needed to hold the ortho coeffs.
1001  ( Teuchos::View,*H_,lclDim,blockSize_,0,curDim_ ) );
1003  AsubH.append( subH );
1004 
1005  // Add the auxiliary vectors to the current basis vectors if any exist
1006  if (auxVecs_.size() > 0) {
1007  for (Array_size_type i=0; i<auxVecs_.size(); i++) {
1008  AVprev.append( auxVecs_[i] );
1009  AsubH.append( Teuchos::null );
1010  }
1011  }
1012 
1013  // Get a view of the part of the Hessenberg matrix needed to hold the norm coeffs.
1014  // om_->stream(Debug) << "V before ortho: " << std::endl;
1015  // MVT::MvPrint(*Vprev,om_->stream(Debug));
1018  ( Teuchos::View,*H_,blockSize_,blockSize_,lclDim,curDim_ ) );
1019  int rank = orthman_->projectAndNormalize(*Vnext,AVprev,AsubH,subR);
1020  // om_->stream(Debug) << "Vnext after ortho: " << std::endl;
1021  // MVT::MvPrint(*Vnext,om_->stream(Debug));
1022  // om_->stream(Debug) << "subH: " << std::endl << *subH << std::endl;
1023  // om_->stream(Debug) << "subR: " << std::endl << *subR << std::endl;
1024  // om_->stream(Debug) << "H: " << std::endl << *H_ << std::endl;
1026  "Anasazi::BlockKrylovSchur::iterate(): couldn't generate basis of full rank.");
1027  }
1028  //
1029  // V has been extended, and H has been extended.
1030  //
1031  // Update basis dim and release all pointers.
1032  Vnext = Teuchos::null;
1033  curDim_ += blockSize_;
1034  // The Ritz vectors/values and Schur form are no longer current.
1035  ritzVecsCurrent_ = false;
1036  ritzValsCurrent_ = false;
1037  schurCurrent_ = false;
1038  //
1039  // Update Ritz values and residuals if needed
1040  if (!(iter_%stepSize_)) {
1041  computeRitzValues();
1042  }
1043 
1044  // When required, monitor some orthogonalities
1045  if (om_->isVerbosity( Debug ) ) {
1046  // Check almost everything here
1047  CheckList chk;
1048  chk.checkV = true;
1049  chk.checkArn = true;
1050  om_->print( Debug, accuracyCheck(chk, ": after local update") );
1051  }
1052  else if (om_->isVerbosity( OrthoDetails ) ) {
1053  CheckList chk;
1054  chk.checkV = true;
1055  om_->print( OrthoDetails, accuracyCheck(chk, ": after local update") );
1056  }
1057 
1058  // Print information on current iteration
1059  if (om_->isVerbosity(Debug)) {
1060  currentStatus( om_->stream(Debug) );
1061  }
1062  else if (om_->isVerbosity(IterationDetails)) {
1063  currentStatus( om_->stream(IterationDetails) );
1064  }
1065 
1066  } // end while (statusTest == false)
1067 
1068  } // end of iterate()
1069 
1070 
1072  // Check accuracy, orthogonality, and other debugging stuff
1073  //
1074  // bools specify which tests we want to run (instead of running more than we actually care about)
1075  //
1076  // checkV : V orthonormal
1077  // orthogonal to auxvecs
1078  // checkAux: check that auxiliary vectors are actually orthonormal
1079  //
1080  // checkArn: check the Arnoldi factorization
1081  //
1082  // NOTE: This method needs to check the current dimension of the subspace, since it is possible to
1083  // call this method when curDim_ = 0 (after initialization).
1084  template <class ScalarType, class MV, class OP>
1085  std::string BlockKrylovSchur<ScalarType,MV,OP>::accuracyCheck( const CheckList &chk, const std::string &where ) const
1086  {
1087  std::stringstream os;
1088  os.precision(2);
1089  os.setf(std::ios::scientific, std::ios::floatfield);
1090  MagnitudeType tmp;
1091 
1092  os << " Debugging checks: iteration " << iter_ << where << std::endl;
1093 
1094  // index vectors for V and F
1095  std::vector<int> lclind(curDim_);
1096  for (int i=0; i<curDim_; i++) lclind[i] = i;
1097  std::vector<int> bsind(blockSize_);
1098  for (int i=0; i<blockSize_; i++) { bsind[i] = curDim_ + i; }
1099 
1100  Teuchos::RCP<const MV> lclV,lclF;
1101  Teuchos::RCP<MV> lclAV;
1102  if (curDim_)
1103  lclV = MVT::CloneView(*V_,lclind);
1104  lclF = MVT::CloneView(*V_,bsind);
1105 
1106  if (chk.checkV) {
1107  if (curDim_) {
1108  tmp = orthman_->orthonormError(*lclV);
1109  os << " >> Error in V^H M V == I : " << tmp << std::endl;
1110  }
1111  tmp = orthman_->orthonormError(*lclF);
1112  os << " >> Error in F^H M F == I : " << tmp << std::endl;
1113  if (curDim_) {
1114  tmp = orthman_->orthogError(*lclV,*lclF);
1115  os << " >> Error in V^H M F == 0 : " << tmp << std::endl;
1116  }
1117  for (Array_size_type i=0; i<auxVecs_.size(); i++) {
1118  if (curDim_) {
1119  tmp = orthman_->orthogError(*lclV,*auxVecs_[i]);
1120  os << " >> Error in V^H M Aux[" << i << "] == 0 : " << tmp << std::endl;
1121  }
1122  tmp = orthman_->orthogError(*lclF,*auxVecs_[i]);
1123  os << " >> Error in F^H M Aux[" << i << "] == 0 : " << tmp << std::endl;
1124  }
1125  }
1126 
1127  if (chk.checkArn) {
1128 
1129  if (curDim_) {
1130  // Compute AV
1131  lclAV = MVT::Clone(*V_,curDim_);
1132  {
1133 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1134  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1135 #endif
1136  OPT::Apply(*Op_,*lclV,*lclAV);
1137  }
1138 
1139  // Compute AV - VH
1141  MVT::MvTimesMatAddMv( -ST_ONE, *lclV, subH, ST_ONE, *lclAV );
1142 
1143  // Compute FB_k^T - (AV-VH)
1145  blockSize_,curDim_, curDim_ );
1146  MVT::MvTimesMatAddMv( -ST_ONE, *lclF, curB, ST_ONE, *lclAV );
1147 
1148  // Compute || FE_k^T - (AV-VH) ||
1149  std::vector<MagnitudeType> arnNorms( curDim_ );
1150  orthman_->norm( *lclAV, arnNorms );
1151 
1152  for (int i=0; i<curDim_; i++) {
1153  os << " >> Error in Krylov-Schur factorization (R = AV-VS-FB^H), ||R[" << i << "]|| : " << arnNorms[i] << std::endl;
1154  }
1155  }
1156  }
1157 
1158  if (chk.checkAux) {
1159  for (Array_size_type i=0; i<auxVecs_.size(); i++) {
1160  tmp = orthman_->orthonormError(*auxVecs_[i]);
1161  os << " >> Error in Aux[" << i << "]^H M Aux[" << i << "] == I : " << tmp << std::endl;
1162  for (Array_size_type j=i+1; j<auxVecs_.size(); j++) {
1163  tmp = orthman_->orthogError(*auxVecs_[i],*auxVecs_[j]);
1164  os << " >> Error in Aux[" << i << "]^H M Aux[" << j << "] == 0 : " << tmp << std::endl;
1165  }
1166  }
1167  }
1168 
1169  os << std::endl;
1170 
1171  return os.str();
1172  }
1173 
1175  /* Get the current approximate eigenvalues, i.e. Ritz values.
1176  *
1177  * POST-CONDITIONS:
1178  *
1179  * ritzValues_ contains Ritz w.r.t. V, H
1180  * Q_ contains the Schur vectors w.r.t. H
1181  * schurH_ contains the Schur matrix w.r.t. H
1182  * ritzOrder_ contains the current ordering from sort manager
1183  */
1184 
1185  template <class ScalarType, class MV, class OP>
1187  {
1188  // Can only call this if the solver is initialized
1189  if (initialized_) {
1190 
1191  // This just updates the Ritz values and residuals.
1192  // --> ritzValsCurrent_ will be set to 'true' by this method.
1193  if (!ritzValsCurrent_) {
1194  // Compute the current Ritz values, through computing the Schur form
1195  // without updating the current projection matrix or sorting the Schur form.
1196  computeSchurForm( false );
1197  }
1198  }
1199  }
1200 
1202  /* Get the current approximate eigenvectors, i.e. Ritz vectors.
1203  *
1204  * POST-CONDITIONS:
1205  *
1206  * ritzValues_ contains Ritz w.r.t. V, H
1207  * ritzVectors_ is first blockSize_ Ritz vectors w.r.t. V, H
1208  * Q_ contains the Schur vectors w.r.t. H
1209  * schurH_ contains the Schur matrix w.r.t. H
1210  * ritzOrder_ contains the current ordering from sort manager
1211  */
1212 
1213  template <class ScalarType, class MV, class OP>
1215  {
1216 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1217  Teuchos::TimeMonitor LocalTimer(*timerCompRitzVec_);
1218 #endif
1219 
1220  TEUCHOS_TEST_FOR_EXCEPTION(numRitzVecs_==0, std::invalid_argument,
1221  "Anasazi::BlockKrylovSchur::computeRitzVectors(): no Ritz vectors were required from this solver.");
1222 
1223  TEUCHOS_TEST_FOR_EXCEPTION(curDim_ < numRitzVecs_, std::invalid_argument,
1224  "Anasazi::BlockKrylovSchur::computeRitzVectors(): the current subspace is not large enough to compute the number of requested Ritz vectors.");
1225 
1226 
1227  // Check to see if the current subspace dimension is non-trivial and the solver is initialized
1228  if (curDim_ && initialized_) {
1229 
1230  // Check to see if the Ritz vectors are current.
1231  if (!ritzVecsCurrent_) {
1232 
1233  // Check to see if the Schur factorization of H (schurH_, Q) is current and sorted.
1234  if (!schurCurrent_) {
1235  // Compute the Schur factorization of the current H, which will not directly change H,
1236  // the factorization will be sorted and placed in (schurH_, Q)
1237  computeSchurForm( true );
1238  }
1239 
1240  // After the Schur form is computed, then the Ritz values are current.
1241  // Thus, I can check the Ritz index vector to see if I have enough space for the Ritz vectors requested.
1242  TEUCHOS_TEST_FOR_EXCEPTION(ritzIndex_[numRitzVecs_-1]==1, std::logic_error,
1243  "Anasazi::BlockKrylovSchur::computeRitzVectors(): the number of required Ritz vectors splits a complex conjugate pair.");
1244 
1247 
1248  // Compute the Ritz vectors.
1249  // --> For a Hermitian problem this is simply the current basis times the first numRitzVecs_ Schur vectors
1250  //
1251  // --> For a non-Hermitian problem, this involves solving the projected eigenproblem, then
1252  // placing the product of the current basis times the first numRitzVecs_ Schur vectors times the
1253  // eigenvectors of interest into the Ritz vectors.
1254 
1255  // Get a view of the current Krylov-Schur basis vectors and Schur vectors
1256  std::vector<int> curind( curDim_ );
1257  for (int i=0; i<curDim_; i++) { curind[i] = i; }
1258  Teuchos::RCP<const MV> Vtemp = MVT::CloneView( *V_, curind );
1259  if (problem_->isHermitian()) {
1260  // Get a view into the current Schur vectors
1261  Teuchos::SerialDenseMatrix<int,ScalarType> subQ( Teuchos::View, *Q_, curDim_, numRitzVecs_ );
1262 
1263  // Compute the current Ritz vectors
1264  MVT::MvTimesMatAddMv( ST_ONE, *Vtemp, subQ, ST_ZERO, *ritzVectors_ );
1265 
1266  // Double check that no complex Ritz values have snuck into the set of converged nev.
1267  bool complexRitzVals = false;
1268  for (int i=0; i<numRitzVecs_; i++) {
1269  if (ritzIndex_[i] != 0) {
1270  complexRitzVals = true;
1271  break;
1272  }
1273  }
1274  if (complexRitzVals)
1275  om_->stream(Warnings) << " Eigenproblem is Hermitian and complex eigenvalues have converged, corresponding eigenvectors will be incorrect!!!"
1276  << std::endl;
1277  } else {
1278 
1279  // Get a view into the current Schur vectors.
1280  Teuchos::SerialDenseMatrix<int,ScalarType> subQ( Teuchos::View, *Q_, curDim_, curDim_ );
1281 
1282  // Get a set of work vectors to hold the current Ritz vectors.
1283  Teuchos::RCP<MV> tmpritzVectors_ = MVT::Clone( *V_, curDim_ );
1284 
1285  // Compute the current Krylov-Schur vectors.
1286  MVT::MvTimesMatAddMv( ST_ONE, *Vtemp, subQ, ST_ZERO, *tmpritzVectors_ );
1287 
1288  // Now compute the eigenvectors of the Schur form
1289  // Reset the dense matrix and compute the eigenvalues of the Schur form.
1290  //
1291  // Allocate the work space. This space will be used below for calls to:
1292  // * TREVC (requires 3*N for real, 2*N for complex)
1293 
1294  int lwork = 3*curDim_;
1295  std::vector<ScalarType> work( lwork );
1296  std::vector<MagnitudeType> rwork( curDim_ );
1297  char side = 'R';
1298  int mm, info = 0;
1299  const int ldvl = 1;
1300  ScalarType vl[ ldvl ];
1301  Teuchos::SerialDenseMatrix<int,ScalarType> copyQ( Teuchos::Copy, *Q_, curDim_, curDim_ );
1302  lapack.TREVC( side, curDim_, schurH_->values(), schurH_->stride(), vl, ldvl,
1303  copyQ.values(), copyQ.stride(), curDim_, &mm, &work[0], &rwork[0], &info );
1304  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::logic_error,
1305  "Anasazi::BlockKrylovSchur::computeRitzVectors(): TREVC(n==" << curDim_ << ") returned info " << info << " != 0.");
1306 
1307  // Get a view into the eigenvectors of the Schur form
1308  Teuchos::SerialDenseMatrix<int,ScalarType> subCopyQ( Teuchos::View, copyQ, curDim_, numRitzVecs_ );
1309 
1310  // Convert back to Ritz vectors of the operator.
1311  curind.resize( numRitzVecs_ ); // This is already initialized above
1312  Teuchos::RCP<MV> view_ritzVectors = MVT::CloneViewNonConst( *ritzVectors_, curind );
1313  MVT::MvTimesMatAddMv( ST_ONE, *tmpritzVectors_, subCopyQ, ST_ZERO, *view_ritzVectors );
1314 
1315  // Compute the norm of the new Ritz vectors
1316  std::vector<MagnitudeType> ritzNrm( numRitzVecs_ );
1317  MVT::MvNorm( *view_ritzVectors, ritzNrm );
1318 
1319  // Release memory used to compute Ritz vectors before scaling the current vectors.
1320  tmpritzVectors_ = Teuchos::null;
1321  view_ritzVectors = Teuchos::null;
1322 
1323  // Scale the Ritz vectors to have Euclidean norm.
1324  ScalarType ritzScale = ST_ONE;
1325  for (int i=0; i<numRitzVecs_; i++) {
1326 
1327  // If this is a conjugate pair then normalize by the real and imaginary parts.
1328  if (ritzIndex_[i] == 1 ) {
1329  ritzScale = ST_ONE/lapack_mag.LAPY2(ritzNrm[i],ritzNrm[i+1]);
1330  std::vector<int> newind(2);
1331  newind[0] = i; newind[1] = i+1;
1332  tmpritzVectors_ = MVT::CloneCopy( *ritzVectors_, newind );
1333  view_ritzVectors = MVT::CloneViewNonConst( *ritzVectors_, newind );
1334  MVT::MvAddMv( ritzScale, *tmpritzVectors_, ST_ZERO, *tmpritzVectors_, *view_ritzVectors );
1335 
1336  // Increment counter for imaginary part
1337  i++;
1338  } else {
1339 
1340  // This is a real Ritz value, normalize the vector
1341  std::vector<int> newind(1);
1342  newind[0] = i;
1343  tmpritzVectors_ = MVT::CloneCopy( *ritzVectors_, newind );
1344  view_ritzVectors = MVT::CloneViewNonConst( *ritzVectors_, newind );
1345  MVT::MvAddMv( ST_ONE/ritzNrm[i], *tmpritzVectors_, ST_ZERO, *tmpritzVectors_, *view_ritzVectors );
1346  }
1347  }
1348 
1349  } // if (problem_->isHermitian())
1350 
1351  // The current Ritz vectors have been computed.
1352  ritzVecsCurrent_ = true;
1353 
1354  } // if (!ritzVecsCurrent_)
1355  } // if (curDim_)
1356  } // computeRitzVectors()
1357 
1358 
1360  // Set a new StatusTest for the solver.
1361  template <class ScalarType, class MV, class OP>
1363  TEUCHOS_TEST_FOR_EXCEPTION(test == Teuchos::null,std::invalid_argument,
1364  "Anasazi::BlockKrylovSchur::setStatusTest() was passed a null StatusTest.");
1365  tester_ = test;
1366  }
1367 
1369  // Get the current StatusTest used by the solver.
1370  template <class ScalarType, class MV, class OP>
1372  return tester_;
1373  }
1374 
1375 
1377  /* Get the current approximate eigenvalues, i.e. Ritz values.
1378  *
1379  * POST-CONDITIONS:
1380  *
1381  * ritzValues_ contains Ritz w.r.t. V, H
1382  * Q_ contains the Schur vectors w.r.t. H
1383  * schurH_ contains the Schur matrix w.r.t. H
1384  * ritzOrder_ contains the current ordering from sort manager
1385  * schurCurrent_ = true if sort = true; i.e. the Schur form is sorted according to the index
1386  * vector returned by the sort manager.
1387  */
1388  template <class ScalarType, class MV, class OP>
1390  {
1391  // local timer
1392 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1393  Teuchos::TimeMonitor LocalTimer(*timerCompSF_);
1394 #endif
1395 
1396  // Check to see if the dimension of the factorization is greater than zero.
1397  if (curDim_) {
1398 
1399  // Check to see if the Schur factorization is current.
1400  if (!schurCurrent_) {
1401 
1402  // Check to see if the Ritz values are current
1403  // --> If they are then the Schur factorization is current but not sorted.
1404  if (!ritzValsCurrent_) {
1409 
1410  // Get a view into Q, the storage for H's Schur vectors.
1411  Teuchos::SerialDenseMatrix<int,ScalarType> subQ( Teuchos::View, *Q_, curDim_, curDim_ );
1412 
1413  // Get a copy of H to compute/sort the Schur form.
1414  schurH_ = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>( Teuchos::Copy, *H_, curDim_, curDim_ ) );
1415  //
1416  //---------------------------------------------------
1417  // Compute the Schur factorization of subH
1418  // ---> Use driver GEES to first reduce to upper Hessenberg
1419  // form and then compute Schur form, outputting Ritz values
1420  //---------------------------------------------------
1421  //
1422  // Allocate the work space. This space will be used below for calls to:
1423  // * GEES (requires 3*N for real, 2*N for complex)
1424  // * TREVC (requires 3*N for real, 2*N for complex)
1425  // * TREXC (requires N for real, none for complex)
1426  // Furthermore, GEES requires a real array of length curDim_ (for complex datatypes)
1427  //
1428  int lwork = 3*curDim_;
1429  std::vector<ScalarType> work( lwork );
1430  std::vector<MagnitudeType> rwork( curDim_ );
1431  std::vector<MagnitudeType> tmp_rRitzValues( curDim_ );
1432  std::vector<MagnitudeType> tmp_iRitzValues( curDim_ );
1433  std::vector<int> bwork( curDim_ );
1434  int info = 0, sdim = 0;
1435  char jobvs = 'V';
1436  lapack.GEES( jobvs,curDim_, schurH_->values(), schurH_->stride(), &sdim, &tmp_rRitzValues[0],
1437  &tmp_iRitzValues[0], subQ.values(), subQ.stride(), &work[0], lwork,
1438  &rwork[0], &bwork[0], &info );
1439 
1440  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::logic_error,
1441  "Anasazi::BlockKrylovSchur::computeSchurForm(): GEES(n==" << curDim_ << ") returned info " << info << " != 0.");
1442 
1443  // Check if imaginary part is detected by the Schur factorization of subH for Hermitian eigenproblems
1444  // NOTE: Because of full orthogonalization, there will be small entries above the block tridiagonal in the block Hessenberg matrix.
1445  // The spectrum of this matrix may include imaginary eigenvalues with small imaginary part, which will mess up the Schur
1446  // form sorting below.
1447  bool hermImagDetected = false;
1448  if (problem_->isHermitian()) {
1449  for (int i=0; i<curDim_; i++)
1450  {
1451  if (tmp_iRitzValues[i] != MT_ZERO)
1452  {
1453  hermImagDetected = true;
1454  break;
1455  }
1456  }
1457  if (hermImagDetected) {
1458  // Warn the user that complex eigenvalues have been detected.
1459  om_->stream(Warnings) << " Eigenproblem is Hermitian and complex eigenvalues have been detected!!!" << std::endl;
1460  // Compute || H - H' || to indicate how bad the symmetry is in the projected eigenproblem
1461  Teuchos::SerialDenseMatrix<int,ScalarType> localH( Teuchos::View, *H_, curDim_, curDim_ );
1465  else
1467  (*tLocalH) -= localH;
1468  MagnitudeType normF = tLocalH->normFrobenius();
1469  MagnitudeType norm1 = tLocalH->normOne();
1470  om_->stream(Warnings) << " Symmetry error in projected eigenproblem: || S - S' ||_F = " << normF
1471  << ", || S - S' ||_1 = " << norm1 << std::endl;
1472  }
1473  }
1474  //
1475  //---------------------------------------------------
1476  // Use the Krylov-Schur factorization to compute the current Ritz residuals
1477  // for ALL the eigenvalues estimates (Ritz values)
1478  // || Ax - x\theta || = || U_m+1*B_m+1^H*Q*s ||
1479  // = || B_m+1^H*Q*s ||
1480  //
1481  // where U_m+1 is the current Krylov-Schur basis, Q are the Schur vectors, and x = U_m+1*Q*s
1482  // NOTE: This means that s = e_i if the problem is hermitian, else the eigenvectors
1483  // of the Schur form need to be computed.
1484  //
1485  // First compute H_{m+1,m}*B_m^T, then determine what 's' is.
1486  //---------------------------------------------------
1487  //
1488  // Get current B_m+1
1490  blockSize_, curDim_, curDim_ );
1491  //
1492  // Compute B_m+1^H*Q
1493  Teuchos::SerialDenseMatrix<int,ScalarType> subB( blockSize_, curDim_ );
1494  blas.GEMM( Teuchos::NO_TRANS, Teuchos::NO_TRANS, blockSize_, curDim_, curDim_, ST_ONE,
1495  curB.values(), curB.stride(), subQ.values(), subQ.stride(),
1496  ST_ZERO, subB.values(), subB.stride() );
1497  //
1498  // Determine what 's' is and compute Ritz residuals.
1499  //
1500  ScalarType* b_ptr = subB.values();
1501  if (problem_->isHermitian() && !hermImagDetected) {
1502  //
1503  // 's' is the i-th canonical basis vector.
1504  //
1505  for (int i=0; i<curDim_ ; i++) {
1506  ritzResiduals_[i] = blas.NRM2(blockSize_, b_ptr + i*blockSize_, 1);
1507  }
1508  } else {
1509  //
1510  // Compute S: the eigenvectors of the block upper triangular, Schur matrix.
1511  //
1512  char side = 'R';
1513  int mm;
1514  const int ldvl = 1;
1515  ScalarType vl[ ldvl ];
1516  Teuchos::SerialDenseMatrix<int,ScalarType> S( curDim_, curDim_ );
1517  lapack.TREVC( side, curDim_, schurH_->values(), schurH_->stride(), vl, ldvl,
1518  S.values(), S.stride(), curDim_, &mm, &work[0], &rwork[0], &info );
1519 
1520  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::logic_error,
1521  "Anasazi::BlockKrylovSchur::computeSchurForm(): TREVC(n==" << curDim_ << ") returned info " << info << " != 0.");
1522  //
1523  // Scale the eigenvectors so that their Euclidean norms are all one.
1524  //
1525  HelperTraits<ScalarType>::scaleRitzVectors( tmp_iRitzValues, &S );
1526  //
1527  // Compute ritzRes = *B_m+1^H*Q*S where the i-th column of S is 's' for the i-th Ritz-value
1528  //
1529  Teuchos::SerialDenseMatrix<int,ScalarType> ritzRes( blockSize_, curDim_ );
1530  blas.GEMM( Teuchos::NO_TRANS, Teuchos::NO_TRANS, blockSize_, curDim_, curDim_, ST_ONE,
1531  subB.values(), subB.stride(), S.values(), S.stride(),
1532  ST_ZERO, ritzRes.values(), ritzRes.stride() );
1533 
1534  /* TO DO: There's be an incorrect assumption made in the computation of the Ritz residuals.
1535  This assumption is that the next vector in the Krylov subspace is Euclidean orthonormal.
1536  It may not be normalized using Euclidean norm.
1537  Teuchos::RCP<MV> ritzResVecs = MVT::Clone( *V_, curDim_ );
1538  std::vector<int> curind(blockSize_);
1539  for (int i=0; i<blockSize_; i++) { curind[i] = curDim_ + i; }
1540  Teuchos::RCP<MV> Vtemp = MVT::CloneView(*V_,curind);
1541 
1542  MVT::MvTimesMatAddMv( ST_ONE, *Vtemp, ritzRes, ST_ZERO, *ritzResVecs );
1543  std::vector<MagnitudeType> ritzResNrms(curDim_);
1544  MVT::MvNorm( *ritzResVecs, ritzResNrms );
1545  i = 0;
1546  while( i < curDim_ ) {
1547  if ( tmp_ritzValues[curDim_+i] != MT_ZERO ) {
1548  ritzResiduals_[i] = lapack_mag.LAPY2( ritzResNrms[i], ritzResNrms[i+1] );
1549  ritzResiduals_[i+1] = ritzResiduals_[i];
1550  i = i+2;
1551  } else {
1552  ritzResiduals_[i] = ritzResNrms[i];
1553  i++;
1554  }
1555  }
1556  */
1557  //
1558  // Compute the Ritz residuals for each Ritz value.
1559  //
1560  HelperTraits<ScalarType>::computeRitzResiduals( tmp_iRitzValues, ritzRes, &ritzResiduals_ );
1561  }
1562  //
1563  // Sort the Ritz values.
1564  //
1565  {
1566 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1567  Teuchos::TimeMonitor LocalTimer2(*timerSortRitzVal_);
1568 #endif
1569  int i=0;
1570  if (problem_->isHermitian() && !hermImagDetected) {
1571  //
1572  // Sort using just the real part of the Ritz values.
1573  sm_->sort(tmp_rRitzValues, Teuchos::rcpFromRef(ritzOrder_), curDim_); // don't catch exception
1574  ritzIndex_.clear();
1575  while ( i < curDim_ ) {
1576  // The Ritz value is not complex.
1577  ritzValues_[i].set(tmp_rRitzValues[i], MT_ZERO);
1578  ritzIndex_.push_back(0);
1579  i++;
1580  }
1581  }
1582  else {
1583  //
1584  // Sort using both the real and imaginary parts of the Ritz values.
1585  sm_->sort(tmp_rRitzValues, tmp_iRitzValues, Teuchos::rcpFromRef(ritzOrder_) , curDim_);
1586  HelperTraits<ScalarType>::sortRitzValues( tmp_rRitzValues, tmp_iRitzValues, &ritzValues_, &ritzOrder_, &ritzIndex_ );
1587  }
1588  //
1589  // Sort the ritzResiduals_ based on the ordering from the Sort Manager.
1590  std::vector<MagnitudeType> ritz2( curDim_ );
1591  for (i=0; i<curDim_; i++) { ritz2[i] = ritzResiduals_[ ritzOrder_[i] ]; }
1592  blas_mag.COPY( curDim_, &ritz2[0], 1, &ritzResiduals_[0], 1 );
1593 
1594  // The Ritz values have now been updated.
1595  ritzValsCurrent_ = true;
1596  }
1597 
1598  } // if (!ritzValsCurrent_) ...
1599  //
1600  //---------------------------------------------------
1601  // * The Ritz values and residuals have been updated at this point.
1602  //
1603  // * The Schur factorization of the projected matrix has been computed,
1604  // and is stored in (schurH_, Q_).
1605  //
1606  // Now the Schur factorization needs to be sorted.
1607  //---------------------------------------------------
1608  //
1609  // Sort the Schur form using the ordering from the Sort Manager.
1610  if (sort) {
1611  sortSchurForm( *schurH_, *Q_, ritzOrder_ );
1612  //
1613  // Indicate the Schur form in (schurH_, Q_) is current and sorted
1614  schurCurrent_ = true;
1615  }
1616  } // if (!schurCurrent_) ...
1617 
1618  } // if (curDim_) ...
1619 
1620  } // computeSchurForm( ... )
1621 
1622 
1624  // Sort the Schur form of H stored in (H,Q) using the ordering vector.
1625  template <class ScalarType, class MV, class OP>
1628  std::vector<int>& order )
1629  {
1630  // local timer
1631 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1632  Teuchos::TimeMonitor LocalTimer(*timerSortSF_);
1633 #endif
1634  //
1635  //---------------------------------------------------
1636  // Reorder real Schur factorization, remember to add one to the indices for the
1637  // fortran call and determine offset. The offset is necessary since the TREXC
1638  // method reorders in a nonsymmetric fashion, thus we use the reordering in
1639  // a stack-like fashion. Also take into account conjugate pairs, which may mess
1640  // up the reordering, since the pair is moved if one of the pair is moved.
1641  //---------------------------------------------------
1642  //
1643  int i = 0, nevtemp = 0;
1644  char compq = 'V';
1645  std::vector<int> offset2( curDim_ );
1646  std::vector<int> order2( curDim_ );
1647 
1648  // LAPACK objects.
1650  int lwork = 3*curDim_;
1651  std::vector<ScalarType> work( lwork );
1652 
1653  while (i < curDim_) {
1654  if ( ritzIndex_[i] != 0 ) { // This is the first value of a complex conjugate pair
1655  offset2[nevtemp] = 0;
1656  for (int j=i; j<curDim_; j++) {
1657  if (order[j] > order[i]) { offset2[nevtemp]++; }
1658  }
1659  order2[nevtemp] = order[i];
1660  i = i+2;
1661  } else {
1662  offset2[nevtemp] = 0;
1663  for (int j=i; j<curDim_; j++) {
1664  if (order[j] > order[i]) { offset2[nevtemp]++; }
1665  }
1666  order2[nevtemp] = order[i];
1667  i++;
1668  }
1669  nevtemp++;
1670  }
1671  ScalarType *ptr_h = H.values();
1672  ScalarType *ptr_q = Q.values();
1673  int ldh = H.stride(), ldq = Q.stride();
1674  int info = 0;
1675  for (i=nevtemp-1; i>=0; i--) {
1676  int ifst = order2[i]+1+offset2[i];
1677  int ilst = 1;
1678  lapack.TREXC( compq, curDim_, ptr_h, ldh, ptr_q, ldq, &ifst,
1679  &ilst, &work[0], &info );
1680  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::logic_error,
1681  "Anasazi::BlockKrylovSchur::computeSchurForm(): TREXC(n==" << curDim_ << ") returned info " << info << " != 0.");
1682  }
1683  }
1684 
1686  // Print the current status of the solver
1687  template <class ScalarType, class MV, class OP>
1689  {
1690  using std::endl;
1691 
1692  os.setf(std::ios::scientific, std::ios::floatfield);
1693  os.precision(6);
1694  os <<"================================================================================" << endl;
1695  os << endl;
1696  os <<" BlockKrylovSchur Solver Status" << endl;
1697  os << endl;
1698  os <<"The solver is "<<(initialized_ ? "initialized." : "not initialized.") << endl;
1699  os <<"The number of iterations performed is " <<iter_<<endl;
1700  os <<"The block size is " << blockSize_<<endl;
1701  os <<"The number of blocks is " << numBlocks_<<endl;
1702  os <<"The current basis size is " << curDim_<<endl;
1703  os <<"The number of auxiliary vectors is " << numAuxVecs_ << endl;
1704  os <<"The number of operations Op*x is "<<count_ApplyOp_<<endl;
1705 
1706  os.setf(std::ios_base::right, std::ios_base::adjustfield);
1707 
1708  os << endl;
1709  if (initialized_) {
1710  os <<"CURRENT RITZ VALUES "<<endl;
1711  if (ritzIndex_.size() != 0) {
1712  int numPrint = (curDim_ < numRitzPrint_? curDim_: numRitzPrint_);
1713  if (problem_->isHermitian()) {
1714  os << std::setw(20) << "Ritz Value"
1715  << std::setw(20) << "Ritz Residual"
1716  << endl;
1717  os <<"--------------------------------------------------------------------------------"<<endl;
1718  for (int i=0; i<numPrint; i++) {
1719  os << std::setw(20) << ritzValues_[i].realpart
1720  << std::setw(20) << ritzResiduals_[i]
1721  << endl;
1722  }
1723  } else {
1724  os << std::setw(24) << "Ritz Value"
1725  << std::setw(30) << "Ritz Residual"
1726  << endl;
1727  os <<"--------------------------------------------------------------------------------"<<endl;
1728  for (int i=0; i<numPrint; i++) {
1729  // Print out the real eigenvalue.
1730  os << std::setw(15) << ritzValues_[i].realpart;
1731  if (ritzValues_[i].imagpart < MT_ZERO) {
1732  os << " - i" << std::setw(15) << Teuchos::ScalarTraits<MagnitudeType>::magnitude(ritzValues_[i].imagpart);
1733  } else {
1734  os << " + i" << std::setw(15) << ritzValues_[i].imagpart;
1735  }
1736  os << std::setw(20) << ritzResiduals_[i] << endl;
1737  }
1738  }
1739  } else {
1740  os << std::setw(20) << "[ NONE COMPUTED ]" << endl;
1741  }
1742  }
1743  os << endl;
1744  os <<"================================================================================" << endl;
1745  os << endl;
1746  }
1747 
1748 } // End of namespace Anasazi
1749 
1750 #endif
1751 
1752 // End of file AnasaziBlockKrylovSchur.hpp
ScalarType * values() const
void computeRitzValues()
Compute the Ritz values using the current Krylov factorization.
BlockKrylovSchurOrthoFailure is thrown when the orthogonalization manager is unable to generate ortho...
Teuchos::Array< Teuchos::RCP< const MV > > getAuxVecs() const
Get the auxiliary vectors for the solver.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > H
The current Hessenberg matrix.
BlockKrylovSchurState< ScalarType, MV > getState() const
Get the current state of the eigensolver.
void setStepSize(int stepSize)
Set the step size.
bool isRitzValsCurrent() const
Get the status of the Ritz values currently stored in the eigensolver.
const Eigenproblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the eigenvalue problem.
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getRitzRes2Norms()
Get the current Ritz residual 2-norms.
Array< T > & append(const T &x)
Teuchos::RCP< const MV > getRitzVectors()
Get the Ritz vectors.
This class defines the interface required by an eigensolver and status test class to compute solution...
virtual ~BlockKrylovSchur()
BlockKrylovSchur destructor.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this eigenproblem.
Declaration of basic traits for the multivector type.
bool isInitialized() const
Indicates whether the solver has been initialized or not.
T & get(const std::string &name, T def_value)
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getResNorms()
Get the current residual norms.
void setAuxVecs(const Teuchos::Array< Teuchos::RCP< const MV > > &auxvecs)
Set the auxiliary vectors for the solver.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void iterate()
This method performs Block Krylov-Schur iterations until the status test indicates the need to stop o...
Virtual base class which defines basic traits for the operator type.
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
static void scaleRitzVectors(const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &iRV, Teuchos::SerialDenseMatrix< int, ScalarType > *S)
Helper function for correctly scaling the eigenvectors of the projected eigenproblem.
Structure to contain pointers to BlockKrylovSchur state variables.
void GEMM(ETransp transa, ETransp transb, const OrdinalType &m, const OrdinalType &n, const OrdinalType &k, const alpha_type alpha, const A_type *A, const OrdinalType &lda, const B_type *B, const OrdinalType &ldb, const beta_type beta, ScalarType *C, const OrdinalType &ldc) const
An exception class parent to all Anasazi exceptions.
void setStatusTest(Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test)
Set a new StatusTest for the solver.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues...
void TREVC(const char &SIDE, const char &HOWMNY, OrdinalType *select, const OrdinalType &n, const ScalarType *T, const OrdinalType &ldt, ScalarType *VL, const OrdinalType &ldvl, ScalarType *VR, const OrdinalType &ldvr, const OrdinalType &mm, OrdinalType *m, ScalarType *WORK, OrdinalType *info) const
int getNumIters() const
Get the current iteration count.
Pure virtual base class which describes the basic interface to the iterative eigensolver.
bool isRitzVecsCurrent() const
Get the status of the Ritz vectors currently stored in the eigensolver.
void initialize()
Initialize the solver with the initial vectors from the eigenproblem or random data.
void computeRitzVectors()
Compute the Ritz vectors using the current Krylov factorization.
int getStepSize() const
Get the step size.
bool isParameter(const std::string &name) const
void setNumRitzVectors(int numRitzVecs)
Set the number of Ritz vectors to compute.
void COPY(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx, ScalarType *y, const OrdinalType &incy) const
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > getStatusTest() const
Get the current StatusTest used by the solver.
void computeSchurForm(const bool sort=true)
Compute the Schur form of the projected eigenproblem from the current Krylov factorization.
static void sortRitzValues(const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &rRV, const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &iRV, std::vector< Value< ScalarType > > *RV, std::vector< int > *RO, std::vector< int > *RI)
Helper function for correctly storing the Ritz values when the eigenproblem is non-Hermitian.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
bool isSchurCurrent() const
Get the status of the Schur form currently stored in the eigensolver.
void GEES(const char &JOBVS, const char &SORT, OrdinalType &(*ptr2func)(ScalarType *, ScalarType *), const OrdinalType &n, ScalarType *A, const OrdinalType &lda, OrdinalType *sdim, ScalarType *WR, ScalarType *WI, ScalarType *VS, const OrdinalType &ldvs, ScalarType *WORK, const OrdinalType &lwork, OrdinalType *BWORK, OrdinalType *info) const
BlockKrylovSchur(const Teuchos::RCP< Eigenproblem< ScalarType, MV, OP > > &problem, const Teuchos::RCP< SortManager< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > &sorter, const Teuchos::RCP< OutputManager< ScalarType > > &printer, const Teuchos::RCP< StatusTest< ScalarType, MV, OP > > &tester, const Teuchos::RCP< OrthoManager< ScalarType, MV > > &ortho, Teuchos::ParameterList &params)
BlockKrylovSchur constructor with eigenproblem, solver utilities, and parameter list of solver option...
Output managers remove the need for the eigensolver to know any information about the required output...
void resetNumIters()
Reset the iteration count.
Templated virtual class for providing orthogonalization/orthonormalization methods.
Traits class which defines basic operations on multivectors.
Virtual base class which defines basic traits for the operator type.
void TREXC(const char &COMPQ, const OrdinalType &n, ScalarType *T, const OrdinalType &ldt, ScalarType *Q, const OrdinalType &ldq, OrdinalType *ifst, OrdinalType *ilst, ScalarType *WORK, OrdinalType *info) const
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace.
int getNumRitzVectors() const
Get the number of Ritz vectors to compute.
ScalarType LAPY2(const ScalarType &x, const ScalarType &y) const
std::vector< Value< ScalarType > > getRitzValues()
Get the Ritz values.
std::vector< int > getRitzIndex()
Get the Ritz index vector.
static magnitudeType magnitude(T a)
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getRes2Norms()
Get the current residual 2-norms.
Teuchos::RCP< const MulVec > V
The current Krylov basis.
BlockKrylovSchurInitFailure is thrown when the BlockKrylovSchur solver is unable to generate an initi...
Types and exceptions used within Anasazi solvers and interfaces.
void currentStatus(std::ostream &os)
This method requests that the solver print out its current status to screen.
Anasazi&#39;s templated virtual class for providing routines for orthogonalization and orthonormalization...
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > Q
The current Schur vectors of the valid part of H.
static void computeRitzResiduals(const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > &iRV, const Teuchos::SerialDenseMatrix< int, ScalarType > &S, std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > *RR)
Helper function for correctly computing the Ritz residuals of the projected eigenproblem.
void setBlockSize(int blockSize)
Set the blocksize.
This class implements the block Krylov-Schur iteration, for solving linear eigenvalue problems...
Anasazi&#39;s templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
void setSize(int blockSize, int numBlocks)
Set the blocksize and number of blocks to be used by the iterative solver in solving this eigenproble...
int curDim
The current dimension of the reduction.
Common interface of stopping criteria for Anasazi&#39;s solvers.
The Eigensolver is a templated virtual base class that defines the basic interface that any eigensolv...
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > S
The current Schur form reduction of the valid part of H.
OrdinalType stride() const