Anasazi  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
AnasaziRTRBase.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_RTRBASE_HPP
15 #define ANASAZI_RTRBASE_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 
82 namespace Anasazi {
83 
85 
86 
91  template <class ScalarType, class MV>
92  struct RTRState {
107  RTRState() : X(Teuchos::null),AX(Teuchos::null),BX(Teuchos::null),
108  R(Teuchos::null),T(Teuchos::null),rho(0) {};
109  };
110 
112 
114 
115 
129  class RTRRitzFailure : public AnasaziError {public:
130  RTRRitzFailure(const std::string& what_arg) : AnasaziError(what_arg)
131  {}};
132 
141  class RTRInitFailure : public AnasaziError {public:
142  RTRInitFailure(const std::string& what_arg) : AnasaziError(what_arg)
143  {}};
144 
161  class RTROrthoFailure : public AnasaziError {public:
162  RTROrthoFailure(const std::string& what_arg) : AnasaziError(what_arg)
163  {}};
164 
165 
167 
168 
169  template <class ScalarType, class MV, class OP>
170  class RTRBase : public Eigensolver<ScalarType,MV,OP> {
171  public:
172 
174 
175 
183  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
186  Teuchos::ParameterList &params,
187  const std::string &solverLabel, bool skinnySolver
188  );
189 
191  virtual ~RTRBase() {};
192 
194 
196 
197 
219  virtual void iterate() = 0;
220 
245  void initialize(RTRState<ScalarType,MV>& newstate);
246 
250  void initialize();
251 
264  bool isInitialized() const;
265 
274 
276 
278 
279 
281  int getNumIters() const;
282 
284  void resetNumIters();
285 
294 
300  std::vector<Value<ScalarType> > getRitzValues();
301 
309  std::vector<int> getRitzIndex();
310 
316  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getResNorms();
317 
318 
324  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRes2Norms();
325 
326 
331  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> getRitzRes2Norms();
332 
333 
342  int getCurSubspaceDim() const;
343 
346  int getMaxSubspaceDim() const;
347 
349 
351 
352 
355 
358 
361 
362 
371  void setBlockSize(int blockSize);
372 
373 
375  int getBlockSize() const;
376 
377 
398  void setAuxVecs(const Teuchos::Array<Teuchos::RCP<const MV> > &auxvecs);
399 
402 
404 
406 
407 
409  virtual void currentStatus(std::ostream &os);
410 
412 
413  protected:
414  //
415  // Convenience typedefs
416  //
417  typedef SolverUtils<ScalarType,MV,OP> Utils;
418  typedef MultiVecTraits<ScalarType,MV> MVT;
421  typedef typename SCT::magnitudeType MagnitudeType;
423  const MagnitudeType ONE;
424  const MagnitudeType ZERO;
425  const MagnitudeType NANVAL;
426  typedef typename std::vector<MagnitudeType>::iterator vecMTiter;
427  typedef typename std::vector<ScalarType>::iterator vecSTiter;
428  //
429  // Internal structs
430  //
431  struct CheckList {
432  bool checkX, checkAX, checkBX;
433  bool checkEta, checkAEta, checkBEta;
434  bool checkR, checkQ, checkBR;
435  bool checkZ, checkPBX;
436  CheckList() : checkX(false),checkAX(false),checkBX(false),
437  checkEta(false),checkAEta(false),checkBEta(false),
438  checkR(false),checkQ(false),checkBR(false),
439  checkZ(false), checkPBX(false) {};
440  };
441  //
442  // Internal methods
443  //
444  std::string accuracyCheck(const CheckList &chk, const std::string &where) const;
445  // solves the model minimization
446  virtual void solveTRSubproblem() = 0;
447  // Riemannian metric
448  typename Teuchos::ScalarTraits<ScalarType>::magnitudeType ginner(const MV &xi) const;
449  typename Teuchos::ScalarTraits<ScalarType>::magnitudeType ginner(const MV &xi, const MV &zeta) const;
450  void ginnersep(const MV &xi, std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType > &d) const;
451  void ginnersep(const MV &xi, const MV &zeta, std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType > &d) const;
452  //
453  // Classes input through constructor that define the eigenproblem to be solved.
454  //
460  //
461  // Information obtained from the eigenproblem
462  //
466  bool hasBOp_, hasPrec_, olsenPrec_;
467  //
468  // Internal timers
469  //
470  Teuchos::RCP<Teuchos::Time> timerAOp_, timerBOp_, timerPrec_,
471  timerSort_,
472  timerLocalProj_, timerDS_,
473  timerLocalUpdate_, timerCompRes_,
474  timerOrtho_, timerInit_;
475  //
476  // Counters
477  //
478  // Number of operator applications
479  int counterAOp_, counterBOp_, counterPrec_;
480 
481  //
482  // Algorithmic parameters.
483  //
484  // blockSize_ is the solver block size
485  int blockSize_;
486  //
487  // Current solver state
488  //
489  // initialized_ specifies that the basis vectors have been initialized and the iterate() routine
490  // is capable of running; _initialize is controlled by the initialize() member method
491  // For the implications of the state of initialized_, please see documentation for initialize()
492  bool initialized_;
493  //
494  // nevLocal_ reflects how much of the current basis is valid (0 <= nevLocal_ <= blockSize_)
495  // this tells us how many of the values in theta_ are valid Ritz values
496  int nevLocal_;
497  //
498  // are we implementing a skinny solver? (SkinnyIRTR)
499  bool isSkinny_;
500  //
501  // are we computing leftmost or rightmost eigenvalue?
502  bool leftMost_;
503  //
504  // State Multivecs
505  //
506  // if we are implementing a skinny solver (SkinnyIRTR),
507  // then some of these will never be allocated
508  //
509  // In order to handle auxiliary vectors, we need to handle the projector
510  // P_{[BQ BX],[BQ BX]}
511  // Using an orthomanager with B-inner product, this requires calling with multivectors
512  // [BQ,BX] and [Q,X].
513  // These multivectors must be combined because <[BQ,BX],[Q,X]>_B != I
514  // Hence, we will create two multivectors V and BV, which store
515  // V = [Q,X]
516  // BV = [BQ,BX]
517  //
518  // In the context of preconditioning, we may need to apply the projector
519  // P_{prec*[BQ,BX],[BQ,BX]}
520  // Because [BQ,BX] do not change during the supproblem solver, we will apply
521  // the preconditioner to [BQ,BX] only once. This result is stored in PBV.
522  //
523  // X,BX are views into V,BV
524  // We don't need views for Q,BQ
525  // Inside the subproblem solver, X,BX are constant, so we can allow these
526  // views to exist alongside the full view of V,BV without violating
527  // view semantics.
528  //
529  // Skinny solver allocates 6/7/8 multivectors:
530  // V_, BV_ (if hasB)
531  // PBV_ (if hasPrec and olsenPrec)
532  // R_, Z_ (regardless of hasPrec)
533  // eta_, delta_, Hdelta_
534  //
535  // Hefty solver allocates 10/11/12/13 multivectors:
536  // V_, AX_, BV_ (if hasB)
537  // PBV_ (if hasPrec and olsenPrec)
538  // R_, Z_ (if hasPrec)
539  // eta_, Aeta_, Beta_
540  // delta_, Adelta_, Bdelta_, Hdelta_
541  //
542  Teuchos::RCP<MV> V_, BV_, PBV_, // V = [Q,X]; B*V; Prec*B*V
543  AX_, R_, // A*X_; residual, gradient, and residual of model minimization
544  eta_, Aeta_, Beta_, // update vector from model minimization
545  delta_, Adelta_, Bdelta_, Hdelta_, // search direction in model minimization
546  Z_; // preconditioned residual
547  Teuchos::RCP<const MV> X_, BX_;
548  //
549  // auxiliary vectors
551  int numAuxVecs_;
552  //
553  // Number of iterations that have been performed.
554  int iter_;
555  //
556  // Current eigenvalues, residual norms
557  std::vector<MagnitudeType> theta_, Rnorms_, R2norms_, ritz2norms_;
558  //
559  // are the residual norms current with the residual?
560  bool Rnorms_current_, R2norms_current_;
561  //
562  // parameters solver and inner solver
563  MagnitudeType conv_kappa_, conv_theta_;
564  MagnitudeType rho_;
565  //
566  // current objective function value
567  MagnitudeType fx_;
568  };
569 
570 
571 
572 
574  // Constructor
575  template <class ScalarType, class MV, class OP>
579  const Teuchos::RCP<OutputManager<ScalarType> > &printer,
582  Teuchos::ParameterList &params,
583  const std::string &solverLabel, bool skinnySolver
584  ) :
585  ONE(Teuchos::ScalarTraits<MagnitudeType>::one()),
586  ZERO(Teuchos::ScalarTraits<MagnitudeType>::zero()),
587  NANVAL(Teuchos::ScalarTraits<MagnitudeType>::nan()),
588  // problem, tools
589  problem_(problem),
590  sm_(sorter),
591  om_(printer),
592  tester_(tester),
593  orthman_(ortho),
594 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
595  timerAOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Operation A*x")),
596  timerBOp_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Operation B*x")),
597  timerPrec_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Operation Prec*x")),
598  timerSort_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Sorting eigenvalues")),
599  timerLocalProj_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Local projection")),
600  timerDS_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Direct solve")),
601  timerLocalUpdate_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Local update")),
602  timerCompRes_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Computing residuals")),
603  timerOrtho_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Orthogonalization")),
604  timerInit_(Teuchos::TimeMonitor::getNewTimer("Anasazi: "+solverLabel+"::Initialization")),
605 #endif
606  counterAOp_(0),
607  counterBOp_(0),
608  counterPrec_(0),
609  // internal data
610  blockSize_(0),
611  initialized_(false),
612  nevLocal_(0),
613  isSkinny_(skinnySolver),
614  leftMost_(true),
615  auxVecs_( Teuchos::Array<Teuchos::RCP<const MV> >(0) ),
616  numAuxVecs_(0),
617  iter_(0),
618  Rnorms_current_(false),
619  R2norms_current_(false),
620  conv_kappa_(.1),
621  conv_theta_(1),
622  rho_( MAT::nan() ),
623  fx_( MAT::nan() )
624  {
625  TEUCHOS_TEST_FOR_EXCEPTION(problem_ == Teuchos::null,std::invalid_argument,
626  "Anasazi::RTRBase::constructor: user passed null problem pointer.");
627  TEUCHOS_TEST_FOR_EXCEPTION(sm_ == Teuchos::null,std::invalid_argument,
628  "Anasazi::RTRBase::constructor: user passed null sort manager pointer.");
629  TEUCHOS_TEST_FOR_EXCEPTION(om_ == Teuchos::null,std::invalid_argument,
630  "Anasazi::RTRBase::constructor: user passed null output manager pointer.");
631  TEUCHOS_TEST_FOR_EXCEPTION(tester_ == Teuchos::null,std::invalid_argument,
632  "Anasazi::RTRBase::constructor: user passed null status test pointer.");
633  TEUCHOS_TEST_FOR_EXCEPTION(orthman_ == Teuchos::null,std::invalid_argument,
634  "Anasazi::RTRBase::constructor: user passed null orthogonalization manager pointer.");
635  TEUCHOS_TEST_FOR_EXCEPTION(problem_->isProblemSet() == false, std::invalid_argument,
636  "Anasazi::RTRBase::constructor: problem is not set.");
637  TEUCHOS_TEST_FOR_EXCEPTION(problem_->isHermitian() == false, std::invalid_argument,
638  "Anasazi::RTRBase::constructor: problem is not Hermitian.");
639 
640  // get the problem operators
641  AOp_ = problem_->getOperator();
642  TEUCHOS_TEST_FOR_EXCEPTION(AOp_ == Teuchos::null, std::invalid_argument,
643  "Anasazi::RTRBase::constructor: problem provides no A matrix.");
644  BOp_ = problem_->getM();
645  Prec_ = problem_->getPrec();
646  hasBOp_ = (BOp_ != Teuchos::null);
647  hasPrec_ = (Prec_ != Teuchos::null);
648  olsenPrec_ = params.get<bool>("Olsen Prec", true);
649 
650  TEUCHOS_TEST_FOR_EXCEPTION(orthman_->getOp() != BOp_,std::invalid_argument,
651  "Anasazi::RTRBase::constructor: orthogonalization manager must use mass matrix for inner product.");
652 
653  // set the block size and allocate data
654  int bs = params.get("Block Size", problem_->getNEV());
655  setBlockSize(bs);
656 
657  // leftmost or rightmost?
658  leftMost_ = params.get("Leftmost",leftMost_);
659 
660  conv_kappa_ = params.get("Kappa Convergence",conv_kappa_);
661  TEUCHOS_TEST_FOR_EXCEPTION(conv_kappa_ <= 0 || conv_kappa_ >= 1,std::invalid_argument,
662  "Anasazi::RTRBase::constructor: kappa must be in (0,1).");
663  conv_theta_ = params.get("Theta Convergence",conv_theta_);
664  TEUCHOS_TEST_FOR_EXCEPTION(conv_theta_ <= 0,std::invalid_argument,
665  "Anasazi::RTRBase::constructor: theta must be strictly postitive.");
666  }
667 
668 
670  // Set the block size and make necessary adjustments.
671  template <class ScalarType, class MV, class OP>
673  {
674  // time spent here counts towards timerInit_
675 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
676  Teuchos::TimeMonitor lcltimer( *timerInit_ );
677 #endif
678 
679  // This routine only allocates space; it doesn't not perform any computation
680  // if solver is initialized and size is to be decreased, take the first blockSize vectors of all to preserve state
681  // otherwise, shrink/grow/allocate space and set solver to unitialized
682 
684  // grab some Multivector to Clone
685  // in practice, getInitVec() should always provide this, but it is possible to use a
686  // Eigenproblem with nothing in getInitVec() by manually initializing with initialize();
687  // in case of that strange scenario, we will try to Clone from R_
688  // we like R_ for this, because it has minimal size (blockSize_), as opposed to V_ (blockSize_+numAuxVecs_)
689  if (blockSize_ > 0) {
690  tmp = R_;
691  }
692  else {
693  tmp = problem_->getInitVec();
694  TEUCHOS_TEST_FOR_EXCEPTION(tmp == Teuchos::null,std::logic_error,
695  "Anasazi::RTRBase::setBlockSize(): Eigenproblem did not specify initial vectors to clone from");
696  }
697 
698  TEUCHOS_TEST_FOR_EXCEPTION(blockSize <= 0 || blockSize > MVT::GetGlobalLength(*tmp), std::invalid_argument,
699  "Anasazi::RTRBase::setBlockSize was passed a non-positive block size");
700 
701  // last chance to quit before causing side-effects
702  if (blockSize == blockSize_) {
703  // do nothing
704  return;
705  }
706 
707  // clear views
708  X_ = Teuchos::null;
709  BX_ = Teuchos::null;
710 
711  // regardless of whether we preserve any data, any content in V, BV and PBV corresponding to the
712  // auxiliary vectors must be retained
713  // go ahead and do these first
714  //
715  // two cases here:
716  // * we are growing (possibly, from empty)
717  // any aux data must be copied over, nothing from X need copying
718  // * we are shrinking
719  // any aux data must be copied over, go ahead and copy X material (initialized or not)
720  //
721  if (blockSize > blockSize_)
722  {
723  // GROWING
724  // get a pointer for Q-related material, and an index for extracting/setting it
726  std::vector<int> indQ(numAuxVecs_);
727  for (int i=0; i<numAuxVecs_; i++) indQ[i] = i;
728  // if numAuxVecs_ > 0, then necessarily blockSize_ > 0 (we have already been allocated once)
729  TEUCHOS_TEST_FOR_EXCEPTION(numAuxVecs_ > 0 && blockSize_ == 0, std::logic_error,
730  "Anasazi::RTRBase::setSize(): logic error. Please contact Anasazi team.");
731  // V
732  if (numAuxVecs_ > 0) Q = MVT::CloneView(*V_,indQ);
733  V_ = MVT::Clone(*tmp,numAuxVecs_ + blockSize);
734  if (numAuxVecs_ > 0) MVT::SetBlock(*Q,indQ,*V_);
735  // BV
736  if (hasBOp_) {
737  if (numAuxVecs_ > 0) Q = MVT::CloneView(*BV_,indQ);
738  BV_ = MVT::Clone(*tmp,numAuxVecs_ + blockSize);
739  if (numAuxVecs_ > 0) MVT::SetBlock(*Q,indQ,*BV_);
740  }
741  else {
742  BV_ = V_;
743  }
744  // PBV
745  if (hasPrec_ && olsenPrec_) {
746  if (numAuxVecs_ > 0) Q = MVT::CloneView(*PBV_,indQ);
747  PBV_ = MVT::Clone(*tmp,numAuxVecs_ + blockSize);
748  if (numAuxVecs_ > 0) MVT::SetBlock(*Q,indQ,*PBV_);
749  }
750  }
751  else
752  {
753  // SHRINKING
754  std::vector<int> indV(numAuxVecs_+blockSize);
755  for (int i=0; i<numAuxVecs_+blockSize; i++) indV[i] = i;
756  // V
757  V_ = MVT::CloneCopy(*V_,indV);
758  // BV
759  if (hasBOp_) {
760  BV_ = MVT::CloneCopy(*BV_,indV);
761  }
762  else {
763  BV_ = V_;
764  }
765  // PBV
766  if (hasPrec_ && olsenPrec_) {
767  PBV_ = MVT::CloneCopy(*PBV_,indV);
768  }
769  }
770 
771  if (blockSize < blockSize_) {
772  // shrink vectors
773  blockSize_ = blockSize;
774 
775  theta_.resize(blockSize_);
776  ritz2norms_.resize(blockSize_);
777  Rnorms_.resize(blockSize_);
778  R2norms_.resize(blockSize_);
779 
780  if (initialized_) {
781  // shrink multivectors with copy
782  std::vector<int> ind(blockSize_);
783  for (int i=0; i<blockSize_; i++) ind[i] = i;
784 
785  // Z can be deleted, no important info there
786  Z_ = Teuchos::null;
787 
788  // we will not use tmp for cloning, clear it and free some space
789  tmp = Teuchos::null;
790 
791  R_ = MVT::CloneCopy(*R_ ,ind);
792  eta_ = MVT::CloneCopy(*eta_ ,ind);
793  delta_ = MVT::CloneCopy(*delta_ ,ind);
794  Hdelta_ = MVT::CloneCopy(*Hdelta_,ind);
795  if (!isSkinny_) {
796  AX_ = MVT::CloneCopy(*AX_ ,ind);
797  Aeta_ = MVT::CloneCopy(*Aeta_ ,ind);
798  Adelta_ = MVT::CloneCopy(*Adelta_,ind);
799  }
800  else {
801  AX_ = Teuchos::null;
802  Aeta_ = Teuchos::null;
803  Adelta_ = Teuchos::null;
804  }
805 
806  if (hasBOp_) {
807  if (!isSkinny_) {
808  Beta_ = MVT::CloneCopy(*Beta_,ind);
809  Bdelta_ = MVT::CloneCopy(*Bdelta_,ind);
810  }
811  else {
812  Beta_ = Teuchos::null;
813  Bdelta_ = Teuchos::null;
814  }
815  }
816  else {
817  Beta_ = eta_;
818  Bdelta_ = delta_;
819  }
820 
821  // we need Z if we have a preconditioner
822  // we also use Z for temp storage in the skinny solvers
823  if (hasPrec_ || isSkinny_) {
824  Z_ = MVT::Clone(*V_,blockSize_);
825  }
826  else {
827  Z_ = R_;
828  }
829 
830  }
831  else {
832  // release previous multivectors
833  // shrink multivectors without copying
834  AX_ = Teuchos::null;
835  R_ = Teuchos::null;
836  eta_ = Teuchos::null;
837  Aeta_ = Teuchos::null;
838  delta_ = Teuchos::null;
839  Adelta_ = Teuchos::null;
840  Hdelta_ = Teuchos::null;
841  Beta_ = Teuchos::null;
842  Bdelta_ = Teuchos::null;
843  Z_ = Teuchos::null;
844 
845  R_ = MVT::Clone(*tmp,blockSize_);
846  eta_ = MVT::Clone(*tmp,blockSize_);
847  delta_ = MVT::Clone(*tmp,blockSize_);
848  Hdelta_ = MVT::Clone(*tmp,blockSize_);
849  if (!isSkinny_) {
850  AX_ = MVT::Clone(*tmp,blockSize_);
851  Aeta_ = MVT::Clone(*tmp,blockSize_);
852  Adelta_ = MVT::Clone(*tmp,blockSize_);
853  }
854 
855  if (hasBOp_) {
856  if (!isSkinny_) {
857  Beta_ = MVT::Clone(*tmp,blockSize_);
858  Bdelta_ = MVT::Clone(*tmp,blockSize_);
859  }
860  }
861  else {
862  Beta_ = eta_;
863  Bdelta_ = delta_;
864  }
865 
866  // we need Z if we have a preconditioner
867  // we also use Z for temp storage in the skinny solvers
868  if (hasPrec_ || isSkinny_) {
869  Z_ = MVT::Clone(*tmp,blockSize_);
870  }
871  else {
872  Z_ = R_;
873  }
874  } // if initialized_
875  } // if blockSize is shrinking
876  else { // if blockSize > blockSize_
877  // this is also the scenario for our initial call to setBlockSize(), in the constructor
878  initialized_ = false;
879 
880  // grow/allocate vectors
881  theta_.resize(blockSize,NANVAL);
882  ritz2norms_.resize(blockSize,NANVAL);
883  Rnorms_.resize(blockSize,NANVAL);
884  R2norms_.resize(blockSize,NANVAL);
885 
886  // deallocate old multivectors
887  AX_ = Teuchos::null;
888  R_ = Teuchos::null;
889  eta_ = Teuchos::null;
890  Aeta_ = Teuchos::null;
891  delta_ = Teuchos::null;
892  Adelta_ = Teuchos::null;
893  Hdelta_ = Teuchos::null;
894  Beta_ = Teuchos::null;
895  Bdelta_ = Teuchos::null;
896  Z_ = Teuchos::null;
897 
898  // clone multivectors off of tmp
899  R_ = MVT::Clone(*tmp,blockSize);
900  eta_ = MVT::Clone(*tmp,blockSize);
901  delta_ = MVT::Clone(*tmp,blockSize);
902  Hdelta_ = MVT::Clone(*tmp,blockSize);
903  if (!isSkinny_) {
904  AX_ = MVT::Clone(*tmp,blockSize);
905  Aeta_ = MVT::Clone(*tmp,blockSize);
906  Adelta_ = MVT::Clone(*tmp,blockSize);
907  }
908 
909  if (hasBOp_) {
910  if (!isSkinny_) {
911  Beta_ = MVT::Clone(*tmp,blockSize);
912  Bdelta_ = MVT::Clone(*tmp,blockSize);
913  }
914  }
915  else {
916  Beta_ = eta_;
917  Bdelta_ = delta_;
918  }
919  if (hasPrec_ || isSkinny_) {
920  Z_ = MVT::Clone(*tmp,blockSize);
921  }
922  else {
923  Z_ = R_;
924  }
925  blockSize_ = blockSize;
926  }
927 
928  // get view of X from V, BX from BV
929  // these are located after the first numAuxVecs columns
930  {
931  std::vector<int> indX(blockSize_);
932  for (int i=0; i<blockSize_; i++) indX[i] = numAuxVecs_+i;
933  X_ = MVT::CloneView(*V_,indX);
934  if (hasBOp_) {
935  BX_ = MVT::CloneView(*BV_,indX);
936  }
937  else {
938  BX_ = X_;
939  }
940  }
941  }
942 
943 
945  // Set a new StatusTest for the solver.
946  template <class ScalarType, class MV, class OP>
948  TEUCHOS_TEST_FOR_EXCEPTION(test == Teuchos::null,std::invalid_argument,
949  "Anasazi::RTRBase::setStatusTest() was passed a null StatusTest.");
950  tester_ = test;
951  }
952 
953 
955  // Get the current StatusTest used by the solver.
956  template <class ScalarType, class MV, class OP>
958  return tester_;
959  }
960 
961 
963  // Set the auxiliary vectors
964  template <class ScalarType, class MV, class OP>
966  typedef typename Teuchos::Array<Teuchos::RCP<const MV> >::const_iterator tarcpmv;
967 
968  // set new auxiliary vectors
969  auxVecs_.resize(0);
970  auxVecs_.reserve(auxvecs.size());
971 
972  numAuxVecs_ = 0;
973  for (tarcpmv v=auxvecs.begin(); v != auxvecs.end(); ++v) {
974  numAuxVecs_ += MVT::GetNumberVecs(**v);
975  }
976 
977  // If the solver has been initialized, X is not necessarily orthogonal to new auxiliary vectors
978  if (numAuxVecs_ > 0 && initialized_) {
979  initialized_ = false;
980  }
981 
982  // clear X,BX views
983  X_ = Teuchos::null;
984  BX_ = Teuchos::null;
985  // deallocate, we'll clone off R_ below
986  V_ = Teuchos::null;
987  BV_ = Teuchos::null;
988  PBV_ = Teuchos::null;
989 
990  // put auxvecs into V, update BV and PBV if necessary
991  if (numAuxVecs_ > 0) {
992  V_ = MVT::Clone(*R_,numAuxVecs_ + blockSize_);
993  int numsofar = 0;
994  for (tarcpmv v=auxvecs.begin(); v != auxvecs.end(); ++v) {
995  std::vector<int> ind(MVT::GetNumberVecs(**v));
996  for (size_t j=0; j<ind.size(); j++) ind[j] = numsofar++;
997  MVT::SetBlock(**v,ind,*V_);
998  auxVecs_.push_back(MVT::CloneView(*Teuchos::rcp_static_cast<const MV>(V_),ind));
999  }
1000  TEUCHOS_TEST_FOR_EXCEPTION(numsofar != numAuxVecs_, std::logic_error,
1001  "Anasazi::RTRBase::setAuxVecs(): Logic error. Please contact Anasazi team.");
1002  // compute B*V, Prec*B*V
1003  if (hasBOp_) {
1004  BV_ = MVT::Clone(*R_,numAuxVecs_ + blockSize_);
1005  OPT::Apply(*BOp_,*V_,*BV_);
1006  }
1007  else {
1008  BV_ = V_;
1009  }
1010  if (hasPrec_ && olsenPrec_) {
1011  PBV_ = MVT::Clone(*R_,numAuxVecs_ + blockSize_);
1012  OPT::Apply(*Prec_,*BV_,*V_);
1013  }
1014  }
1015  //
1016 
1017  if (om_->isVerbosity( Debug ) ) {
1018  // Check almost everything here
1019  CheckList chk;
1020  chk.checkQ = true;
1021  om_->print( Debug, accuracyCheck(chk, "in setAuxVecs()") );
1022  }
1023  }
1024 
1025 
1027  /* Initialize the state of the solver
1028  *
1029  * POST-CONDITIONS:
1030  *
1031  * initialized_ == true
1032  * X is orthonormal, orthogonal to auxVecs_
1033  * AX = A*X if not skinnny
1034  * BX = B*X if hasBOp_
1035  * theta_ contains Ritz values of X
1036  * R = AX - BX*diag(theta_)
1037  */
1038  template <class ScalarType, class MV, class OP>
1040  {
1041  // NOTE: memory has been allocated by setBlockSize(). Use SetBlock below; do not Clone
1042  // NOTE: Overall time spent in this routine is counted to timerInit_; portions will also be counted towards other primitives
1043 
1044  // clear const views to X,BX in V,BV
1045  // work with temporary non-const views
1046  X_ = Teuchos::null;
1047  BX_ = Teuchos::null;
1048  Teuchos::RCP<MV> X, BX;
1049  {
1050  std::vector<int> ind(blockSize_);
1051  for (int i=0; i<blockSize_; ++i) ind[i] = numAuxVecs_+i;
1052  X = MVT::CloneViewNonConst(*V_,ind);
1053  if (hasBOp_) {
1054  BX = MVT::CloneViewNonConst(*BV_,ind);
1055  }
1056  else {
1057  BX = X;
1058  }
1059  }
1060 
1061 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1062  Teuchos::TimeMonitor inittimer( *timerInit_ );
1063 #endif
1064 
1065  std::vector<int> bsind(blockSize_);
1066  for (int i=0; i<blockSize_; i++) bsind[i] = i;
1067 
1068  // in RTR, X (the subspace iterate) is primary
1069  // the order of dependence follows like so.
1070  // --init-> X
1071  // --op apply-> AX,BX
1072  // --ritz analysis-> theta
1073  //
1074  // if the user specifies all data for a level, we will accept it.
1075  // otherwise, we will generate the whole level, and all subsequent levels.
1076  //
1077  // the data members are ordered based on dependence, and the levels are
1078  // partitioned according to the amount of work required to produce the
1079  // items in a level.
1080  //
1081  // inconsitent multivectors widths and lengths will not be tolerated, and
1082  // will be treated with exceptions.
1083 
1084  // set up X, AX, BX: get them from "state" if user specified them
1085  if (newstate.X != Teuchos::null) {
1086  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.X) != MVT::GetGlobalLength(*X),
1087  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): vector length of newstate.X not correct." );
1088  // newstate.X must have blockSize_ vectors; any more will be ignored
1089  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.X) < blockSize_,
1090  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): newstate.X must have at least block size vectors.");
1091 
1092  // put data in X
1093  MVT::SetBlock(*newstate.X,bsind,*X);
1094 
1095  // put data in AX
1096  // if we are implementing a skinny solver, then we don't have memory allocated for AX
1097  // in this case, point AX at Z (skinny solvers allocate Z) and use it for temporary storage
1098  // we will clear the pointer later
1099  if (isSkinny_) {
1100  AX_ = Z_;
1101  }
1102  if (newstate.AX != Teuchos::null) {
1103  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.AX) != MVT::GetGlobalLength(*X),
1104  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): vector length of newstate.AX not correct." );
1105  // newstate.AX must have blockSize_ vectors; any more will be ignored
1106  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.AX) < blockSize_,
1107  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): newstate.AX must have at least block size vectors.");
1108  MVT::SetBlock(*newstate.AX,bsind,*AX_);
1109  }
1110  else {
1111  {
1112 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1113  Teuchos::TimeMonitor lcltimer( *timerAOp_ );
1114 #endif
1115  OPT::Apply(*AOp_,*X,*AX_);
1116  counterAOp_ += blockSize_;
1117  }
1118  // we generated AX; we will generate R as well
1119  newstate.R = Teuchos::null;
1120  }
1121 
1122  // put data in BX
1123  // skinny solvers always allocate BX if hasB, so this is unconditionally appropriate
1124  if (hasBOp_) {
1125  if (newstate.BX != Teuchos::null) {
1126  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.BX) != MVT::GetGlobalLength(*X),
1127  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): vector length of newstate.BX not correct." );
1128  // newstate.BX must have blockSize_ vectors; any more will be ignored
1129  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.BX) < blockSize_,
1130  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): newstate.BX must have at least block size vectors.");
1131  MVT::SetBlock(*newstate.BX,bsind,*BX);
1132  }
1133  else {
1134  {
1135 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1136  Teuchos::TimeMonitor lcltimer( *timerBOp_ );
1137 #endif
1138  OPT::Apply(*BOp_,*X,*BX);
1139  counterBOp_ += blockSize_;
1140  }
1141  // we generated BX; we will generate R as well
1142  newstate.R = Teuchos::null;
1143  }
1144  }
1145  else {
1146  // the assignment BX_==X_ would be redundant; take advantage of this opportunity to debug a little
1147  TEUCHOS_TEST_FOR_EXCEPTION(BX != X, std::logic_error, "Anasazi::RTRBase::initialize(): solver invariant not satisfied (BX==X).");
1148  }
1149 
1150  }
1151  else {
1152  // user did not specify X
1153 
1154  // clear state so we won't use any data from it below
1155  newstate.R = Teuchos::null;
1156  newstate.T = Teuchos::null;
1157 
1158  // generate something and projectAndNormalize
1159  Teuchos::RCP<const MV> ivec = problem_->getInitVec();
1160  TEUCHOS_TEST_FOR_EXCEPTION(ivec == Teuchos::null,std::logic_error,
1161  "Anasazi::RTRBase::initialize(): Eigenproblem did not specify initial vectors to clone from.");
1162 
1163  int initSize = MVT::GetNumberVecs(*ivec);
1164  if (initSize > blockSize_) {
1165  // we need only the first blockSize_ vectors from ivec; get a view of them
1166  initSize = blockSize_;
1167  std::vector<int> ind(blockSize_);
1168  for (int i=0; i<blockSize_; i++) ind[i] = i;
1169  ivec = MVT::CloneView(*ivec,ind);
1170  }
1171 
1172  // assign ivec to first part of X
1173  if (initSize > 0) {
1174  std::vector<int> ind(initSize);
1175  for (int i=0; i<initSize; i++) ind[i] = i;
1176  MVT::SetBlock(*ivec,ind,*X);
1177  }
1178  // fill the rest of X with random
1179  if (blockSize_ > initSize) {
1180  std::vector<int> ind(blockSize_ - initSize);
1181  for (int i=0; i<blockSize_ - initSize; i++) ind[i] = initSize + i;
1182  Teuchos::RCP<MV> rX = MVT::CloneViewNonConst(*X,ind);
1183  MVT::MvRandom(*rX);
1184  rX = Teuchos::null;
1185  }
1186 
1187  // put data in BX
1188  if (hasBOp_) {
1189 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1190  Teuchos::TimeMonitor lcltimer( *timerBOp_ );
1191 #endif
1192  OPT::Apply(*BOp_,*X,*BX);
1193  counterBOp_ += blockSize_;
1194  }
1195  else {
1196  // the assignment BX==X would be redundant; take advantage of this opportunity to debug a little
1197  TEUCHOS_TEST_FOR_EXCEPTION(BX != X, std::logic_error, "Anasazi::RTRBase::initialize(): solver invariant not satisfied (BX==X).");
1198  }
1199 
1200  // remove auxVecs from X and normalize it
1201  if (numAuxVecs_ > 0) {
1202 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1203  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
1204 #endif
1206  int rank = orthman_->projectAndNormalizeMat(*X,auxVecs_,dummyC,Teuchos::null,BX);
1207  TEUCHOS_TEST_FOR_EXCEPTION(rank != blockSize_, RTRInitFailure,
1208  "Anasazi::RTRBase::initialize(): Couldn't generate initial basis of full rank.");
1209  }
1210  else {
1211 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1212  Teuchos::TimeMonitor lcltimer( *timerOrtho_ );
1213 #endif
1214  int rank = orthman_->normalizeMat(*X,Teuchos::null,BX);
1215  TEUCHOS_TEST_FOR_EXCEPTION(rank != blockSize_, RTRInitFailure,
1216  "Anasazi::RTRBase::initialize(): Couldn't generate initial basis of full rank.");
1217  }
1218 
1219  // put data in AX
1220  if (isSkinny_) {
1221  AX_ = Z_;
1222  }
1223  {
1224 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1225  Teuchos::TimeMonitor lcltimer( *timerAOp_ );
1226 #endif
1227  OPT::Apply(*AOp_,*X,*AX_);
1228  counterAOp_ += blockSize_;
1229  }
1230 
1231  } // end if (newstate.X != Teuchos::null)
1232 
1233 
1234  // set up Ritz values
1235  if (newstate.T != Teuchos::null) {
1236  TEUCHOS_TEST_FOR_EXCEPTION( (signed int)(newstate.T->size()) < blockSize_,
1237  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): newstate.T must contain at least block size Ritz values.");
1238  for (int i=0; i<blockSize_; i++) {
1239  theta_[i] = (*newstate.T)[i];
1240  }
1241  }
1242  else {
1243  // get ritz vecs/vals
1244  Teuchos::SerialDenseMatrix<int,ScalarType> AA(blockSize_,blockSize_),
1245  BB(blockSize_,blockSize_),
1246  S(blockSize_,blockSize_);
1247  // project A
1248  {
1249 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1250  Teuchos::TimeMonitor lcltimer( *timerLocalProj_ );
1251 #endif
1252  MVT::MvTransMv(ONE,*X,*AX_,AA);
1253  if (hasBOp_) {
1254  MVT::MvTransMv(ONE,*X,*BX,BB);
1255  }
1256  }
1257  nevLocal_ = blockSize_;
1258 
1259  // solve the projected problem
1260  // project B
1261  int ret;
1262  if (hasBOp_) {
1263 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1264  Teuchos::TimeMonitor lcltimer( *timerDS_ );
1265 #endif
1266  ret = Utils::directSolver(blockSize_,AA,Teuchos::rcpFromRef(BB),S,theta_,nevLocal_,1);
1267  }
1268  else {
1269 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1270  Teuchos::TimeMonitor lcltimer( *timerDS_ );
1271 #endif
1272  ret = Utils::directSolver(blockSize_,AA,Teuchos::null,S,theta_,nevLocal_,10);
1273  }
1275  "Anasazi::RTRBase::initialize(): failure solving projected eigenproblem after retraction. LAPACK returns " << ret);
1276  TEUCHOS_TEST_FOR_EXCEPTION(nevLocal_ != blockSize_,RTRInitFailure,"Anasazi::RTRBase::initialize(): retracted iterate failed in Ritz analysis.");
1277 
1278  // We only have blockSize_ ritz pairs, ergo we do not need to select.
1279  // However, we still require them to be ordered correctly
1280  {
1281 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1282  Teuchos::TimeMonitor lcltimer( *timerSort_ );
1283 #endif
1284 
1285  std::vector<int> order(blockSize_);
1286  //
1287  // sort the first blockSize_ values in theta_
1288  sm_->sort(theta_, Teuchos::rcpFromRef(order), blockSize_); // don't catch exception
1289  //
1290  // apply the same ordering to the primitive ritz vectors
1291  Utils::permuteVectors(order,S);
1292  }
1293 
1294  // compute ritz residual norms
1295  {
1297  Teuchos::SerialDenseMatrix<int,ScalarType> RR(blockSize_,blockSize_);
1298  // RR = AA*S - BB*S*diag(theta)
1299  int info;
1300  if (hasBOp_) {
1301  info = RR.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,ONE,BB,S,ZERO);
1302  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::logic_error, "Anasazi::RTRBase::initialize(): Logic error calling SerialDenseMatrix::multiply.");
1303  }
1304  else {
1305  RR.assign(S);
1306  }
1307  for (int i=0; i<blockSize_; i++) {
1308  blas.SCAL(blockSize_,theta_[i],RR[i],1);
1309  }
1310  info = RR.multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,ONE,AA,S,-ONE);
1311  TEUCHOS_TEST_FOR_EXCEPTION(info != 0, std::logic_error, "Anasazi::RTRBase::initialize(): Logic error calling SerialDenseMatrix::multiply.");
1312  for (int i=0; i<blockSize_; i++) {
1313  ritz2norms_[i] = blas.NRM2(blockSize_,RR[i],1);
1314  }
1315  }
1316 
1317  // update the solution
1318  // use R as local work space
1319  // Z may already be in use as work space (holding AX if isSkinny)
1320  {
1321 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1322  Teuchos::TimeMonitor lcltimer( *timerLocalUpdate_ );
1323 #endif
1324  // X <- X*S
1325  MVT::MvAddMv( ONE, *X, ZERO, *X, *R_ );
1326  MVT::MvTimesMatAddMv( ONE, *R_, S, ZERO, *X );
1327  // AX <- AX*S
1328  MVT::MvAddMv( ONE, *AX_, ZERO, *AX_, *R_ );
1329  MVT::MvTimesMatAddMv( ONE, *R_, S, ZERO, *AX_ );
1330  if (hasBOp_) {
1331  // BX <- BX*S
1332  MVT::MvAddMv( ONE, *BX, ZERO, *BX, *R_ );
1333  MVT::MvTimesMatAddMv( ONE, *R_, S, ZERO, *BX );
1334  }
1335  }
1336  }
1337 
1338  // done modifying X,BX; clear temp views and setup const views
1339  X = Teuchos::null;
1340  BX = Teuchos::null;
1341  {
1342  std::vector<int> ind(blockSize_);
1343  for (int i=0; i<blockSize_; ++i) ind[i] = numAuxVecs_+i;
1344  this->X_ = MVT::CloneView(static_cast<const MV&>(*V_),ind);
1345  if (hasBOp_) {
1346  this->BX_ = MVT::CloneView(static_cast<const MV&>(*BV_),ind);
1347  }
1348  else {
1349  this->BX_ = this->X_;
1350  }
1351  }
1352 
1353 
1354  // get objective function value
1355  fx_ = std::accumulate(theta_.begin(),theta_.end(),ZERO);
1356 
1357  // set up R
1358  if (newstate.R != Teuchos::null) {
1359  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetNumberVecs(*newstate.R) < blockSize_,
1360  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): newstate.R must have blockSize number of vectors." );
1361  TEUCHOS_TEST_FOR_EXCEPTION( MVT::GetGlobalLength(*newstate.R) != MVT::GetGlobalLength(*R_),
1362  std::invalid_argument, "Anasazi::RTRBase::initialize(newstate): vector length of newstate.R not correct." );
1363  MVT::SetBlock(*newstate.R,bsind,*R_);
1364  }
1365  else {
1366 #ifdef ANASAZI_TEUCHOS_TIME_MONITOR
1367  Teuchos::TimeMonitor lcltimer( *timerCompRes_ );
1368 #endif
1369  // form R <- AX - BX*T
1370  MVT::MvAddMv(ZERO,*AX_,ONE,*AX_,*R_);
1371  Teuchos::SerialDenseMatrix<int,ScalarType> T(blockSize_,blockSize_);
1372  T.putScalar(ZERO);
1373  for (int i=0; i<blockSize_; i++) T(i,i) = theta_[i];
1374  MVT::MvTimesMatAddMv(-ONE,*BX_,T,ONE,*R_);
1375  }
1376 
1377  // R has been updated; mark the norms as out-of-date
1378  Rnorms_current_ = false;
1379  R2norms_current_ = false;
1380 
1381  // if isSkinny, then AX currently points to Z for temp storage
1382  // set AX back to null
1383  if (isSkinny_) {
1384  AX_ = Teuchos::null;
1385  }
1386 
1387  // finally, we are initialized
1388  initialized_ = true;
1389 
1390  if (om_->isVerbosity( Debug ) ) {
1391  // Check almost everything here
1392  CheckList chk;
1393  chk.checkX = true;
1394  chk.checkAX = true;
1395  chk.checkBX = true;
1396  chk.checkR = true;
1397  chk.checkQ = true;
1398  om_->print( Debug, accuracyCheck(chk, "after initialize()") );
1399  }
1400  }
1401 
1402  template <class ScalarType, class MV, class OP>
1404  {
1406  initialize(empty);
1407  }
1408 
1409 
1410 
1411 
1413  // compute/return residual B-norms
1414  template <class ScalarType, class MV, class OP>
1415  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
1417  if (Rnorms_current_ == false) {
1418  // Update the residual norms
1419  orthman_->norm(*R_,Rnorms_);
1420  Rnorms_current_ = true;
1421  }
1422  return Rnorms_;
1423  }
1424 
1425 
1427  // compute/return residual 2-norms
1428  template <class ScalarType, class MV, class OP>
1429  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
1431  if (R2norms_current_ == false) {
1432  // Update the residual 2-norms
1433  MVT::MvNorm(*R_,R2norms_);
1434  R2norms_current_ = true;
1435  }
1436  return R2norms_;
1437  }
1438 
1439 
1440 
1441 
1443  // Check accuracy, orthogonality, and other debugging stuff
1444  //
1445  // bools specify which tests we want to run (instead of running more than we actually care about)
1446  //
1447  // we don't bother checking the following because they are computed explicitly:
1448  // AH == A*H
1449  //
1450  //
1451  // checkX : X orthonormal
1452  // orthogonal to auxvecs
1453  // checkAX : check AX == A*X
1454  // checkBX : check BX == B*X
1455  // checkEta : check that Eta'*B*X == 0
1456  // orthogonal to auxvecs
1457  // checkAEta : check that AEta = A*Eta
1458  // checkBEta : check that BEta = B*Eta
1459  // checkR : check R orthogonal to X
1460  // checkBR : check R B-orthogonal to X, valid in and immediately after solveTRSubproblem
1461  // checkQ : check that auxiliary vectors are actually orthonormal
1462  //
1463  // TODO:
1464  // add checkTheta
1465  //
1466  template <class ScalarType, class MV, class OP>
1467  std::string RTRBase<ScalarType,MV,OP>::accuracyCheck( const CheckList &chk, const std::string &where ) const
1468  {
1469  using std::setprecision;
1470  using std::scientific;
1471  using std::setw;
1472  using std::endl;
1473  std::stringstream os;
1474  MagnitudeType tmp;
1475 
1476  os << " Debugging checks: " << where << endl;
1477 
1478  // X and friends
1479  if (chk.checkX && initialized_) {
1480  tmp = orthman_->orthonormError(*X_);
1481  os << " >> Error in X^H B X == I : " << scientific << setprecision(10) << tmp << endl;
1482  for (Array_size_type i=0; i<auxVecs_.size(); i++) {
1483  tmp = orthman_->orthogError(*X_,*auxVecs_[i]);
1484  os << " >> Error in X^H B Q[" << i << "] == 0 : " << scientific << setprecision(10) << tmp << endl;
1485  }
1486  }
1487  if (chk.checkAX && !isSkinny_ && initialized_) {
1488  tmp = Utils::errorEquality(*X_, *AX_, AOp_);
1489  os << " >> Error in AX == A*X : " << scientific << setprecision(10) << tmp << endl;
1490  }
1491  if (chk.checkBX && hasBOp_ && initialized_) {
1492  Teuchos::RCP<MV> BX = MVT::Clone(*this->X_,this->blockSize_);
1493  OPT::Apply(*BOp_,*this->X_,*BX);
1494  tmp = Utils::errorEquality(*BX, *BX_);
1495  os << " >> Error in BX == B*X : " << scientific << setprecision(10) << tmp << endl;
1496  }
1497 
1498  // Eta and friends
1499  if (chk.checkEta && initialized_) {
1500  tmp = orthman_->orthogError(*X_,*eta_);
1501  os << " >> Error in X^H B Eta == 0 : " << scientific << setprecision(10) << tmp << endl;
1502  for (Array_size_type i=0; i<auxVecs_.size(); i++) {
1503  tmp = orthman_->orthogError(*eta_,*auxVecs_[i]);
1504  os << " >> Error in Eta^H B Q[" << i << "]==0 : " << scientific << setprecision(10) << tmp << endl;
1505  }
1506  }
1507  if (chk.checkAEta && !isSkinny_ && initialized_) {
1508  tmp = Utils::errorEquality(*eta_, *Aeta_, AOp_);
1509  os << " >> Error in AEta == A*Eta : " << scientific << setprecision(10) << tmp << endl;
1510  }
1511  if (chk.checkBEta && !isSkinny_ && hasBOp_ && initialized_) {
1512  tmp = Utils::errorEquality(*eta_, *Beta_, BOp_);
1513  os << " >> Error in BEta == B*Eta : " << scientific << setprecision(10) << tmp << endl;
1514  }
1515 
1516  // R: this is not B-orthogonality, but standard euclidean orthogonality
1517  if (chk.checkR && initialized_) {
1518  Teuchos::SerialDenseMatrix<int,ScalarType> xTx(blockSize_,blockSize_);
1519  MVT::MvTransMv(ONE,*X_,*R_,xTx);
1520  tmp = xTx.normFrobenius();
1521  os << " >> Error in X^H R == 0 : " << scientific << setprecision(10) << tmp << endl;
1522  }
1523 
1524  // BR: this is B-orthogonality: this is only valid inside and immediately after solveTRSubproblem
1525  // only check if B != I, otherwise, it is equivalent to the above test
1526  if (chk.checkBR && hasBOp_ && initialized_) {
1527  Teuchos::SerialDenseMatrix<int,ScalarType> xTx(blockSize_,blockSize_);
1528  MVT::MvTransMv(ONE,*BX_,*R_,xTx);
1529  tmp = xTx.normFrobenius();
1530  os << " >> Error in X^H B R == 0 : " << scientific << setprecision(10) << tmp << endl;
1531  }
1532 
1533  // Z: Z is preconditioned R, should be on tangent plane
1534  if (chk.checkZ && initialized_) {
1535  Teuchos::SerialDenseMatrix<int,ScalarType> xTx(blockSize_,blockSize_);
1536  MVT::MvTransMv(ONE,*BX_,*Z_,xTx);
1537  tmp = xTx.normFrobenius();
1538  os << " >> Error in X^H B Z == 0 : " << scientific << setprecision(10) << tmp << endl;
1539  }
1540 
1541  // Q
1542  if (chk.checkQ) {
1543  for (Array_size_type i=0; i<auxVecs_.size(); i++) {
1544  tmp = orthman_->orthonormError(*auxVecs_[i]);
1545  os << " >> Error in Q[" << i << "]^H B Q[" << i << "]==I: " << scientific << setprecision(10) << tmp << endl;
1546  for (Array_size_type j=i+1; j<auxVecs_.size(); j++) {
1547  tmp = orthman_->orthogError(*auxVecs_[i],*auxVecs_[j]);
1548  os << " >> Error in Q[" << i << "]^H B Q[" << j << "]==0: " << scientific << setprecision(10) << tmp << endl;
1549  }
1550  }
1551  }
1552  os << endl;
1553  return os.str();
1554  }
1555 
1556 
1558  // Print the current status of the solver
1559  template <class ScalarType, class MV, class OP>
1560  void
1562  {
1563  using std::setprecision;
1564  using std::scientific;
1565  using std::setw;
1566  using std::endl;
1567 
1568  os <<endl;
1569  os <<"================================================================================" << endl;
1570  os << endl;
1571  os <<" RTRBase Solver Status" << endl;
1572  os << endl;
1573  os <<"The solver is "<<(initialized_ ? "initialized." : "not initialized.") << endl;
1574  os <<"The number of iterations performed is " << iter_ << endl;
1575  os <<"The current block size is " << blockSize_ << endl;
1576  os <<"The number of auxiliary vectors is " << numAuxVecs_ << endl;
1577  os <<"The number of operations A*x is " << counterAOp_ << endl;
1578  os <<"The number of operations B*x is " << counterBOp_ << endl;
1579  os <<"The number of operations Prec*x is " << counterPrec_ << endl;
1580  os <<"The most recent rho was " << scientific << setprecision(10) << rho_ << endl;
1581  os <<"The current value of f(x) is " << scientific << setprecision(10) << fx_ << endl;
1582 
1583  if (initialized_) {
1584  os << endl;
1585  os <<"CURRENT EIGENVALUE ESTIMATES "<<endl;
1586  os << setw(20) << "Eigenvalue"
1587  << setw(20) << "Residual(B)"
1588  << setw(20) << "Residual(2)"
1589  << endl;
1590  os <<"--------------------------------------------------------------------------------"<<endl;
1591  for (int i=0; i<blockSize_; i++) {
1592  os << scientific << setprecision(10) << setw(20) << theta_[i];
1593  if (Rnorms_current_) os << scientific << setprecision(10) << setw(20) << Rnorms_[i];
1594  else os << scientific << setprecision(10) << setw(20) << "not current";
1595  if (R2norms_current_) os << scientific << setprecision(10) << setw(20) << R2norms_[i];
1596  else os << scientific << setprecision(10) << setw(20) << "not current";
1597  os << endl;
1598  }
1599  }
1600  os <<"================================================================================" << endl;
1601  os << endl;
1602  }
1603 
1604 
1606  // Inner product 1
1607  template <class ScalarType, class MV, class OP>
1609  RTRBase<ScalarType,MV,OP>::ginner(const MV &xi) const
1610  {
1611  std::vector<MagnitudeType> d(MVT::GetNumberVecs(xi));
1612  MVT::MvNorm(xi,d);
1613  MagnitudeType ret = 0;
1614  for (vecMTiter i=d.begin(); i != d.end(); ++i) {
1615  ret += (*i)*(*i);
1616  }
1617  return ret;
1618  }
1619 
1620 
1622  // Inner product 2
1623  template <class ScalarType, class MV, class OP>
1625  RTRBase<ScalarType,MV,OP>::ginner(const MV &xi, const MV &zeta) const
1626  {
1627  std::vector<ScalarType> d(MVT::GetNumberVecs(xi));
1628  MVT::MvDot(xi,zeta,d);
1629  return SCT::real(std::accumulate(d.begin(),d.end(),SCT::zero()));
1630  }
1631 
1632 
1634  // Inner product 1 without trace accumulation
1635  template <class ScalarType, class MV, class OP>
1636  void RTRBase<ScalarType,MV,OP>::ginnersep(
1637  const MV &xi,
1638  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> &d) const
1639  {
1640  MVT::MvNorm(xi,d);
1641  for (vecMTiter i=d.begin(); i != d.end(); ++i) {
1642  (*i) = (*i)*(*i);
1643  }
1644  }
1645 
1646 
1648  // Inner product 2 without trace accumulation
1649  template <class ScalarType, class MV, class OP>
1650  void RTRBase<ScalarType,MV,OP>::ginnersep(
1651  const MV &xi, const MV &zeta,
1652  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType> &d) const
1653  {
1654  std::vector<ScalarType> dC(MVT::GetNumberVecs(xi));
1655  MVT::MvDot(xi,zeta,dC);
1656  vecMTiter iR=d.begin();
1657  vecSTiter iS=dC.begin();
1658  for (; iR != d.end(); iR++, iS++) {
1659  (*iR) = SCT::real(*iS);
1660  }
1661  }
1662 
1663  template <class ScalarType, class MV, class OP>
1665  return auxVecs_;
1666  }
1667 
1668  template <class ScalarType, class MV, class OP>
1670  return(blockSize_);
1671  }
1672 
1673  template <class ScalarType, class MV, class OP>
1675  return(*problem_);
1676  }
1677 
1678  template <class ScalarType, class MV, class OP>
1680  return blockSize_;
1681  }
1682 
1683  template <class ScalarType, class MV, class OP>
1685  {
1686  if (!initialized_) return 0;
1687  return nevLocal_;
1688  }
1689 
1690  template <class ScalarType, class MV, class OP>
1691  std::vector<typename Teuchos::ScalarTraits<ScalarType>::magnitudeType>
1693  {
1694  std::vector<MagnitudeType> ret = ritz2norms_;
1695  ret.resize(nevLocal_);
1696  return ret;
1697  }
1698 
1699  template <class ScalarType, class MV, class OP>
1700  std::vector<Value<ScalarType> >
1702  {
1703  std::vector<Value<ScalarType> > ret(nevLocal_);
1704  for (int i=0; i<nevLocal_; i++) {
1705  ret[i].realpart = theta_[i];
1706  ret[i].imagpart = ZERO;
1707  }
1708  return ret;
1709  }
1710 
1711  template <class ScalarType, class MV, class OP>
1714  {
1715  return X_;
1716  }
1717 
1718  template <class ScalarType, class MV, class OP>
1720  {
1721  iter_=0;
1722  }
1723 
1724  template <class ScalarType, class MV, class OP>
1726  {
1727  return iter_;
1728  }
1729 
1730  template <class ScalarType, class MV, class OP>
1732  {
1734  state.X = X_;
1735  state.AX = AX_;
1736  if (hasBOp_) {
1737  state.BX = BX_;
1738  }
1739  else {
1740  state.BX = Teuchos::null;
1741  }
1742  state.rho = rho_;
1743  state.R = R_;
1744  state.T = Teuchos::rcp(new std::vector<MagnitudeType>(theta_));
1745  return state;
1746  }
1747 
1748  template <class ScalarType, class MV, class OP>
1750  {
1751  return initialized_;
1752  }
1753 
1754  template <class ScalarType, class MV, class OP>
1756  {
1757  std::vector<int> ret(nevLocal_,0);
1758  return ret;
1759  }
1760 
1761 
1762 } // end Anasazi namespace
1763 
1764 #endif // ANASAZI_RTRBASE_HPP
RTROrthoFailure is thrown when an orthogonalization attempt fails.
int getBlockSize() const
Get the blocksize to be used by the iterative solver in solving this eigenproblem.
int getCurSubspaceDim() const
Get the dimension of the search subspace used to generate the current eigenvectors and eigenvalues...
void resetNumIters()
Reset the iteration count.
std::vector< Value< ScalarType > > getRitzValues()
Get the Ritz values from the previous iteration.
Structure to contain pointers to RTR state variables.
virtual void currentStatus(std::ostream &os)
This method requests that the solver print out its current status to screen.
Teuchos::RCP< const MV > AX
The image of the current eigenvectors under A, or Teuchos::null is we implement a skinny solver...
This class defines the interface required by an eigensolver and status test class to compute solution...
Teuchos::RCP< const MV > BX
The image of the current eigenvectors under B, or Teuchos::null if B was not specified.
Declaration of basic traits for the multivector type.
T & get(const std::string &name, T def_value)
int multiply(ETransp transa, ETransp transb, ScalarType alpha, const SerialDenseMatrix< OrdinalType, ScalarType > &A, const SerialDenseMatrix< OrdinalType, ScalarType > &B, ScalarType beta)
Teuchos::RCP< const MV > getRitzVectors()
Get the Ritz vectors from the previous iteration.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getRitzRes2Norms()
Get the 2-norms of the Ritz residuals.
This class is an abstract base class for Implicit Riemannian Trust-Region based eigensolvers. The class provides the interfaces shared by the IRTR solvers (e.g., getState() and initialize()) as well as the shared implementations (e.g., inner products).
Virtual base class which defines basic traits for the operator type.
int getMaxSubspaceDim() const
Get the maximum dimension allocated for the search subspace. For RTR, this always returns getBlockSiz...
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getRes2Norms()
Get the current residual 2-norms.
Templated virtual class for providing orthogonalization/orthonormalization methods with matrix-based ...
virtual void iterate()=0
This method performs RTR iterations until the status test indicates the need to stop or an error occu...
ScalarTraits< ScalarType >::magnitudeType NRM2(const OrdinalType &n, const ScalarType *x, const OrdinalType &incx) const
Teuchos::RCP< const MV > X
The current eigenvectors.
An exception class parent to all Anasazi exceptions.
void setStatusTest(Teuchos::RCP< StatusTest< ScalarType, MV, OP > > test)
Set a new StatusTest for the solver.
Pure virtual base class which describes the basic interface to the iterative eigensolver.
virtual ~RTRBase()
RTRBase destructor
Anasazi&#39;s templated, static class providing utilities for the solvers.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Output managers remove the need for the eigensolver to know any information about the required output...
void setBlockSize(int blockSize)
Set the blocksize to be used by the iterative solver in solving this eigenproblem.
int putScalar(const ScalarType value=Teuchos::ScalarTraits< ScalarType >::zero())
std::vector< int > getRitzIndex()
Get the index used for extracting Ritz vectors from getRitzVectors().
Traits class which defines basic operations on multivectors.
Virtual base class which defines basic traits for the operator type.
Teuchos::Array< Teuchos::RCP< const MV > > getAuxVecs() const
Get the current auxiliary vectors.
bool isInitialized() const
Indicates whether the solver has been initialized or not.
void resize(size_type new_size, const value_type &x=value_type())
RTRInitFailure is thrown when the RTR solver is unable to generate an initial iterate in the RTRBase:...
Teuchos::ScalarTraits< ScalarType >::magnitudeType rho
The current rho value. This is only valid if the debugging level of verbosity is enabled.
std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > getResNorms()
Get the current residual norms.
int getNumIters() const
Get the current iteration count.
Teuchos::RCP< const std::vector< typename Teuchos::ScalarTraits< ScalarType >::magnitudeType > > T
The current Ritz values.
Teuchos::RCP< const MV > R
The current residual vectors.
void setAuxVecs(const Teuchos::Array< Teuchos::RCP< const MV > > &auxvecs)
Set the auxiliary vectors for the solver.
void SCAL(const OrdinalType &n, const ScalarType &alpha, ScalarType *x, const OrdinalType &incx) const
Types and exceptions used within Anasazi solvers and interfaces.
void initialize()
Initialize the solver with the initial vectors from the eigenproblem or random data.
const Eigenproblem< ScalarType, MV, OP > & getProblem() const
Get a constant reference to the eigenvalue problem.
Anasazi&#39;s templated pure virtual class for managing the sorting of approximate eigenvalues computed b...
RTRState< ScalarType, MV > getState() const
Get the current state of the eigensolver.
Common interface of stopping criteria for Anasazi&#39;s solvers.
RTRRitzFailure is thrown when the RTR solver is unable to continue a call to RTRBase::iterate() due t...
SerialDenseMatrix< OrdinalType, ScalarType > & assign(const SerialDenseMatrix< OrdinalType, ScalarType > &Source)
The Eigensolver is a templated virtual base class that defines the basic interface that any eigensolv...
RTRBase(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< GenOrthoManager< ScalarType, MV, OP > > &ortho, Teuchos::ParameterList &params, const std::string &solverLabel, bool skinnySolver)
RTRBase constructor with eigenproblem, solver utilities, and parameter list of solver options...
Teuchos::RCP< StatusTest< ScalarType, MV, OP > > getStatusTest() const
Get the current StatusTest used by the solver.
Class which provides internal utilities for the Anasazi solvers.