Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnasaziBlockDavidson.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_DAVIDSON_HPP
15 #define ANASAZI_BLOCK_DAVIDSON_HPP
16 
17 #include "AnasaziTypes.hpp"
18 
19 #include "AnasaziEigensolver.hpp"
22 #include "Teuchos_ScalarTraits.hpp"
23 
25 #include "AnasaziSolverUtils.hpp"
26 
27 #include "Teuchos_LAPACK.hpp"
28 #include "Teuchos_BLAS.hpp"
31 #include "Teuchos_TimeMonitor.hpp"
32 
48 namespace Anasazi {
49 
51 
52 
57  template <class ScalarType, class MV>
63  int curDim;
90  BlockDavidsonState() : curDim(0), V(Teuchos::null),
91  X(Teuchos::null), KX(Teuchos::null), MX(Teuchos::null),
92  R(Teuchos::null), H(Teuchos::null),
93  T(Teuchos::null), KK(Teuchos::null) {}
94  };
95 
97 
99 
100 
113  class BlockDavidsonInitFailure : public AnasaziError {public:
114  BlockDavidsonInitFailure(const std::string& what_arg) : AnasaziError(what_arg)
115  {}};
116 
124  class BlockDavidsonOrthoFailure : public AnasaziError {public:
125  BlockDavidsonOrthoFailure(const std::string& what_arg) : AnasaziError(what_arg)
126  {}};
127 
129 
130 
131  template <class ScalarType, class MV, class OP>
132  class BlockDavidson : public Eigensolver<ScalarType,MV,OP> {
133  public:
135 
136 
146  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
149  Teuchos::ParameterList &params
150  );
151 
153  virtual ~BlockDavidson();
155 
156 
158 
159 
183  void iterate();
184 
219 
223  void initialize();
224 
240  bool isInitialized() const;
241 
255 
257 
258 
260 
261 
263  int getNumIters() const;
264 
266  void resetNumIters();
267 
276 
282  std::vector<Value<ScalarType> > getRitzValues();
283 
284 
292  std::vector<int> getRitzIndex();
293 
294 
300  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getResNorms();
301 
302 
308  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRes2Norms();
309 
310 
318  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRitzRes2Norms();
319 
325  int getCurSubspaceDim() const;
326 
328  int getMaxSubspaceDim() const;
329 
331 
332 
334 
335 
338 
341 
344 
354  void setBlockSize(int blockSize);
355 
357  int getBlockSize() const;
358 
371  void setAuxVecs(const Teuchos::Array<Teuchos::RCP<const MV> > &auxvecs);
372 
375 
377 
379 
380 
390  void setSize(int blockSize, int numBlocks);
391 
393 
395 
396 
398  void currentStatus(std::ostream &os);
399 
401 
402  private:
403  //
404  // Convenience typedefs
405  //
410  typedef typename SCT::magnitudeType MagnitudeType;
411  const MagnitudeType ONE;
412  const MagnitudeType ZERO;
413  const MagnitudeType NANVAL;
414  //
415  // Internal structs
416  //
417  struct CheckList {
418  bool checkV;
419  bool checkX, checkMX, checkKX;
420  bool checkH, checkMH, checkKH;
421  bool checkR, checkQ;
422  bool checkKK;
423  CheckList() : checkV(false),
424  checkX(false),checkMX(false),checkKX(false),
425  checkH(false),checkMH(false),checkKH(false),
426  checkR(false),checkQ(false),checkKK(false) {};
427  };
428  //
429  // Internal methods
430  //
431  std::string accuracyCheck(const CheckList &chk, const std::string &where) const;
432  //
433  // Classes inputed through constructor that define the eigenproblem to be solved.
434  //
440  //
441  // Information obtained from the eigenproblem
442  //
446  bool hasM_;
447  //
448  // Internal timers
449  //
450  Teuchos::RCP<Teuchos::Time> timerOp_, timerMOp_, timerPrec_,
451  timerSortEval_, timerDS_,
452  timerLocal_, timerCompRes_,
453  timerOrtho_, timerInit_;
454  //
455  // Counters
456  //
457  int count_ApplyOp_, count_ApplyM_, count_ApplyPrec_;
458 
459  //
460  // Algorithmic parameters.
461  //
462  // blockSize_ is the solver block size; it controls the number of eigenvectors that
463  // we compute, the number of residual vectors that we compute, and therefore the number
464  // of vectors added to the basis on each iteration.
465  int blockSize_;
466  // numBlocks_ is the size of the allocated space for the Krylov basis, in blocks.
467  int numBlocks_;
468 
469  //
470  // Current solver state
471  //
472  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
473  // is capable of running; _initialize is controlled by the initialize() member method
474  // For the implications of the state of initialized_, please see documentation for initialize()
475  bool initialized_;
476  //
477  // curDim_ reflects how much of the current basis is valid
478  // NOTE: 0 <= curDim_ <= blockSize_*numBlocks_
479  // this also tells us how many of the values in theta_ are valid Ritz values
480  int curDim_;
481  //
482  // State Multivecs
483  // H_,KH_,MH_ will not own any storage
484  // H_ will occasionally point at the current block of vectors in the basis V_
485  // MH_,KH_ will occasionally point at MX_,KX_ when they are used as temporary storage
486  Teuchos::RCP<MV> X_, KX_, MX_, R_,
487  H_, KH_, MH_,
488  V_;
489  //
490  // Projected matrices
491  //
493  //
494  // auxiliary vectors
496  int numAuxVecs_;
497  //
498  // Number of iterations that have been performed.
499  int iter_;
500  //
501  // Current eigenvalues, residual norms
502  std::vector<MagnitudeType> theta_, Rnorms_, R2norms_;
503  //
504  // are the residual norms current with the residual?
505  bool Rnorms_current_, R2norms_current_;
506 
507  };
508 
511  //
512  // Implementations
513  //
516 
517 
519  // Constructor
520  template <class ScalarType, class MV, class OP>
524  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
527  Teuchos::ParameterList &params
528  ) :
529  ONE(Teuchos::ScalarTraits<MagnitudeType>::one()),
530  ZERO(Teuchos::ScalarTraits<MagnitudeType>::zero()),
531  NANVAL(Teuchos::ScalarTraits<MagnitudeType>::nan()),
532  // problem, tools
533  problem_(problem),
534  sm_(sorter),
535  om_(printer),
536  tester_(tester),
537  orthman_(ortho),
538  // timers, counters
539 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
540  timerOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Operation Op*x")),
541  timerMOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Operation M*x")),
542  timerPrec_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Operation Prec*x")),
543  timerSortEval_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Sorting eigenvalues")),
544  timerDS_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Direct solve")),
545  timerLocal_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Local update")),
546  timerCompRes_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Computing residuals")),
547  timerOrtho_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Orthogonalization")),
548  timerInit_(Teuchos::TimeMonitor::getNewTimer("Anasazi: BlockDavidson::Initialization")),
549 #endif
550  count_ApplyOp_(0),
551  count_ApplyM_(0),
552  count_ApplyPrec_(0),
553  // internal data
554  blockSize_(0),
555  numBlocks_(0),
556  initialized_(false),
557  curDim_(0),
558  auxVecs_( Teuchos::Array<Teuchos::RCP<const MV> >(0) ),
559  numAuxVecs_(0),
560  iter_(0),
561  Rnorms_current_(false),
562  R2norms_current_(false)
563  {
564  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,std::invalid_argument,
565  "Anasazi::BlockDavidson::constructor: user passed null problem pointer.");
566  TEUCHOS_TEST_FOR_EXCEPTION(sm_ == Teuchos::null,std::invalid_argument,
567  "Anasazi::BlockDavidson::constructor: user passed null sort manager pointer.");
568  TEUCHOS_TEST_FOR_EXCEPTION(om_ == Teuchos::null,std::invalid_argument,
569  "Anasazi::BlockDavidson::constructor: user passed null output manager pointer.");
570  TEUCHOS_TEST_FOR_EXCEPTION(tester_ == Teuchos::null,std::invalid_argument,
571  "Anasazi::BlockDavidson::constructor: user passed null status test pointer.");
572  TEUCHOS_TEST_FOR_EXCEPTION(orthman_ == Teuchos::null,std::invalid_argument,
573  "Anasazi::BlockDavidson::constructor: user passed null orthogonalization manager pointer.");
574  TEUCHOS_TEST_FOR_EXCEPTION(problem_->isProblemSet() == false, std::invalid_argument,
575  "Anasazi::BlockDavidson::constructor: problem is not set.");
576  TEUCHOS_TEST_FOR_EXCEPTION(problem_->isHermitian() == false, std::invalid_argument,
577  "Anasazi::BlockDavidson::constructor: problem is not hermitian.");
578 
579  // get the problem operators
580  Op_ = problem_->getOperator();
581  TEUCHOS_TEST_FOR_EXCEPTION(Op_ == Teuchos::null, std::invalid_argument,
582  "Anasazi::BlockDavidson::constructor: problem provides no operator.");
583  MOp_ = problem_->getM();
584  Prec_ = problem_->getPrec();
585  hasM_ = (MOp_ != Teuchos::null);
586 
587  // set the block size and allocate data
588  int bs = params.get("Block Size", problem_->getNEV());
589  int nb = params.get("Num Blocks", 2);
590  setSize(bs,nb);
591  }
592 
593 
595  // Destructor
596  template <class ScalarType, class MV, class OP>
598 
599 
601  // Set the block size
602  // This simply calls setSize(), modifying the block size while retaining the number of blocks.
603  template <class ScalarType, class MV, class OP>
605  {
606  setSize(blockSize,numBlocks_);
607  }
608 
609 
611  // Return the current auxiliary vectors
612  template <class ScalarType, class MV, class OP>
614  return auxVecs_;
615  }
616 
617 
618 
620  // return the current block size
621  template <class ScalarType, class MV, class OP>
623  return(blockSize_);
624  }
625 
626 
628  // return eigenproblem
629  template <class ScalarType, class MV, class OP>
631  return(*problem_);
632  }
633 
634 
636  // return max subspace dim
637  template <class ScalarType, class MV, class OP>
639  return blockSize_*numBlocks_;
640  }
641 
643  // return current subspace dim
644  template <class ScalarType, class MV, class OP>
646  if (!initialized_) return 0;
647  return curDim_;
648  }
649 
650 
652  // return ritz residual 2-norms
653  template <class ScalarType, class MV, class OP>
654  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
656  return this->getRes2Norms();
657  }
658 
659 
661  // return ritz index
662  template <class ScalarType, class MV, class OP>
664  std::vector<int> ret(curDim_,0);
665  return ret;
666  }
667 
668 
670  // return ritz values
671  template <class ScalarType, class MV, class OP>
672  std::vector<Value<ScalarType> > BlockDavidson<ScalarType,MV,OP>::getRitzValues() {
673  std::vector<Value<ScalarType> > ret(curDim_);
674  for (int i=0; i<curDim_; ++i) {
675  ret[i].realpart = theta_[i];
676  ret[i].imagpart = ZERO;
677  }
678  return ret;
679  }
680 
681 
683  // return pointer to ritz vectors
684  template <class ScalarType, class MV, class OP>
686  return X_;
687  }
688 
689 
691  // reset number of iterations
692  template <class ScalarType, class MV, class OP>
694  iter_=0;
695  }
696 
697 
699  // return number of iterations
700  template <class ScalarType, class MV, class OP>
702  return(iter_);
703  }
704 
705 
707  // return state pointers
708  template <class ScalarType, class MV, class OP>
711  state.curDim = curDim_;
712  state.V = V_;
713  state.X = X_;
714  state.KX = KX_;
715  if (hasM_) {
716  state.MX = MX_;
717  }
718  else {
719  state.MX = Teuchos::null;
720  }
721  state.R = R_;
722  state.H = H_;
723  state.KK = KK_;
724  if (curDim_ > 0) {
725  state.T = Teuchos::rcp(new std::vector<MagnitudeType>(theta_.begin(),theta_.begin()+curDim_) );
726  } else {
727  state.T = Teuchos::rcp(new std::vector<MagnitudeType>(0));
728  }
729  return state;
730  }
731 
732 
734  // Return initialized state
735  template <class ScalarType, class MV, class OP>
736  bool BlockDavidson<ScalarType,MV,OP>::isInitialized() const { return initialized_; }
737 
738 
740  // Set the block size and make necessary adjustments.
741  template <class ScalarType, class MV, class OP>
742  void BlockDavidson<ScalarType,MV,OP>::setSize (int blockSize, int numBlocks)
743  {
744  // time spent here counts towards timerInit_
745 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
746  Teuchos::TimeMonitor initimer( *timerInit_ );
747 #endif
748 
749  // This routine only allocates space; it doesn't not perform any computation
750  // any change in size will invalidate the state of the solver.
751 
752  TEUCHOS_TEST_FOR_EXCEPTION(blockSize < 1, std::invalid_argument, "Anasazi::BlockDavidson::setSize(blocksize,numblocks): blocksize must be strictly positive.");
753  TEUCHOS_TEST_FOR_EXCEPTION(numBlocks < 2, std::invalid_argument, "Anasazi::BlockDavidson::setSize(blocksize,numblocks): numblocks must be greater than one.");
754  if (blockSize == blockSize_ && numBlocks == numBlocks_) {
755  // do nothing
756  return;
757  }
758 
759  blockSize_ = blockSize;
760  numBlocks_ = numBlocks;
761 
763  // grab some Multivector to Clone
764  // in practice, getInitVec() should always provide this, but it is possible to use a
765  // Eigenproblem with nothing in getInitVec() by manually initializing with initialize();
766  // in case of that strange scenario, we will try to Clone from X_ first, then resort to getInitVec()
767  if (X_ != Teuchos::null) { // this is equivalent to blockSize_ > 0
768  tmp = X_;
769  }
770  else {
771  tmp = problem_->getInitVec();
772  TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::invalid_argument,
773  "Anasazi::BlockDavidson::setSize(): eigenproblem did not specify initial vectors to clone from.");
774  }
775 
776  TEUCHOS_TEST_FOR_EXCEPTION(numAuxVecs_+blockSize*static_cast<ptrdiff_t>(numBlocks) > MVT::GetGlobalLength(*tmp),std::invalid_argument,
777  "Anasazi::BlockDavidson::setSize(): max subspace dimension and auxilliary subspace too large.");
778 
779 
781  // blockSize dependent
782  //
783  // grow/allocate vectors
784  Rnorms_.resize(blockSize_,NANVAL);
785  R2norms_.resize(blockSize_,NANVAL);
786  //
787  // clone multivectors off of tmp
788  //
789  // free current allocation first, to make room for new allocation
790  X_ = Teuchos::null;
791  KX_ = Teuchos::null;
792  MX_ = Teuchos::null;
793  R_ = Teuchos::null;
794  V_ = Teuchos::null;
795 
796  om_->print(Debug," >> Allocating X_\n");
797  X_ = MVT::Clone(*tmp,blockSize_);
798  om_->print(Debug," >> Allocating KX_\n");
799  KX_ = MVT::Clone(*tmp,blockSize_);
800  if (hasM_) {
801  om_->print(Debug," >> Allocating MX_\n");
802  MX_ = MVT::Clone(*tmp,blockSize_);
803  }
804  else {
805  MX_ = X_;
806  }
807  om_->print(Debug," >> Allocating R_\n");
808  R_ = MVT::Clone(*tmp,blockSize_);
809 
811  // blockSize*numBlocks dependent
812  //
813  int newsd = blockSize_*numBlocks_;
814  theta_.resize(blockSize_*numBlocks_,NANVAL);
815  om_->print(Debug," >> Allocating V_\n");
816  V_ = MVT::Clone(*tmp,newsd);
818 
819  om_->print(Debug," >> done allocating.\n");
820 
821  initialized_ = false;
822  curDim_ = 0;
823  }
824 
825 
827  // Set the auxiliary vectors
828  template <class ScalarType, class MV, class OP>
830  typedef typename Teuchos::Array<Teuchos::RCP<const MV> >::iterator tarcpmv;
831 
832  // set new auxiliary vectors
833  auxVecs_ = auxvecs;
834  numAuxVecs_ = 0;
835  for (tarcpmv i=auxVecs_.begin(); i != auxVecs_.end(); ++i) {
836  numAuxVecs_ += MVT::GetNumberVecs(**i);
837  }
838 
839  // If the solver has been initialized, V is not necessarily orthogonal to new auxiliary vectors
840  if (numAuxVecs_ > 0 && initialized_) {
841  initialized_ = false;
842  }
843 
844  if (om_->isVerbosity( Debug ) ) {
845  CheckList chk;
846  chk.checkQ = true;
847  om_->print( Debug, accuracyCheck(chk, ": in setAuxVecs()") );
848  }
849  }
850 
851 
853  /* Initialize the state of the solver
854  *
855  * POST-CONDITIONS:
856  *
857  * V_ is orthonormal, orthogonal to auxVecs_, for first curDim_ vectors
858  * theta_ contains Ritz w.r.t. V_(1:curDim_)
859  * X is Ritz vectors w.r.t. V_(1:curDim_)
860  * KX = Op*X
861  * MX = M*X if hasM_
862  * R = KX - MX*diag(theta_)
863  *
864  */
865  template <class ScalarType, class MV, class OP>
867  {
868  // NOTE: memory has been allocated by setBlockSize(). Use setBlock below; do not Clone
869  // NOTE: Overall time spent in this routine is counted to timerInit_; portions will also be counted towards other primitives
870 
871 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
872  Teuchos::TimeMonitor inittimer( *timerInit_ );
873 #endif
874 
875  std::vector<int> bsind(blockSize_);
876  for (int i=0; i<blockSize_; ++i) bsind[i] = i;
877 
879 
880  // in BlockDavidson, V is primary
881  // the order of dependence follows like so.
882  // --init-> V,KK
883  // --ritz analysis-> theta,X
884  // --op apply-> KX,MX
885  // --compute-> R
886  //
887  // if the user specifies all data for a level, we will accept it.
888  // otherwise, we will generate the whole level, and all subsequent levels.
889  //
890  // the data members are ordered based on dependence, and the levels are
891  // partitioned according to the amount of work required to produce the
892  // items in a level.
893  //
894  // inconsistent multivectors widths and lengths will not be tolerated, and
895  // will be treated with exceptions.
896  //
897  // for multivector pointers in newstate which point directly (as opposed to indirectly, via a view) to
898  // multivectors in the solver, the copy will not be affected.
899 
900  // set up V and KK: get them from newstate if user specified them
901  // otherwise, set them manually
902  Teuchos::RCP<MV> lclV;
904 
905  if (newstate.V != Teuchos::null && newstate.KK != Teuchos::null) {
906  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.V) != MVT::GetGlobalLength(*V_), std::invalid_argument,
907  "Anasazi::BlockDavidson::initialize(newstate): Vector length of V not correct." );
908  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim < blockSize_, std::invalid_argument,
909  "Anasazi::BlockDavidson::initialize(newstate): Rank of new state must be at least blockSize().");
910  TEUCHOS_TEST_FOR_EXCEPTION( newstate.curDim > blockSize_*numBlocks_, std::invalid_argument,
911  "Anasazi::BlockDavidson::initialize(newstate): Rank of new state must be less than getMaxSubspaceDim().");
912  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.V) < newstate.curDim, std::invalid_argument,
913  "Anasazi::BlockDavidson::initialize(newstate): Multivector for basis in new state must be as large as specified state rank.");
914 
915  curDim_ = newstate.curDim;
916  // pick an integral amount
917  curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
918 
919  TEUCHOS_TEST_FOR_EXCEPTION( curDim_ != newstate.curDim, std::invalid_argument,
920  "Anasazi::BlockDavidson::initialize(newstate): Rank of new state must be a multiple of getBlockSize().");
921 
922  // check size of KK
923  TEUCHOS_TEST_FOR_EXCEPTION( newstate.KK->numRows() < curDim_ || newstate.KK->numCols() < curDim_, std::invalid_argument,
924  "Anasazi::BlockDavidson::initialize(newstate): Projected matrix in new state must be as large as specified state rank.");
925 
926  // put data in V
927  std::vector<int> nevind(curDim_);
928  for (int i=0; i<curDim_; ++i) nevind[i] = i;
929  if (newstate.V != V_) {
930  if (curDim_ < MVT::GetNumberVecs(*newstate.V)) {
931  newstate.V = MVT::CloneView(*newstate.V,nevind);
932  }
933  MVT::SetBlock(*newstate.V,nevind,*V_);
934  }
935  lclV = MVT::CloneViewNonConst(*V_,nevind);
936 
937  // put data into KK_
938  lclKK = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
939  if (newstate.KK != KK_) {
940  if (newstate.KK->numRows() > curDim_ || newstate.KK->numCols() > curDim_) {
941  newstate.KK = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*newstate.KK,curDim_,curDim_) );
942  }
943  lclKK->assign(*newstate.KK);
944  }
945  //
946  // make lclKK Hermitian in memory (copy the upper half to the lower half)
947  for (int j=0; j<curDim_-1; ++j) {
948  for (int i=j+1; i<curDim_; ++i) {
949  (*lclKK)(i,j) = SCT::conjugate((*lclKK)(j,i));
950  }
951  }
952  }
953  else {
954  // user did not specify one of V or KK
955  // get vectors from problem or generate something, projectAndNormalize
956  Teuchos::RCP<const MV> ivec = problem_->getInitVec();
957  TEUCHOS_TEST_FOR_EXCEPTION(ivec == Teuchos::null,std::invalid_argument,
958  "Anasazi::BlockDavdison::initialize(newstate): Eigenproblem did not specify initial vectors to clone from.");
959  // clear newstate so we won't use any data from it below
960  newstate.X = Teuchos::null;
961  newstate.MX = Teuchos::null;
962  newstate.KX = Teuchos::null;
963  newstate.R = Teuchos::null;
964  newstate.H = Teuchos::null;
965  newstate.T = Teuchos::null;
966  newstate.KK = Teuchos::null;
967  newstate.V = Teuchos::null;
968  newstate.curDim = 0;
969 
970  curDim_ = MVT::GetNumberVecs(*ivec);
971  // pick the largest multiple of blockSize_
972  curDim_ = (int)(curDim_ / blockSize_)*blockSize_;
973  if (curDim_ > blockSize_*numBlocks_) {
974  // user specified too many vectors... truncate
975  // this produces a full subspace, but that is okay
976  curDim_ = blockSize_*numBlocks_;
977  }
978  bool userand = false;
979  if (curDim_ == 0) {
980  // we need at least blockSize_ vectors
981  // use a random multivec: ignore everything from InitVec
982  userand = true;
983  curDim_ = blockSize_;
984  }
985 
986  // get pointers into V,KV,MV
987  // tmpVecs will be used below for M*V and K*V (not simultaneously)
988  // lclV has curDim vectors
989  // if there is space for lclV and tmpVecs in V_, point tmpVecs into V_
990  // otherwise, we must allocate space for these products
991  //
992  // get pointer to first curDim vector in V_
993  std::vector<int> dimind(curDim_);
994  for (int i=0; i<curDim_; ++i) dimind[i] = i;
995  lclV = MVT::CloneViewNonConst(*V_,dimind);
996  if (userand) {
997  // generate random vector data
998  MVT::MvRandom(*lclV);
999  }
1000  else {
1001  if (MVT::GetNumberVecs(*ivec) > curDim_) {
1002  ivec = MVT::CloneView(*ivec,dimind);
1003  }
1004  // assign ivec to first part of lclV
1005  MVT::SetBlock(*ivec,dimind,*lclV);
1006  }
1007  Teuchos::RCP<MV> tmpVecs;
1008  if (curDim_*2 <= blockSize_*numBlocks_) {
1009  // partition V_ = [lclV tmpVecs _leftover_]
1010  std::vector<int> block2(curDim_);
1011  for (int i=0; i<curDim_; ++i) block2[i] = curDim_+i;
1012  tmpVecs = MVT::CloneViewNonConst(*V_,block2);
1013  }
1014  else {
1015  // allocate space for tmpVecs
1016  tmpVecs = MVT::Clone(*V_,curDim_);
1017  }
1018 
1019  // compute M*lclV if hasM_
1020  if (hasM_) {
1021 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1022  Teuchos::TimeMonitor lcltimer( *timerMOp_ );
1023 #endif
1024  OPT::Apply(*MOp_,*lclV,*tmpVecs);
1025  count_ApplyM_ += curDim_;
1026  }
1027 
1028  // remove auxVecs from lclV and normalize it
1029  if (auxVecs_.size() > 0) {
1030 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1031  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
1032 #endif
1033 
1035  int rank = orthman_->projectAndNormalizeMat(*lclV,auxVecs_,dummyC,Teuchos::null,tmpVecs);
1037  "Anasazi::BlockDavidson::initialize(): Couldn't generate initial basis of full rank.");
1038  }
1039  else {
1040 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1041  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
1042 #endif
1043 
1044  int rank = orthman_->normalizeMat(*lclV,Teuchos::null,tmpVecs);
1046  "Anasazi::BlockDavidson::initialize(): Couldn't generate initial basis of full rank.");
1047  }
1048 
1049  // compute K*lclV: we are re-using tmpVecs to store the result
1050  {
1051 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1052  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1053 #endif
1054  OPT::Apply(*Op_,*lclV,*tmpVecs);
1055  count_ApplyOp_ += curDim_;
1056  }
1057 
1058  // generate KK
1059  lclKK = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,curDim_) );
1060  MVT::MvTransMv(ONE,*lclV,*tmpVecs,*lclKK);
1061 
1062  // clear tmpVecs
1063  tmpVecs = Teuchos::null;
1064  }
1065 
1066  // X,theta require Ritz analysis; if we have to generate one of these, we might as well generate both
1067  if (newstate.X != Teuchos::null && newstate.T != Teuchos::null) {
1068  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.X) != blockSize_ || MVT::GetGlobalLength(*newstate.X) != MVT::GetGlobalLength(*X_),
1069  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): Size of X must be consistent with block size and length of V.");
1070  TEUCHOS_TEST_FOR_EXCEPTION((signed int)(newstate.T->size()) != curDim_,
1071  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): Size of T must be consistent with dimension of V.");
1072 
1073  if (newstate.X != X_) {
1074  MVT::SetBlock(*newstate.X,bsind,*X_);
1075  }
1076 
1077  std::copy(newstate.T->begin(),newstate.T->end(),theta_.begin());
1078  }
1079  else {
1080  // compute ritz vecs/vals
1081  Teuchos::SerialDenseMatrix<int,ScalarType> S(curDim_,curDim_);
1082  {
1083 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1084  Teuchos::TimeMonitor lcltimer( *timerDS_ );
1085 #endif
1086  int rank = curDim_;
1087  Utils::directSolver(curDim_, *lclKK, Teuchos::null, S, theta_, rank, 10);
1088  // we want all ritz values back
1090  "Anasazi::BlockDavidson::initialize(newstate): Not enough Ritz vectors to initialize algorithm.");
1091  }
1092  // sort ritz pairs
1093  {
1094 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1095  Teuchos::TimeMonitor lcltimer( *timerSortEval_ );
1096 #endif
1097 
1098  std::vector<int> order(curDim_);
1099  //
1100  // sort the first curDim_ values in theta_
1101  sm_->sort(theta_, Teuchos::rcpFromRef(order), curDim_); // don't catch exception
1102  //
1103  // apply the same ordering to the primitive ritz vectors
1104  Utils::permuteVectors(order,S);
1105  }
1106 
1107  // compute eigenvectors
1109  {
1110 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1111  Teuchos::TimeMonitor lcltimer( *timerLocal_ );
1112 #endif
1113 
1114  // X <- lclV*S
1115  MVT::MvTimesMatAddMv( ONE, *lclV, S1, ZERO, *X_ );
1116  }
1117  // we generated theta,X so we don't want to use the user's KX,MX
1118  newstate.KX = Teuchos::null;
1119  newstate.MX = Teuchos::null;
1120  }
1121 
1122  // done with local pointers
1123  lclV = Teuchos::null;
1124  lclKK = Teuchos::null;
1125 
1126  // set up KX
1127  if ( newstate.KX != Teuchos::null ) {
1128  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.KX) != blockSize_,
1129  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): vector length of newstate.KX not correct." );
1130  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetGlobalLength(*newstate.KX) != MVT::GetGlobalLength(*X_),
1131  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): newstate.KX must have at least block size vectors." );
1132  if (newstate.KX != KX_) {
1133  MVT::SetBlock(*newstate.KX,bsind,*KX_);
1134  }
1135  }
1136  else {
1137  // generate KX
1138  {
1139 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1140  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1141 #endif
1142  OPT::Apply(*Op_,*X_,*KX_);
1143  count_ApplyOp_ += blockSize_;
1144  }
1145  // we generated KX; we will generate R as well
1146  newstate.R = Teuchos::null;
1147  }
1148 
1149  // set up MX
1150  if (hasM_) {
1151  if ( newstate.MX != Teuchos::null ) {
1152  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.MX) != blockSize_,
1153  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): vector length of newstate.MX not correct." );
1154  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetGlobalLength(*newstate.MX) != MVT::GetGlobalLength(*X_),
1155  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): newstate.MX must have at least block size vectors." );
1156  if (newstate.MX != MX_) {
1157  MVT::SetBlock(*newstate.MX,bsind,*MX_);
1158  }
1159  }
1160  else {
1161  // generate MX
1162  {
1163 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1164  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1165 #endif
1166  OPT::Apply(*MOp_,*X_,*MX_);
1167  count_ApplyOp_ += blockSize_;
1168  }
1169  // we generated MX; we will generate R as well
1170  newstate.R = Teuchos::null;
1171  }
1172  }
1173  else {
1174  // the assignment MX_==X_ would be redundant; take advantage of this opportunity to debug a little
1175  TEUCHOS_TEST_FOR_EXCEPTION(MX_ != X_, std::logic_error, "Anasazi::BlockDavidson::initialize(): solver invariant not satisfied (MX==X).");
1176  }
1177 
1178  // set up R
1179  if (newstate.R != Teuchos::null) {
1180  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetNumberVecs(*newstate.R) != blockSize_,
1181  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): vector length of newstate.R not correct." );
1182  TEUCHOS_TEST_FOR_EXCEPTION(MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*X_),
1183  std::invalid_argument, "Anasazi::BlockDavidson::initialize(newstate): newstate.R must have at least block size vectors." );
1184  if (newstate.R != R_) {
1185  MVT::SetBlock(*newstate.R,bsind,*R_);
1186  }
1187  }
1188  else {
1189 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1190  Teuchos::TimeMonitor lcltimer( *timerCompRes_ );
1191 #endif
1192 
1193  // form R <- KX - MX*T
1194  MVT::MvAddMv(ZERO,*KX_,ONE,*KX_,*R_);
1195  Teuchos::SerialDenseMatrix<int,ScalarType> T(blockSize_,blockSize_);
1196  T.putScalar(ZERO);
1197  for (int i=0; i<blockSize_; ++i) T(i,i) = theta_[i];
1198  MVT::MvTimesMatAddMv(-ONE,*MX_,T,ONE,*R_);
1199 
1200  }
1201 
1202  // R has been updated; mark the norms as out-of-date
1203  Rnorms_current_ = false;
1204  R2norms_current_ = false;
1205 
1206  // finally, we are initialized
1207  initialized_ = true;
1208 
1209  if (om_->isVerbosity( Debug ) ) {
1210  // Check almost everything here
1211  CheckList chk;
1212  chk.checkV = true;
1213  chk.checkX = true;
1214  chk.checkKX = true;
1215  chk.checkMX = true;
1216  chk.checkR = true;
1217  chk.checkQ = true;
1218  chk.checkKK = true;
1219  om_->print( Debug, accuracyCheck(chk, ": after initialize()") );
1220  }
1221 
1222  // Print information on current status
1223  if (om_->isVerbosity(Debug)) {
1224  currentStatus( om_->stream(Debug) );
1225  }
1226  else if (om_->isVerbosity(IterationDetails)) {
1227  currentStatus( om_->stream(IterationDetails) );
1228  }
1229  }
1230 
1231 
1233  // initialize the solver with default state
1234  template <class ScalarType, class MV, class OP>
1236  {
1238  initialize(empty);
1239  }
1240 
1241 
1243  // Perform BlockDavidson iterations until the StatusTest tells us to stop.
1244  template <class ScalarType, class MV, class OP>
1246  {
1247  //
1248  // Initialize solver state
1249  if (initialized_ == false) {
1250  initialize();
1251  }
1252 
1253  // as a data member, this would be redundant and require synchronization with
1254  // blockSize_ and numBlocks_; we'll use a constant here.
1255  const int searchDim = blockSize_*numBlocks_;
1256 
1258 
1259  //
1260  // The projected matrices are part of the state, but the eigenvectors are defined locally.
1261  // S = Local eigenvectors (size: searchDim * searchDim
1262  Teuchos::SerialDenseMatrix<int,ScalarType> S( searchDim, searchDim );
1263 
1264 
1266  // iterate until the status test tells us to stop.
1267  // also break if our basis is full
1268  while (tester_->checkStatus(this) != Passed && curDim_ < searchDim) {
1269 
1270  // Print information on current iteration
1271  if (om_->isVerbosity(Debug)) {
1272  currentStatus( om_->stream(Debug) );
1273  }
1274  else if (om_->isVerbosity(IterationDetails)) {
1275  currentStatus( om_->stream(IterationDetails) );
1276  }
1277 
1278  ++iter_;
1279 
1280  // get the current part of the basis
1281  std::vector<int> curind(blockSize_);
1282  for (int i=0; i<blockSize_; ++i) curind[i] = curDim_ + i;
1283  H_ = MVT::CloneViewNonConst(*V_,curind);
1284 
1285  // Apply the preconditioner on the residuals: H <- Prec*R
1286  // H = Prec*R
1287  if (Prec_ != Teuchos::null) {
1288 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1289  Teuchos::TimeMonitor lcltimer( *timerPrec_ );
1290 #endif
1291  OPT::Apply( *Prec_, *R_, *H_ ); // don't catch the exception
1292  count_ApplyPrec_ += blockSize_;
1293  }
1294  else {
1295  std::vector<int> bsind(blockSize_);
1296  for (int i=0; i<blockSize_; ++i) { bsind[i] = i; }
1297  MVT::SetBlock(*R_,bsind,*H_);
1298  }
1299 
1300  // Apply the mass matrix on H
1301  if (hasM_) {
1302  // use memory at MX_ for temporary storage
1303  MH_ = MX_;
1304 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1305  Teuchos::TimeMonitor lcltimer( *timerMOp_ );
1306 #endif
1307  OPT::Apply( *MOp_, *H_, *MH_); // don't catch the exception
1308  count_ApplyM_ += blockSize_;
1309  }
1310  else {
1311  MH_ = H_;
1312  }
1313 
1314  // Get a view of the previous vectors
1315  // this is used for orthogonalization and for computing V^H K H
1316  std::vector<int> prevind(curDim_);
1317  for (int i=0; i<curDim_; ++i) prevind[i] = i;
1318  Teuchos::RCP<const MV> Vprev = MVT::CloneView(*V_,prevind);
1319 
1320  // Orthogonalize H against the previous vectors and the auxiliary vectors, and normalize
1321  {
1322 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1323  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
1324 #endif
1325 
1326  Teuchos::Array<Teuchos::RCP<const MV> > against = auxVecs_;
1327  against.push_back(Vprev);
1329  int rank = orthman_->projectAndNormalizeMat(*H_,against,
1330  dummyC,
1331  Teuchos::null,MH_);
1333  "Anasazi::BlockDavidson::iterate(): unable to compute orthonormal basis for H.");
1334  }
1335 
1336  // Apply the stiffness matrix to H
1337  {
1338  // use memory at KX_ for temporary storage
1339  KH_ = KX_;
1340 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1341  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1342 #endif
1343  OPT::Apply( *Op_, *H_, *KH_); // don't catch the exception
1344  count_ApplyOp_ += blockSize_;
1345  }
1346 
1347  if (om_->isVerbosity( Debug ) ) {
1348  CheckList chk;
1349  chk.checkH = true;
1350  chk.checkMH = true;
1351  chk.checkKH = true;
1352  om_->print( Debug, accuracyCheck(chk, ": after ortho H") );
1353  }
1354  else if (om_->isVerbosity( OrthoDetails ) ) {
1355  CheckList chk;
1356  chk.checkH = true;
1357  chk.checkMH = true;
1358  chk.checkKH = true;
1359  om_->print( OrthoDetails, accuracyCheck(chk,": after ortho H") );
1360  }
1361 
1362  // compute next part of the projected matrices
1363  // this this in two parts
1365  // Vprev*K*H
1366  nextKK = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,curDim_,blockSize_,0,curDim_) );
1367  MVT::MvTransMv(ONE,*Vprev,*KH_,*nextKK);
1368  // H*K*H
1369  nextKK = Teuchos::rcp( new Teuchos::SerialDenseMatrix<int,ScalarType>(Teuchos::View,*KK_,blockSize_,blockSize_,curDim_,curDim_) );
1370  MVT::MvTransMv(ONE,*H_,*KH_,*nextKK);
1371  //
1372  // make sure that KK_ is Hermitian in memory
1373  nextKK = Teuchos::null;
1374  for (int i=curDim_; i<curDim_+blockSize_; ++i) {
1375  for (int j=0; j<i; ++j) {
1376  (*KK_)(i,j) = SCT::conjugate((*KK_)(j,i));
1377  }
1378  }
1379 
1380  // V has been extended, and KK has been extended. Update basis dim and release all pointers.
1381  curDim_ += blockSize_;
1382  H_ = KH_ = MH_ = Teuchos::null;
1383  Vprev = Teuchos::null;
1384 
1385  if (om_->isVerbosity( Debug ) ) {
1386  CheckList chk;
1387  chk.checkKK = true;
1388  om_->print( Debug, accuracyCheck(chk, ": after expanding KK") );
1389  }
1390 
1391  // Get pointer to complete basis
1392  curind.resize(curDim_);
1393  for (int i=0; i<curDim_; ++i) curind[i] = i;
1394  Teuchos::RCP<const MV> curV = MVT::CloneView(*V_,curind);
1395 
1396  // Perform spectral decomposition
1397  {
1398 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1399  Teuchos::TimeMonitor lcltimer(*timerDS_);
1400 #endif
1401  int nevlocal = curDim_;
1402  int info = Utils::directSolver(curDim_,*KK_,Teuchos::null,S,theta_,nevlocal,10);
1403  TEUCHOS_TEST_FOR_EXCEPTION(info != 0,std::logic_error,"Anasazi::BlockDavidson::iterate(): direct solve returned error code.");
1404  // we did not ask directSolver to perform deflation, so nevLocal better be curDim_
1405  TEUCHOS_TEST_FOR_EXCEPTION(nevlocal != curDim_,std::logic_error,"Anasazi::BlockDavidson::iterate(): direct solve did not compute all eigenvectors."); // this should never happen
1406  }
1407 
1408  // Sort ritz pairs
1409  {
1410 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1411  Teuchos::TimeMonitor lcltimer( *timerSortEval_ );
1412 #endif
1413 
1414  std::vector<int> order(curDim_);
1415  //
1416  // sort the first curDim_ values in theta_
1417  sm_->sort(theta_, Teuchos::rcp(&order,false), curDim_); // don't catch exception
1418  //
1419  // apply the same ordering to the primitive ritz vectors
1421  Utils::permuteVectors(order,curS);
1422  }
1423 
1424  // Create a view matrix of the first blockSize_ vectors
1425  Teuchos::SerialDenseMatrix<int,ScalarType> S1( Teuchos::View, S, curDim_, blockSize_ );
1426 
1427  // Compute the new Ritz vectors
1428  {
1429 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1430  Teuchos::TimeMonitor lcltimer( *timerLocal_ );
1431 #endif
1432  MVT::MvTimesMatAddMv(ONE,*curV,S1,ZERO,*X_);
1433  }
1434 
1435  // Apply the stiffness matrix for the Ritz vectors
1436  {
1437 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1438  Teuchos::TimeMonitor lcltimer( *timerOp_ );
1439 #endif
1440  OPT::Apply( *Op_, *X_, *KX_); // don't catch the exception
1441  count_ApplyOp_ += blockSize_;
1442  }
1443  // Apply the mass matrix for the Ritz vectors
1444  if (hasM_) {
1445 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1446  Teuchos::TimeMonitor lcltimer( *timerMOp_ );
1447 #endif
1448  OPT::Apply(*MOp_,*X_,*MX_);
1449  count_ApplyM_ += blockSize_;
1450  }
1451  else {
1452  MX_ = X_;
1453  }
1454 
1455  // Compute the residual
1456  // R = KX - MX*diag(theta)
1457  {
1458 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1459  Teuchos::TimeMonitor lcltimer( *timerCompRes_ );
1460 #endif
1461 
1462  MVT::MvAddMv( ONE, *KX_, ZERO, *KX_, *R_ );
1463  Teuchos::SerialDenseMatrix<int,ScalarType> T( blockSize_, blockSize_ );
1464  for (int i = 0; i < blockSize_; ++i) {
1465  T(i,i) = theta_[i];
1466  }
1467  MVT::MvTimesMatAddMv( -ONE, *MX_, T, ONE, *R_ );
1468  }
1469 
1470  // R has been updated; mark the norms as out-of-date
1471  Rnorms_current_ = false;
1472  R2norms_current_ = false;
1473 
1474 
1475  // When required, monitor some orthogonalities
1476  if (om_->isVerbosity( Debug ) ) {
1477  // Check almost everything here
1478  CheckList chk;
1479  chk.checkV = true;
1480  chk.checkX = true;
1481  chk.checkKX = true;
1482  chk.checkMX = true;
1483  chk.checkR = true;
1484  om_->print( Debug, accuracyCheck(chk, ": after local update") );
1485  }
1486  else if (om_->isVerbosity( OrthoDetails )) {
1487  CheckList chk;
1488  chk.checkX = true;
1489  chk.checkKX = true;
1490  chk.checkMX = true;
1491  chk.checkR = true;
1492  om_->print( OrthoDetails, accuracyCheck(chk, ": after local update") );
1493  }
1494  } // end while (statusTest == false)
1495 
1496  } // end of iterate()
1497 
1498 
1499 
1501  // compute/return residual M-norms
1502  template <class ScalarType, class MV, class OP>
1503  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
1505  if (Rnorms_current_ == false) {
1506  // Update the residual norms
1507  orthman_->norm(*R_,Rnorms_);
1508  Rnorms_current_ = true;
1509  }
1510  return Rnorms_;
1511  }
1512 
1513 
1515  // compute/return residual 2-norms
1516  template <class ScalarType, class MV, class OP>
1517  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
1519  if (R2norms_current_ == false) {
1520  // Update the residual 2-norms
1521  MVT::MvNorm(*R_,R2norms_);
1522  R2norms_current_ = true;
1523  }
1524  return R2norms_;
1525  }
1526 
1528  // Set a new StatusTest for the solver.
1529  template <class ScalarType, class MV, class OP>
1531  TEUCHOS_TEST_FOR_EXCEPTION(test == Teuchos::null,std::invalid_argument,
1532  "Anasazi::BlockDavidson::setStatusTest() was passed a null StatusTest.");
1533  tester_ = test;
1534  }
1535 
1537  // Get the current StatusTest used by the solver.
1538  template <class ScalarType, class MV, class OP>
1540  return tester_;
1541  }
1542 
1543 
1545  // Check accuracy, orthogonality, and other debugging stuff
1546  //
1547  // bools specify which tests we want to run (instead of running more than we actually care about)
1548  //
1549  // we don't bother checking the following because they are computed explicitly:
1550  // H == Prec*R
1551  // KH == K*H
1552  //
1553  //
1554  // checkV : V orthonormal
1555  // orthogonal to auxvecs
1556  // checkX : X orthonormal
1557  // orthogonal to auxvecs
1558  // checkMX: check MX == M*X
1559  // checkKX: check KX == K*X
1560  // checkH : H orthonormal
1561  // orthogonal to V and H and auxvecs
1562  // checkMH: check MH == M*H
1563  // checkR : check R orthogonal to X
1564  // checkQ : check that auxiliary vectors are actually orthonormal
1565  // checkKK: check that KK is symmetric in memory
1566  //
1567  // TODO:
1568  // add checkTheta
1569  //
1570  template <class ScalarType, class MV, class OP>
1571  std::string BlockDavidson<ScalarType,MV,OP>::accuracyCheck( const CheckList &chk, const std::string &where ) const
1572  {
1573  using std::endl;
1574 
1575  std::stringstream os;
1576  os.precision(2);
1577  os.setf(std::ios::scientific, std::ios::floatfield);
1578 
1579  os << " Debugging checks: iteration " << iter_ << where << endl;
1580 
1581  // V and friends
1582  std::vector<int> lclind(curDim_);
1583  for (int i=0; i<curDim_; ++i) lclind[i] = i;
1585  if (initialized_) {
1586  lclV = MVT::CloneView(*V_,lclind);
1587  }
1588  if (chk.checkV && initialized_) {
1589  MagnitudeType err = orthman_->orthonormError(*lclV);
1590  os << " >> Error in V^H M V == I : " << err << endl;
1591  for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
1592  err = orthman_->orthogError(*lclV,*auxVecs_[i]);
1593  os << " >> Error in V^H M Q[" << i << "] == 0 : " << err << endl;
1594  }
1595  Teuchos::SerialDenseMatrix<int,ScalarType> curKK(curDim_,curDim_);
1596  Teuchos::RCP<MV> lclKV = MVT::Clone(*V_,curDim_);
1597  OPT::Apply(*Op_,*lclV,*lclKV);
1598  MVT::MvTransMv(ONE,*lclV,*lclKV,curKK);
1599  Teuchos::SerialDenseMatrix<int,ScalarType> subKK(Teuchos::View,*KK_,curDim_,curDim_);
1600  curKK -= subKK;
1601  // dup the lower tri part
1602  for (int j=0; j<curDim_; ++j) {
1603  for (int i=j+1; i<curDim_; ++i) {
1604  curKK(i,j) = curKK(j,i);
1605  }
1606  }
1607  os << " >> Error in V^H K V == KK : " << curKK.normFrobenius() << endl;
1608  }
1609 
1610  // X and friends
1611  if (chk.checkX && initialized_) {
1612  MagnitudeType err = orthman_->orthonormError(*X_);
1613  os << " >> Error in X^H M X == I : " << err << endl;
1614  for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
1615  err = orthman_->orthogError(*X_,*auxVecs_[i]);
1616  os << " >> Error in X^H M Q[" << i << "] == 0 : " << err << endl;
1617  }
1618  }
1619  if (chk.checkMX && hasM_ && initialized_) {
1620  MagnitudeType err = Utils::errorEquality(*X_, *MX_, MOp_);
1621  os << " >> Error in MX == M*X : " << err << endl;
1622  }
1623  if (chk.checkKX && initialized_) {
1624  MagnitudeType err = Utils::errorEquality(*X_, *KX_, Op_);
1625  os << " >> Error in KX == K*X : " << err << endl;
1626  }
1627 
1628  // H and friends
1629  if (chk.checkH && initialized_) {
1630  MagnitudeType err = orthman_->orthonormError(*H_);
1631  os << " >> Error in H^H M H == I : " << err << endl;
1632  err = orthman_->orthogError(*H_,*lclV);
1633  os << " >> Error in H^H M V == 0 : " << err << endl;
1634  err = orthman_->orthogError(*H_,*X_);
1635  os << " >> Error in H^H M X == 0 : " << err << endl;
1636  for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
1637  err = orthman_->orthogError(*H_,*auxVecs_[i]);
1638  os << " >> Error in H^H M Q[" << i << "] == 0 : " << err << endl;
1639  }
1640  }
1641  if (chk.checkKH && initialized_) {
1642  MagnitudeType err = Utils::errorEquality(*H_, *KH_, Op_);
1643  os << " >> Error in KH == K*H : " << err << endl;
1644  }
1645  if (chk.checkMH && hasM_ && initialized_) {
1646  MagnitudeType err = Utils::errorEquality(*H_, *MH_, MOp_);
1647  os << " >> Error in MH == M*H : " << err << endl;
1648  }
1649 
1650  // R: this is not M-orthogonality, but standard Euclidean orthogonality
1651  if (chk.checkR && initialized_) {
1652  Teuchos::SerialDenseMatrix<int,ScalarType> xTx(blockSize_,blockSize_);
1653  MVT::MvTransMv(ONE,*X_,*R_,xTx);
1654  MagnitudeType err = xTx.normFrobenius();
1655  os << " >> Error in X^H R == 0 : " << err << endl;
1656  }
1657 
1658  // KK
1659  if (chk.checkKK && initialized_) {
1660  Teuchos::SerialDenseMatrix<int,ScalarType> SDMerr(curDim_,curDim_), lclKK(Teuchos::View,*KK_,curDim_,curDim_);
1661  for (int j=0; j<curDim_; ++j) {
1662  for (int i=0; i<curDim_; ++i) {
1663  SDMerr(i,j) = lclKK(i,j) - SCT::conjugate(lclKK(j,i));
1664  }
1665  }
1666  os << " >> Error in KK - KK^H == 0 : " << SDMerr.normFrobenius() << endl;
1667  }
1668 
1669  // Q
1670  if (chk.checkQ) {
1671  for (Array_size_type i=0; i<auxVecs_.size(); ++i) {
1672  MagnitudeType err = orthman_->orthonormError(*auxVecs_[i]);
1673  os << " >> Error in Q[" << i << "]^H M Q[" << i << "] == I : " << err << endl;
1674  for (Array_size_type j=i+1; j<auxVecs_.size(); ++j) {
1675  err = orthman_->orthogError(*auxVecs_[i],*auxVecs_[j]);
1676  os << " >> Error in Q[" << i << "]^H M Q[" << j << "] == 0 : " << err << endl;
1677  }
1678  }
1679  }
1680 
1681  os << endl;
1682 
1683  return os.str();
1684  }
1685 
1686 
1688  // Print the current status of the solver
1689  template <class ScalarType, class MV, class OP>
1690  void
1692  {
1693  using std::endl;
1694 
1695  os.setf(std::ios::scientific, std::ios::floatfield);
1696  os.precision(6);
1697  os <<endl;
1698  os <<"================================================================================" << endl;
1699  os << endl;
1700  os <<" BlockDavidson Solver Status" << endl;
1701  os << endl;
1702  os <<"The solver is "<<(initialized_ ? "initialized." : "not initialized.") << endl;
1703  os <<"The number of iterations performed is " <<iter_<<endl;
1704  os <<"The block size is " << blockSize_<<endl;
1705  os <<"The number of blocks is " << numBlocks_<<endl;
1706  os <<"The current basis size is " << curDim_<<endl;
1707  os <<"The number of auxiliary vectors is "<< numAuxVecs_ << endl;
1708  os <<"The number of operations Op*x is "<<count_ApplyOp_<<endl;
1709  os <<"The number of operations M*x is "<<count_ApplyM_<<endl;
1710  os <<"The number of operations Prec*x is "<<count_ApplyPrec_<<endl;
1711 
1712  os.setf(std::ios_base::right, std::ios_base::adjustfield);
1713 
1714  if (initialized_) {
1715  os << endl;
1716  os <<"CURRENT EIGENVALUE ESTIMATES "<<endl;
1717  os << std::setw(20) << "Eigenvalue"
1718  << std::setw(20) << "Residual(M)"
1719  << std::setw(20) << "Residual(2)"
1720  << endl;
1721  os <<"--------------------------------------------------------------------------------"<<endl;
1722  for (int i=0; i<blockSize_; ++i) {
1723  os << std::setw(20) << theta_[i];
1724  if (Rnorms_current_) os << std::setw(20) << Rnorms_[i];
1725  else os << std::setw(20) << "not current";
1726  if (R2norms_current_) os << std::setw(20) << R2norms_[i];
1727  else os << std::setw(20) << "not current";
1728  os << endl;
1729  }
1730  }
1731  os <<"================================================================================" << endl;
1732  os << endl;
1733  }
1734 
1735 } // End of namespace Anasazi
1736 
1737 #endif
1738 
1739 // End of file AnasaziBlockDavidson.hpp
virtual ~BlockDavidson()
Anasazi::BlockDavidson destructor.
bool isInitialized() const
Indicates whether the solver has been initialized or not.
Teuchos::RCP< const MV > H
The current preconditioned residual vectors.
int getNumIters() const
Get the current iteration count.
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getRitzRes2Norms()
Get the 2-norms of the residuals.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
void setSize(int blockSize, int numBlocks)
Set the blocksize and number of blocks to be used by the iterative solver in solving this eigenproble...
This class defines the interface required by an eigensolver and status test class to compute solution...
Teuchos::RCP< const MV > R
The current residual vectors.
Teuchos::RCP< const Teuchos::SerialDenseMatrix< int, ScalarType > > KK
The current projected K matrix.
Declaration of basic traits for the multivector type.
T & get(const std::string &name, T def_value)
BlockDavidsonOrthoFailure is thrown when the orthogonalization manager is unable to orthogonalize the...
void currentStatus(std::ostream &os)
This method requests that the solver print out its current status to the given output stream...
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues...
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void setStatusTest(Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test)
Set a new StatusTest for the solver.
Virtual base class which defines basic traits for the operator type.
Teuchos::RCP< const MV > KX
The image of the current eigenvectors under K.
BlockDavidsonInitFailure is thrown when the BlockDavidson solver is unable to generate an initial ite...
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > getStatusTest() const
Get the current StatusTest used by the solver.
An exception class parent to all Anasazi exceptions.
std::vector< Value< ScalarType > > getRitzValues()
Get the Ritz values for the previous iteration.
This class implements a Block Davidson iteration, a preconditioned iteration for solving linear Hermi...
Teuchos::Array< Teuchos::RCP< const MV > > getAuxVecs() const
Get the auxiliary vectors for the solver.
Pure virtual base class which describes the basic interface to the iterative eigensolver.
Anasazi&#39;s templated, static class providing utilities for the solvers.
Anasazi&#39;s templated virtual class for providing routines for orthogonalization and orthonormalization...
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace. For BlockDavidson, this always returns n...
Output managers remove the need for the eigensolver to know any information about the required output...
Teuchos::RCP< const MV > V
The basis for the Krylov space.
int putScalar(const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
Traits class which defines basic operations on multivectors.
Virtual base class which defines basic traits for the operator type.
int getBlockSize() const
Get the blocksize used by the iterative solver.
void setAuxVecs(const Teuchos::Array< Teuchos::RCP< const MV > > &auxvecs)
Set the auxiliary vectors for the solver.
const Eigenproblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the eigenvalue problem.
int curDim
The current dimension of the solver.
void initialize()
Initialize the solver with the initial vectors from the eigenproblem or random data.
void push_back(const value_type &x)
Structure to contain pointers to BlockDavidson state variables.
Teuchos::RCP< const MV > X
The current eigenvectors.
Teuchos::RCP< const MV > MX
The image of the current eigenvectors under M, or Teuchos::null if M was not specified.
Types and exceptions used within Anasazi solvers and interfaces.
std::vector< int > getRitzIndex()
Get the index used for extracting individual Ritz vectors from getRitzVectors().
void iterate()
This method performs BlockDavidson iterations until the status test indicates the need to stop or an ...
BlockDavidsonState< ScalarType, MV > getState() const
Get access to the current state of the eigensolver.
Anasazi&#39;s templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
Common interface of stopping criteria for Anasazi&#39;s solvers.
BlockDavidson(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< MatOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params)
BlockDavidson constructor with eigenproblem, solver utilities, and parameter list of solver options...
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getRes2Norms()
Get the current residual 2-norms, computing the norms if they are not up-to-date with the current res...
The Eigensolver is a templated virtual base class that defines the basic interface that any eigensolv...
void setBlockSize(int blockSize)
Set the blocksize.
Teuchos::RCP< const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > T
The current Ritz values. This vector is a copy of the internal data.
void resetNumIters()
Reset the iteration count.
Class which provides internal utilities for the Anasazi solvers.
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getResNorms()
Get the current residual norms, computing the norms if they are not up-to-date with the current resid...
Teuchos::RCP< const MV > getRitzVectors()
Get access to the current Ritz vectors.